pi-DoMUS: Parallel Deal.II MUltiphysics Solver
non_conservative.h
Go to the documentation of this file.
1 /*
2  * Non Conservative Interface
3  *
4  * This class is a child of interface.h
5  *
6  * It is used to define a problem in the case whole the problem
7  * could not be stated in terms of energy.
8  *
9  */
10 
11 #ifndef _non_conservative_interface_h
12 #define _non_conservative_interface_h
13 
14 #include "interfaces/interface.h"
15 
16 using namespace deal2lkit;
17 
18 template<int dim, int spacedim, int n_components, class Implementation, typename LAC=LATrilinos>
19 class NonConservativeInterface : public Interface<dim,spacedim,n_components,LAC>
20 {
21 
22  typedef FEValuesCache<dim,spacedim> Scratch;
23  typedef Assembly::CopyData::piDoMUSPreconditioner<dim,spacedim> CopyPreconditioner;
24  typedef Assembly::CopyData::piDoMUSSystem<dim,spacedim> CopySystem;
25 public:
26 
28 
29  NonConservativeInterface(const std::string &name="",
30  const std::string &default_fe="FE_Q(1)",
31  const std::string &default_component_names="u",
32  const std::string &default_coupling="",
33  const std::string &default_preconditioner_coupling="",
34  const std::string &default_differential_components="") :
35  Interface<dim,spacedim,n_components,LAC>(name, default_fe, default_component_names,
36  default_coupling, default_preconditioner_coupling,
37  default_differential_components) {};
38 
40  {
42  }
44  {
46  }
47 
49  Scratch &scratch,
50  CopySystem &data,
51  std::vector<double> &local_residual) const
52  {
53  static_cast<const Implementation *>(this)->system_residual(cell, scratch, data, local_residual);
54  }
55 
57  Scratch &scratch,
58  CopySystem &data,
59  std::vector<Sdouble> &local_residual) const
60  {
61  static_cast<const Implementation *>(this)->system_residual(cell, scratch, data, local_residual);
62  }
63 
65  Scratch &scratch,
66  CopyPreconditioner &data,
67  std::vector<Sdouble> &local_residual) const
68  {
69  static_cast<const Implementation *>(this)->preconditioner_residual(cell, scratch, data, local_residual);
70  }
71 
72 // aux matrices ////////////////////////////////////////////////////////////////
73 
75  Scratch &scratch,
76  CopyPreconditioner &data,
77  std::vector<std::vector<double> > &local_residuals) const
78  {
79  static_cast<const Implementation *>(this)->aux_matrix_residuals(cell, scratch, data, local_residuals);
80  }
81 
83  Scratch &scratch,
84  CopyPreconditioner &data,
85  std::vector<std::vector<Sdouble> > &local_residuals) const
86  {
87  static_cast<const Implementation *>(this)->aux_matrix_residuals(cell, scratch, data, local_residuals);
88  }
89 };
90 
91 #endif
virtual void get_system_residual(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< double > &local_residual) const
Definition: non_conservative.h:48
virtual void get_preconditioner_residual(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< Sdouble > &local_residual) const
Definition: non_conservative.h:64
Definition: lac_type.h:33
NonConservativeInterface(const std::string &name="", const std::string &default_fe="FE_Q(1)", const std::string &default_component_names="u", const std::string &default_coupling="", const std::string &default_preconditioner_coupling="", const std::string &default_differential_components="")
Definition: non_conservative.h:29
virtual void get_aux_matrix_residuals(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< std::vector< double > > &local_residuals) const
Definition: non_conservative.h:74
ActiveSelector::active_cell_iterator active_cell_iterator
Interface.
Definition: old_interface.h:56
virtual ~NonConservativeInterface()
Definition: non_conservative.h:27
virtual void parse_parameters_call_back()
virtual void get_system_residual(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< Sdouble > &local_residual) const
Definition: non_conservative.h:56
virtual void parse_parameters_call_back()
Definition: non_conservative.h:43
virtual void declare_parameters(ParameterHandler &prm)
Definition: non_conservative.h:39
virtual void get_aux_matrix_residuals(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< std::vector< Sdouble > > &local_residuals) const
Definition: non_conservative.h:82
Definition: non_conservative.h:19
virtual void declare_parameters(ParameterHandler &prm)