GRINS-0.7.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-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 #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 
40 //libMesh
41 #include "libmesh/libmesh.h"
42 #include "libmesh/point.h"
43 #include "libmesh/fe_base.h"
44 
45 // GRVY
46 #ifdef GRINS_HAVE_GRVY
47 #include "grvy.h" // GRVY timers
48 #endif
49 
50 // libMesh forward declarations
51 class GetPot;
52 namespace libMesh
53 {
54  template <typename Scalar>
55  class CompositeFunction;
56 
57  class FEMSystem;
58  class Elem;
59 
60  template <typename Scalar>
61  class ParameterMultiAccessor;
62 }
63 
65 namespace GRINS
66 {
67  // GRINS forward declarations
68  class BCHandlingBase;
69  class ICHandlingBase;
70  class NBCContainer;
71  class DBCContainer;
72  class AssemblyContext;
73  class MultiphysicsSystem;
74 
75  template <typename Scalar>
76  class PostProcessedQuantities;
77 
79 
99  //TODO: is it F(u) or F(u_{\theta})?
100 
101  // element* routines work on element interiors
102  // side* routines work on element sides
103 
104  // *_time_derivative correspond to calculating terms for F(u)
105  // *_mass_residual correspond to calculating terms for M(u)\dot{u}
106 
107  class Physics : public ParameterUser
108  {
109 
110  public:
111 
112  Physics( const GRINS::PhysicsName& physics_name, const GetPot& input );
113  virtual ~Physics();
114 
116  virtual void init_variables( libMesh::FEMSystem* system ) = 0;
117 
119  virtual bool enabled_on_elem( const libMesh::Elem* elem );
120 
122 
124  void set_is_steady( bool is_steady );
125 
127  bool is_steady() const;
128 
132 
133  static bool is_axisymmetric()
134  { return _is_axisymmetric; }
135 
137 
143  virtual void set_time_evolving_vars( libMesh::FEMSystem* system );
144 
146 
148  virtual void auxiliary_init( MultiphysicsSystem& system );
149 
151 
156  virtual void register_postprocessing_vars( const GetPot& input,
157  PostProcessedQuantities<libMesh::Real>& postprocessing );
158 
160  virtual void init_context( AssemblyContext& context );
161 
162  // residual and jacobian calculations
163  // element_*, side_* as *time_derivative, *constraint, *mass_residual
164 
166  virtual void element_time_derivative( bool compute_jacobian,
167  AssemblyContext& context,
168  CachedValues& cache );
169 
171  virtual void side_time_derivative( bool compute_jacobian,
172  AssemblyContext& context,
173  CachedValues& cache );
174 
176  virtual void nonlocal_time_derivative( bool compute_jacobian,
177  AssemblyContext& context,
178  CachedValues& cache );
179 
181  virtual void element_constraint( bool compute_jacobian,
182  AssemblyContext& context,
183  CachedValues& cache );
184 
186  virtual void side_constraint( bool compute_jacobian,
187  AssemblyContext& context,
188  CachedValues& cache );
189 
191  virtual void nonlocal_constraint( bool compute_jacobian,
192  AssemblyContext& context,
193  CachedValues& cache );
194 
196  virtual void damping_residual( bool compute_jacobian,
197  AssemblyContext& context,
198  CachedValues& cache );
199 
201  virtual void mass_residual( bool compute_jacobian,
202  AssemblyContext& context,
203  CachedValues& cache );
204 
206  virtual void nonlocal_mass_residual( bool compute_jacobian,
207  AssemblyContext& context,
208  CachedValues& cache );
209 
210  void init_ics( libMesh::FEMSystem* system,
212 
213  virtual void compute_element_time_derivative_cache( const AssemblyContext& context,
214  CachedValues& cache );
215 
216  virtual void compute_side_time_derivative_cache( const AssemblyContext& context,
217  CachedValues& cache );
218 
219  virtual void compute_nonlocal_time_derivative_cache( const AssemblyContext& context,
220  CachedValues& cache );
221 
222  virtual void compute_element_constraint_cache( const AssemblyContext& context,
223  CachedValues& cache );
224 
225  virtual void compute_side_constraint_cache( const AssemblyContext& context,
226  CachedValues& cache );
227 
228  virtual void compute_nonlocal_constraint_cache( const AssemblyContext& context,
229  CachedValues& cache );
230 
231  virtual void compute_damping_residual_cache( const AssemblyContext& context,
232  CachedValues& cache );
233 
234  virtual void compute_mass_residual_cache( const AssemblyContext& context,
235  CachedValues& cache );
236 
237  virtual void compute_nonlocal_mass_residual_cache( const AssemblyContext& context,
238  CachedValues& cache );
239 
240  virtual void compute_postprocessed_quantity( unsigned int quantity_index,
241  const AssemblyContext& context,
242  const libMesh::Point& point,
243  libMesh::Real& value );
244 
246 
247 #ifdef GRINS_USE_GRVY_TIMERS
248  void attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer );
249 #endif
250 
251  protected:
252 
254  libMesh::UniquePtr<libMesh::FEGenericBase<libMesh::Real> > build_new_fe( const libMesh::Elem* elem,
255  const libMesh::FEGenericBase<libMesh::Real>* fe,
256  const libMesh::Point p );
257 
258  void parse_enabled_subdomains( const GetPot& input,
259  const std::string& physics_name );
260 
262 
268 
270 
272  std::set<libMesh::subdomain_id_type> _enabled_subdomains;
273 
275 
277  static bool _is_steady;
278 
280  static bool _is_axisymmetric;
281 
282 #ifdef GRINS_USE_GRVY_TIMERS
283  GRVY::GRVY_Timer_Class* _timer;
284 #endif
285 
286  private:
287  Physics();
288 
289  }; // End Physics class declarations
290 
291  /* ------------------------- Inline Functions -------------------------*/
292  inline
294  {
295  return _ic_handler;
296  }
297 
298 } // End namespace GRINS
299 
300 #endif //GRINS_PHYSICS_H
static bool is_axisymmetric()
Definition: physics.h:133
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
Definition: physics.C:251
virtual void init_variables(libMesh::FEMSystem *system)=0
Initialize variables for this physics.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:96
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:107
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register name of postprocessed quantity with PostProcessedQuantities.
Definition: physics.C:128
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:259
virtual void damping_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Definition: physics.C:230
ICHandlingBase * get_ic_handler()
Definition: physics.h:293
virtual void nonlocal_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for scalar variables.
Definition: physics.C:202
virtual void compute_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:176
virtual ~Physics()
Definition: physics.C:59
virtual void compute_element_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:134
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for element interiors.
Definition: physics.C:209
virtual void compute_side_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:158
GRINS namespace.
virtual void compute_element_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:152
virtual void compute_side_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:140
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:267
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:76
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:216
virtual void auxiliary_init(MultiphysicsSystem &system)
Any auxillary initialization a Physics class may need.
Definition: physics.C:106
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:244
static void set_is_axisymmetric(bool is_axisymmetric)
Set whether we should treat the problem as axisymmetric.
Definition: physics.h:130
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:277
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Set which variables are time evolving.
Definition: physics.C:101
Interface with libMesh for solving Multiphysics problems.
virtual void compute_damping_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:170
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: physics.C:123
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
Definition: physics.C:188
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
Definition: physics.C:64
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:195
virtual void compute_nonlocal_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:182
virtual void compute_nonlocal_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:146
virtual void compute_nonlocal_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:164
static bool _is_axisymmetric
Caches whether we are solving an axisymmetric problem or not.
Definition: physics.h:280
virtual void nonlocal_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for scalar variables.
Definition: physics.C:223
void init_ics(libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Definition: physics.C:112
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:272
void set_is_steady(bool is_steady)
Sets whether this physics is to be solved with a steady solver or not.
Definition: physics.C:90
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:237

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