GRINS-0.8.0
physics.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 #ifndef GRINS_PHYSICS_H
26 #define GRINS_PHYSICS_H
27 
28 // C++
29 #include <string>
30 #include <set>
31 
32 //GRINS
33 #include "grins_config.h"
35 #include "grins/var_typedefs.h"
36 #include "grins/physics_naming.h"
37 #include "grins/cached_values.h"
38 #include "grins/parameter_user.h"
39 #include "grins/assembly_context.h"
40 
41 //libMesh
42 #include "libmesh/libmesh.h"
43 #include "libmesh/point.h"
44 #include "libmesh/fe_base.h"
45 #include "libmesh/system.h"
46 #include "libmesh/mesh_base.h"
47 #include "libmesh/auto_ptr.h"
48 
49 // libMesh forward declarations
50 class GetPot;
51 namespace libMesh
52 {
53  template <typename Scalar>
54  class CompositeFunction;
55 
56  class FEMSystem;
57  class Elem;
58 
59  template <typename Scalar>
60  class ParameterMultiAccessor;
61 }
62 
64 namespace GRINS
65 {
66  // GRINS forward declarations
67  class BCHandlingBase;
68  class ICHandlingBase;
69  class NBCContainer;
70  class DBCContainer;
71  class MultiphysicsSystem;
72  class FEVariablesBase;
73 
74  template <typename Scalar>
75  class PostProcessedQuantities;
76 
78 
98  //TODO: is it F(u) or F(u_{\theta})?
99 
100  // element* routines work on element interiors
101  // side* routines work on element sides
102 
103  // *_time_derivative correspond to calculating terms for F(u)
104  // *_mass_residual correspond to calculating terms for M(u)\dot{u}
105 
106  class Physics : public ParameterUser
107  {
108 
109  public:
110 
111  Physics( const GRINS::PhysicsName& physics_name, const GetPot& input );
112  virtual ~Physics();
113 
115  virtual void init_variables( libMesh::FEMSystem* /*system*/ ){};
116 
118  virtual bool enabled_on_elem( const libMesh::Elem* elem );
119 
121 
123  void set_is_steady( bool is_steady );
124 
126  bool is_steady() const;
127 
131 
132  static bool is_axisymmetric()
133  { return _is_axisymmetric; }
134 
136 
142  virtual void set_time_evolving_vars( libMesh::FEMSystem* system );
143 
145 
147  virtual void auxiliary_init( MultiphysicsSystem& system );
148 
150 
155  virtual void register_postprocessing_vars( const GetPot& input,
156  PostProcessedQuantities<libMesh::Real>& postprocessing );
157 
159  virtual void init_context( AssemblyContext& context );
160 
162  virtual void preassembly( MultiphysicsSystem & /*system*/ ){};
163 
165 
168  virtual void reinit( MultiphysicsSystem & /*system*/ ){};
169 
170  // residual and jacobian calculations
171  // element_*, side_* as *time_derivative, *constraint, *mass_residual
172 
174  virtual void element_time_derivative( bool /*compute_jacobian*/,
175  AssemblyContext & /*context*/ ){}
176 
178  virtual void side_time_derivative( bool /*compute_jacobian*/,
179  AssemblyContext & /*context*/ ){}
180 
182  virtual void nonlocal_time_derivative( bool /*compute_jacobian*/,
183  AssemblyContext & /*context*/ ){}
184 
186  virtual void element_constraint( bool /*compute_jacobian*/,
187  AssemblyContext & /*context*/ ){}
188 
190  virtual void side_constraint( bool /*compute_jacobian*/,
191  AssemblyContext & /*context*/ ){}
192 
194  virtual void nonlocal_constraint( bool /*compute_jacobian*/,
195  AssemblyContext & /*context*/ ){}
196 
198  virtual void damping_residual( bool /*compute_jacobian*/,
199  AssemblyContext & /*context*/ ){}
200 
202  virtual void mass_residual( bool /*compute_jacobian*/,
203  AssemblyContext & /*context*/ ){}
204 
206  virtual void nonlocal_mass_residual( bool /*compute_jacobian*/,
207  AssemblyContext & /*context*/ ){}
208 
209  void init_ics( libMesh::FEMSystem* system,
211 
213 
214  virtual void compute_side_time_derivative_cache( AssemblyContext & /*context*/ ){}
215 
217 
218  virtual void compute_element_constraint_cache( AssemblyContext & /*context*/ ){}
219 
220  virtual void compute_side_constraint_cache( AssemblyContext & /*context*/ ){}
221 
222  virtual void compute_nonlocal_constraint_cache( AssemblyContext & /*context*/ ){}
223 
224  virtual void compute_damping_residual_cache( AssemblyContext & /*context*/ ){}
225 
226  virtual void compute_mass_residual_cache( AssemblyContext & /*context*/ ){}
227 
229 
230  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
231  const AssemblyContext& context,
232  const libMesh::Point& point,
233  libMesh::Real& value );
234 
236 
237  protected:
238 
240  libMesh::UniquePtr<libMesh::FEGenericBase<libMesh::Real> > build_new_fe( const libMesh::Elem* elem,
241  const libMesh::FEGenericBase<libMesh::Real>* fe,
242  const libMesh::Point p );
243 
244  void parse_enabled_subdomains( const GetPot& input,
245  const std::string& physics_name );
246 
248  void check_var_subdomain_consistency( const FEVariablesBase& var ) const;
249 
251 
257 
259 
261  std::set<libMesh::subdomain_id_type> _enabled_subdomains;
262 
264 
266  static bool _is_steady;
267 
269  static bool _is_axisymmetric;
270 
271  private:
272  Physics();
273 
274  }; // End Physics class declarations
275 
276  /* ------------------------- Inline Functions -------------------------*/
277  inline
279  {
280  return _ic_handler;
281  }
282 
283 } // End namespace GRINS
284 
285 #endif //GRINS_PHYSICS_H
virtual void compute_element_constraint_cache(AssemblyContext &)
Definition: physics.h:218
static bool is_axisymmetric()
Definition: physics.h:132
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
Definition: physics.C:135
virtual void compute_mass_residual_cache(AssemblyContext &)
Definition: physics.h:226
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:97
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
virtual void compute_element_time_derivative_cache(AssemblyContext &)
Definition: physics.h:212
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register name of postprocessed quantity with PostProcessedQuantities.
Definition: physics.C:129
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe(const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
Definition: physics.C:143
virtual void compute_damping_residual_cache(AssemblyContext &)
Definition: physics.h:224
virtual void reinit(MultiphysicsSystem &)
Any reinitialization that needs to be done.
Definition: physics.h:168
virtual void side_time_derivative(bool, AssemblyContext &)
Time dependent part(s) of physics for boundaries of elements on the domain boundary.
Definition: physics.h:178
void check_var_subdomain_consistency(const FEVariablesBase &var) const
Check that var is enabled on at least the subdomains this Physics is.
Definition: physics.C:174
ICHandlingBase * get_ic_handler()
Definition: physics.h:278
virtual void element_time_derivative(bool, AssemblyContext &)
Time dependent part(s) of physics for element interiors.
Definition: physics.h:174
virtual void nonlocal_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for scalar variables.
Definition: physics.h:194
virtual ~Physics()
Definition: physics.C:60
virtual void compute_side_constraint_cache(AssemblyContext &)
Definition: physics.h:220
GRINS namespace.
virtual void nonlocal_mass_residual(bool, AssemblyContext &)
Mass matrix part(s) for scalar variables.
Definition: physics.h:206
virtual void side_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for boundaries of elements on the domain boundary.
Definition: physics.h:190
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:256
ParameterUser base class. Utility methods for subclasses.
virtual bool enabled_on_elem(const libMesh::Elem *elem)
Find if current physics is active on supplied element.
Definition: physics.C:77
virtual void element_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for element interiors.
Definition: physics.h:186
std::string PhysicsName
virtual void compute_nonlocal_time_derivative_cache(AssemblyContext &)
Definition: physics.h:216
virtual void mass_residual(bool, AssemblyContext &)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Definition: physics.h:202
virtual void auxiliary_init(MultiphysicsSystem &system)
Any auxillary initialization a Physics class may need.
Definition: physics.C:107
Base class for reading and handling initial conditions for physics classes.
static void set_is_axisymmetric(bool is_axisymmetric)
Set whether we should treat the problem as axisymmetric.
Definition: physics.h:129
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
virtual void compute_nonlocal_constraint_cache(AssemblyContext &)
Definition: physics.h:222
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Set which variables are time evolving.
Definition: physics.C:102
Interface with libMesh for solving Multiphysics problems.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: physics.C:124
virtual void compute_nonlocal_mass_residual_cache(AssemblyContext &)
Definition: physics.h:228
virtual void init_variables(libMesh::FEMSystem *)
Initialize variables for this physics.
Definition: physics.h:115
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
Definition: physics.C:65
virtual void compute_side_time_derivative_cache(AssemblyContext &)
Definition: physics.h:214
virtual void preassembly(MultiphysicsSystem &)
Perform any necessary setup before element assembly begins.
Definition: physics.h:162
virtual void damping_residual(bool, AssemblyContext &)
Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Definition: physics.h:198
virtual void nonlocal_time_derivative(bool, AssemblyContext &)
Time dependent part(s) of physics for scalar variables.
Definition: physics.h:182
static bool _is_axisymmetric
Caches whether we are solving an axisymmetric problem or not.
Definition: physics.h:269
void init_ics(libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Definition: physics.C:113
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
void set_is_steady(bool is_steady)
Sets whether this physics is to be solved with a steady solver or not.
Definition: physics.C:91

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