GRINS-0.6.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-2015 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"
35 
36 // libMesh
37 #include "libmesh/fem_system.h"
38 
39 #ifdef GRINS_HAVE_GRVY
40 // GRVY timers
41 #include "grvy.h"
42 #endif
43 
44 // libMesh forward declartions
45 class GetPot;
46 
47 namespace libMesh
48 {
49  class EquationSystems;
50  class DiffContext;
51 
52  template <typename Scalar>
54 }
55 
56 namespace GRINS
57 {
58  // Forward Declarations
59  template <typename Scalar>
61 
63 
80  //TODO: is it F(u) or F(u_{\theta})?
81  class MultiphysicsSystem : public libMesh::FEMSystem
82  {
83  public:
84 
86  MultiphysicsSystem( libMesh::EquationSystems& es,
87  const std::string& name,
88  const unsigned int number );
89 
92 
94  void attach_physics_list( PhysicsList physics_list );
95 
97 
102  virtual void read_input_options( const GetPot& input );
103 
105  virtual void init_data();
106 
108  void register_postprocessing_vars( const GetPot& input,
109  PostProcessedQuantities<libMesh::Real>& postprocessing );
110 
112  // named in this call.
113  void register_parameter
114  ( const std::string & param_name,
116 
118  virtual libMesh::AutoPtr<libMesh::DiffContext> build_context();
119 
121  virtual void init_context( libMesh::DiffContext &context );
122 
123  // residual and jacobian calculations
124  // element_*, side_* as *time_derivative, *constraint, *mass_residual
125 
127  virtual bool element_time_derivative( bool request_jacobian,
128  libMesh::DiffContext& context );
129 
131  virtual bool side_time_derivative( bool request_jacobian,
132  libMesh::DiffContext& context );
133 
135  virtual bool nonlocal_time_derivative( bool request_jacobian,
136  libMesh::DiffContext& context );
137 
140  virtual bool element_constraint( bool request_jacobian,
141  libMesh::DiffContext& context );
142 
144  virtual bool side_constraint( bool request_jacobian,
145  libMesh::DiffContext& context );
146 
148  virtual bool nonlocal_constraint( bool request_jacobian,
149  libMesh::DiffContext& context );
150 
152  virtual bool mass_residual( bool request_jacobian,
153  libMesh::DiffContext& context );
154 
156  virtual bool nonlocal_mass_residual( bool request_jacobian,
157  libMesh::DiffContext& context );
158 
160  bool has_physics( const std::string physics_name ) const;
161 
162  std::tr1::shared_ptr<GRINS::Physics> get_physics( const std::string physics_name );
163 
164  std::tr1::shared_ptr<GRINS::Physics> get_physics( const std::string physics_name ) const;
165 
166  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
167  const AssemblyContext& context,
168  const libMesh::Point& point,
169  libMesh::Real& value );
170 
171 #ifdef GRINS_USE_GRVY_TIMERS
172  void attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer );
174 #endif
175 
176  private:
177 
179 
182 
184 
185 #ifdef GRINS_USE_GRVY_TIMERS
186  GRVY::GRVY_Timer_Class* _timer;
187 #endif
188 
189  // Useful typedef for refactoring
191  typedef void (GRINS::Physics::*CacheFuncType) (const AssemblyContext&, CachedValues &);
192 
193  // Refactored residual evaluation implementation
194  bool _general_residual( bool request_jacobian,
195  libMesh::DiffContext& context,
196  ResFuncType resfunc,
197  CacheFuncType cachefunc);
198  };
199 
200  inline
201  std::tr1::shared_ptr<GRINS::Physics> MultiphysicsSystem::get_physics( const std::string physics_name ) const
202  {
203  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
204 
205  return _physics_list.find(physics_name)->second;
206  }
207 
208 } //End namespace block
209 
210 #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.
void(GRINS::Physics::* CacheFuncType)(const AssemblyContext &, CachedValues &)
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
Definition: var_typedefs.h:57
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
std::tr1::shared_ptr< GRINS::Physics > get_physics(const std::string physics_name)
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.
void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer)
Each Physics will register its copy(s) of an independent variable.
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.
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.
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::AutoPtr< 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.
Interface with libMesh for solving Multiphysics problems.
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 .
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()
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled.

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1