GRINS-0.7.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-2016 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 #ifdef GRINS_HAVE_GRVY
41 // GRVY timers
42 #include "grvy.h"
43 #endif
44 
45 // libMesh forward declartions
46 class GetPot;
47 
48 namespace libMesh
49 {
50  class EquationSystems;
51  class DiffContext;
52 
53  template <typename Scalar>
55 }
56 
57 namespace GRINS
58 {
59  // Forward Declarations
60  template <typename Scalar>
62 
64 
81  //TODO: is it F(u) or F(u_{\theta})?
82  class MultiphysicsSystem : public libMesh::FEMSystem
83  {
84  public:
85 
87  MultiphysicsSystem( libMesh::EquationSystems& es,
88  const std::string& name,
89  const unsigned int number );
90 
93 
95  void attach_physics_list( PhysicsList physics_list );
96 
98  { return _physics_list; }
99 
101 
106  virtual void read_input_options( const GetPot& input );
107 
109  virtual void init_data();
110 
112  void register_postprocessing_vars( const GetPot& input,
113  PostProcessedQuantities<libMesh::Real>& postprocessing );
114 
116  // named in this call.
117  void register_parameter
118  ( const std::string & param_name,
120 
122  virtual libMesh::UniquePtr<libMesh::DiffContext> build_context();
123 
125  virtual void init_context( libMesh::DiffContext &context );
126 
127  // residual and jacobian calculations
128  // element_*, side_* as *time_derivative, *constraint, *mass_residual
129 
131  virtual bool element_time_derivative( bool request_jacobian,
132  libMesh::DiffContext& context );
133 
135  virtual bool side_time_derivative( bool request_jacobian,
136  libMesh::DiffContext& context );
137 
139  virtual bool nonlocal_time_derivative( bool request_jacobian,
140  libMesh::DiffContext& context );
141 
144  virtual bool element_constraint( bool request_jacobian,
145  libMesh::DiffContext& context );
146 
148  virtual bool side_constraint( bool request_jacobian,
149  libMesh::DiffContext& context );
150 
152  virtual bool nonlocal_constraint( bool request_jacobian,
153  libMesh::DiffContext& context );
154 
156  virtual bool damping_residual( bool request_jacobian,
157  libMesh::DiffContext& context );
158 
160  virtual bool mass_residual( bool request_jacobian,
161  libMesh::DiffContext& context );
162 
164  virtual bool nonlocal_mass_residual( bool request_jacobian,
165  libMesh::DiffContext& context );
166 
168  bool has_physics( const std::string physics_name ) const;
169 
170  SharedPtr<GRINS::Physics> get_physics( const std::string physics_name );
171 
172  SharedPtr<GRINS::Physics> get_physics( const std::string physics_name ) const;
173 
174  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
175  const AssemblyContext& context,
176  const libMesh::Point& point,
177  libMesh::Real& value );
178 
179  std::vector<SharedPtr<NeumannBCContainer> >& get_neumann_bcs()
180  { return _neumann_bcs; }
181 
182  const std::vector<SharedPtr<NeumannBCContainer> >& get_neumann_bcs() const
183  { return _neumann_bcs; }
184 
185 #ifdef GRINS_USE_GRVY_TIMERS
186  void attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer );
188 #endif
189 
190  private:
191 
193 
196 
198 
199  // A list of names of variables who need their own numerical
200  // jacobian deltas
201  std::vector<std::string> _numerical_jacobian_h_variables;
202 
203  // A list of values for per-variable numerical jacobian deltas
204  std::vector<libMesh::Real> _numerical_jacobian_h_values;
205 
207 
210  const GetPot* _input;
211 
213 
217  std::vector<SharedPtr<NeumannBCContainer> > _neumann_bcs;
218 
219 #ifdef GRINS_USE_GRVY_TIMERS
220  GRVY::GRVY_Timer_Class* _timer;
221 #endif
222 
223  // Useful typedef for refactoring
225  typedef void (GRINS::Physics::*CacheFuncType) (const AssemblyContext&, CachedValues &);
226 
227  // Refactored residual evaluation implementation
228  bool _general_residual( bool request_jacobian,
229  libMesh::DiffContext& context,
230  ResFuncType resfunc,
231  CacheFuncType cachefunc);
232 
235  const std::vector<SharedPtr<NeumannBCContainer> >& neumann_bcs,
236  std::vector<SharedPtr<NeumannBCContainer> >& active_neumann_bcs );
237 
239  bool apply_neumann_bcs( bool request_jacobian,
240  libMesh::DiffContext& context );
241  };
242 
243  inline
244  SharedPtr<GRINS::Physics> MultiphysicsSystem::get_physics( const std::string physics_name ) const
245  {
246  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
247 
248  return _physics_list.find(physics_name)->second;
249  }
250 
251 } //End namespace block
252 
253 #endif // GRINS_MULTIPHYSICS_SYS_H
void(GRINS::Physics::* ResFuncType)(bool, AssemblyContext &, CachedValues &)
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(GRINS::Physics::* CacheFuncType)(const AssemblyContext &, CachedValues &)
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:107
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)
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.
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Element interior contributions to which have time varying components.
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
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.
const GetPot * _input
Cached for helping build boundary conditions.

Generated on Thu Jun 2 2016 21:52:27 for GRINS-0.7.0 by  doxygen 1.8.10