GRINS-0.6.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-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 #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"
37 #include "grins/cached_values.h"
38 #include "grins/parameter_user.h"
39 
40 //libMesh
41 #include "libmesh/libmesh.h"
42 #include "libmesh/point.h"
43 
44 // GRVY
45 #ifdef GRINS_HAVE_GRVY
46 #include "grvy.h" // GRVY timers
47 #endif
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 ParameterMultiPointer;
61 }
62 
64 namespace GRINS
65 {
66  // GRINS forward declarations
67  class BCHandlingBase;
68  class ICHandlingBase;
69  class NBCContainer;
70  class DBCContainer;
71  class AssemblyContext;
72  class MultiphysicsSystem;
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 read_input_options( const GetPot& input );
116 
118  virtual void init_variables( libMesh::FEMSystem* system ) = 0;
119 
121  virtual bool enabled_on_elem( const libMesh::Elem* elem );
122 
124 
126  void set_is_steady( bool is_steady );
127 
129  bool is_steady() const;
130 
132 
138  virtual void set_time_evolving_vars( libMesh::FEMSystem* system );
139 
141 
143  virtual void auxiliary_init( MultiphysicsSystem& system );
144 
146 
151  virtual void register_postprocessing_vars( const GetPot& input,
152  PostProcessedQuantities<libMesh::Real>& postprocessing );
153 
155  virtual void init_context( AssemblyContext& context );
156 
157  // residual and jacobian calculations
158  // element_*, side_* as *time_derivative, *constraint, *mass_residual
159 
161  virtual void element_time_derivative( bool compute_jacobian,
162  AssemblyContext& context,
163  CachedValues& cache );
164 
166  virtual void side_time_derivative( bool compute_jacobian,
167  AssemblyContext& context,
168  CachedValues& cache );
169 
171  virtual void nonlocal_time_derivative( bool compute_jacobian,
172  AssemblyContext& context,
173  CachedValues& cache );
174 
176  virtual void element_constraint( bool compute_jacobian,
177  AssemblyContext& context,
178  CachedValues& cache );
179 
181  virtual void side_constraint( bool compute_jacobian,
182  AssemblyContext& context,
183  CachedValues& cache );
184 
186  virtual void nonlocal_constraint( bool compute_jacobian,
187  AssemblyContext& context,
188  CachedValues& cache );
189 
191  virtual void mass_residual( bool compute_jacobian,
192  AssemblyContext& context,
193  CachedValues& cache );
194 
196  virtual void nonlocal_mass_residual( bool compute_jacobian,
197  AssemblyContext& context,
198  CachedValues& cache );
199 
200  void init_bcs( libMesh::FEMSystem* system );
201 
202  void init_ics( libMesh::FEMSystem* system,
204 
205  void attach_neumann_bound_func( GRINS::NBCContainer& neumann_bcs );
206 
207  void attach_dirichlet_bound_func( const GRINS::DBCContainer& dirichlet_bc );
208 
209  virtual void compute_element_time_derivative_cache( const AssemblyContext& context,
210  CachedValues& cache );
211 
212  virtual void compute_side_time_derivative_cache( const AssemblyContext& context,
213  CachedValues& cache );
214 
215  virtual void compute_nonlocal_time_derivative_cache( const AssemblyContext& context,
216  CachedValues& cache );
217 
218  virtual void compute_element_constraint_cache( const AssemblyContext& context,
219  CachedValues& cache );
220 
221  virtual void compute_side_constraint_cache( const AssemblyContext& context,
222  CachedValues& cache );
223 
224  virtual void compute_nonlocal_constraint_cache( const AssemblyContext& context,
225  CachedValues& cache );
226 
227  virtual void compute_mass_residual_cache( const AssemblyContext& context,
228  CachedValues& cache );
229 
230  virtual void compute_nonlocal_mass_residual_cache( const AssemblyContext& context,
231  CachedValues& cache );
232 
233  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
234  const AssemblyContext& context,
235  const libMesh::Point& point,
236  libMesh::Real& value );
237 
239 
241 
242 #ifdef GRINS_USE_GRVY_TIMERS
243  void attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer );
244 #endif
245 
246  protected:
247 
249 
255 
257 
259 
261  std::set<libMesh::subdomain_id_type> _enabled_subdomains;
262 
264 
266  static bool _is_steady;
267 
269 
270 #ifdef GRINS_USE_GRVY_TIMERS
271  GRVY::GRVY_Timer_Class* _timer;
272 #endif
273 
274  private:
275  Physics();
276 
277  }; // End Physics class declarations
278 
279  /* ------------------------- Inline Functions -------------------------*/
280 
281  inline
283  {
284  return _bc_handler;
285  }
286 
287  inline
289  {
290  return _ic_handler;
291  }
292 
293 } // End namespace GRINS
294 
295 #endif //GRINS_PHYSICS_H
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
Definition: physics.C:271
virtual void init_variables(libMesh::FEMSystem *system)=0
Initialize variables for this physics.
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:103
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
Simple helper class to setup general Dirichlet boundary conditions.
Definition: dbc_container.h:49
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register name of postprocessed quantity with PostProcessedQuantities.
Definition: physics.C:161
ICHandlingBase * get_ic_handler()
Definition: physics.h:288
virtual void nonlocal_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for scalar variables.
Definition: physics.C:229
virtual void compute_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:203
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
Definition: physics.C:150
virtual ~Physics()
Definition: physics.C:60
virtual void compute_element_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:167
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for element interiors.
Definition: physics.C:236
virtual void compute_side_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:191
GRINS namespace.
void init_bcs(libMesh::FEMSystem *system)
Definition: physics.C:118
Simple helper class to setup general Neumann boundary conditions.
Definition: nbc_container.h:41
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file. By default, nothing is read.
Definition: physics.C:70
virtual void compute_element_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:185
virtual void compute_side_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:173
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
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:83
std::string PhysicsName
virtual void side_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for boundaries of elements on the domain boundary.
Definition: physics.C:243
virtual void auxiliary_init(MultiphysicsSystem &system)
Any auxillary initialization a Physics class may need.
Definition: physics.C:113
Base class for reading and handling initial conditions for physics classes.
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for scalar variables.
Definition: physics.C:264
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Set which variables are time evolving.
Definition: physics.C:108
Interface with libMesh for solving Multiphysics problems.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: physics.C:156
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
Definition: physics.C:215
virtual void side_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for boundaries of elements on the domain boundary.
Definition: physics.C:222
virtual void compute_nonlocal_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:209
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
Definition: physics.C:144
virtual void compute_nonlocal_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:179
Base class for reading and handling boundary conditions for physics classes.
BCHandlingBase * get_bc_handler()
Definition: physics.h:282
virtual void compute_nonlocal_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:197
bool _is_axisymmetric
Definition: physics.h:268
virtual void nonlocal_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for scalar variables.
Definition: physics.C:250
void init_ics(libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Definition: physics.C:133
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:97
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Definition: physics.C:257

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