GRINS-0.8.0
multiphysics_sys.h
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 #ifndef GRINS_MULTIPHYSICS_SYS_H
27 #define GRINS_MULTIPHYSICS_SYS_H
28 
29 // C++
30 #include <string>
31 
32 // GRINS
33 #include "grins_config.h"
34 #include "grins/physics.h"
36 
37 // libMesh
38 #include "libmesh/fem_system.h"
39 
40 // libMesh forward declartions
41 class GetPot;
42 
43 namespace libMesh
44 {
45  class EquationSystems;
46  class DiffContext;
47 
48  template <typename Scalar>
50 }
51 
52 namespace GRINS
53 {
54  // Forward Declarations
55  template <typename Scalar>
57 
59 
76  //TODO: is it F(u) or F(u_{\theta})?
77  class MultiphysicsSystem : public libMesh::FEMSystem
78  {
79  public:
80 
82  MultiphysicsSystem( libMesh::EquationSystems& es,
83  const std::string& name,
84  const unsigned int number );
85 
88 
90  void attach_physics_list( PhysicsList physics_list );
91 
93  { return _physics_list; }
94 
96 
101  virtual void read_input_options( const GetPot& input );
102 
104  virtual void init_data();
105 
107  void register_postprocessing_vars( const GetPot& input,
108  PostProcessedQuantities<libMesh::Real>& postprocessing );
109 
111  // named in this call.
112  void register_parameter
113  ( const std::string & param_name,
115 
117  virtual libMesh::UniquePtr<libMesh::DiffContext> build_context();
118 
120  virtual void init_context( libMesh::DiffContext &context );
121 
123 
124  virtual void assembly( bool get_residual,
125  bool get_jacobian,
126  bool apply_heterogeneous_constraints = false,
127  bool apply_no_constraints = false );
128 
130 
132  virtual void reinit();
133 
134  // residual and jacobian calculations
135  // element_*, side_* as *time_derivative, *constraint, *mass_residual
136 
138  virtual bool element_time_derivative( bool request_jacobian,
139  libMesh::DiffContext& context );
140 
142  virtual bool side_time_derivative( bool request_jacobian,
143  libMesh::DiffContext& context );
144 
146  virtual bool nonlocal_time_derivative( bool request_jacobian,
147  libMesh::DiffContext& context );
148 
151  virtual bool element_constraint( bool request_jacobian,
152  libMesh::DiffContext& context );
153 
155  virtual bool side_constraint( bool request_jacobian,
156  libMesh::DiffContext& context );
157 
159  virtual bool nonlocal_constraint( bool request_jacobian,
160  libMesh::DiffContext& context );
161 
163  virtual bool damping_residual( bool request_jacobian,
164  libMesh::DiffContext& context );
165 
167  virtual bool mass_residual( bool request_jacobian,
168  libMesh::DiffContext& context );
169 
171  virtual bool nonlocal_mass_residual( bool request_jacobian,
172  libMesh::DiffContext& context );
173 
175  bool has_physics( const std::string physics_name ) const;
176 
177  SharedPtr<GRINS::Physics> get_physics( const std::string physics_name );
178 
179  SharedPtr<GRINS::Physics> get_physics( const std::string physics_name ) const;
180 
181  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
182  const AssemblyContext& context,
183  const libMesh::Point& point,
184  libMesh::Real& value );
185 
186  std::vector<SharedPtr<NeumannBCContainer> >& get_neumann_bcs()
187  { return _neumann_bcs; }
188 
189  const std::vector<SharedPtr<NeumannBCContainer> >& get_neumann_bcs() const
190  { return _neumann_bcs; }
191 
192  private:
193 
195 
198 
200 
201  // A list of names of variables who need their own numerical
202  // jacobian deltas
203  std::vector<std::string> _numerical_jacobian_h_variables;
204 
205  // A list of values for per-variable numerical jacobian deltas
206  std::vector<libMesh::Real> _numerical_jacobian_h_values;
207 
209 
212  const GetPot* _input;
213 
215 
219  std::vector<SharedPtr<NeumannBCContainer> > _neumann_bcs;
220 
222  libMesh::UniquePtr<libMesh::System::Constraint> _constraint;
223 
224  // Useful typedef to pointer-to-member functions so we can call all
225  // residual and caching functions using a single function (_general_residual)
226  typedef void (GRINS::Physics::*ResFuncType) (bool, AssemblyContext &);
228 
229  // Refactored residual evaluation implementation
230  bool _general_residual( bool request_jacobian,
231  libMesh::DiffContext & context,
232  ResFuncType resfunc,
233  CacheFuncType cachefunc);
234 
237  const std::vector<SharedPtr<NeumannBCContainer> >& neumann_bcs,
238  std::vector<SharedPtr<NeumannBCContainer> >& active_neumann_bcs );
239 
241  bool apply_neumann_bcs( bool request_jacobian,
242  libMesh::DiffContext& context );
243  };
244 
245  inline
246  SharedPtr<GRINS::Physics> MultiphysicsSystem::get_physics( const std::string physics_name ) const
247  {
248  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
249 
250  return _physics_list.find(physics_name)->second;
251  }
252 
253 } //End namespace block
254 
255 #endif // GRINS_MULTIPHYSICS_SYS_H
bool has_physics(const std::string physics_name) const
Query to check if a particular physics has been enabled.
const std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs() const
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
std::vector< std::string > _numerical_jacobian_h_variables
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
virtual void init_data()
System initialization. Calls each physics implementation of init_variables()
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
virtual bool side_constraint(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to which do not have time varying components.
libMesh::boundary_id_type BoundaryID
More descriptive name of the type used for boundary ids.
Definition: var_typedefs.h:56
virtual bool damping_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to .
SharedPtr< GRINS::Physics > get_physics(const std::string physics_name)
libMesh::UniquePtr< libMesh::System::Constraint > _constraint
Constraint application object.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
virtual bool nonlocal_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables which have time varying components.
std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs()
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
GRINS namespace.
~MultiphysicsSystem()
Destructor. Clean up all physics allocations.
void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer)
Each Physics will register its copy(s) of an independent variable.
virtual bool nonlocal_mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables.
void(GRINS::Physics::* ResFuncType)(bool, AssemblyContext &)
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Element interior contributions to which have time varying components.
virtual void reinit()
Override FEMSystem::reinit.
virtual libMesh::UniquePtr< libMesh::DiffContext > build_context()
Override FEMSystem::build_context in order to use our own AssemblyContext.
virtual bool nonlocal_constraint(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables which do not have time varying components.
void get_active_neumann_bcs(BoundaryID bc_id, const std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs, std::vector< SharedPtr< NeumannBCContainer > > &active_neumann_bcs)
Extract the bcs from neumann_bcs that are active on bc_id and return them in active_neumann_bcs.
Interface with libMesh for solving Multiphysics problems.
std::map< std::string, SharedPtr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
Definition: var_typedefs.h:59
std::vector< libMesh::Real > _numerical_jacobian_h_values
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false)
Override FEMSystem::assembly.
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual bool element_constraint(bool request_jacobian, libMesh::DiffContext &context)
virtual bool mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to .
bool apply_neumann_bcs(bool request_jacobian, libMesh::DiffContext &context)
Applies the subset of _neumann_bcs that are active on the current element side.
MultiphysicsSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Constructor. Will be called by libMesh only.
virtual bool side_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to which have time varying components.
virtual void init_context(libMesh::DiffContext &context)
Context initialization. Calls each physics implementation of init_context()
const PhysicsList & get_physics_list() const
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled.
void(GRINS::Physics::* CacheFuncType)(AssemblyContext &)
const GetPot * _input
Cached for helping build boundary conditions.

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1