GRINS-0.6.0
Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
GRINS::MultiphysicsSystem Class Reference

Interface with libMesh for solving Multiphysics problems. More...

#include <multiphysics_sys.h>

Inheritance diagram for GRINS::MultiphysicsSystem:
Inheritance graph
[legend]
Collaboration diagram for GRINS::MultiphysicsSystem:
Collaboration graph
[legend]

Public Member Functions

 MultiphysicsSystem (libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
 Constructor. Will be called by libMesh only. More...
 
 ~MultiphysicsSystem ()
 Destructor. Clean up all physics allocations. More...
 
void attach_physics_list (PhysicsList physics_list)
 PhysicsList gets built by GRINS::PhysicsFactory and attached here. More...
 
virtual void read_input_options (const GetPot &input)
 Reads input options for this class and all physics that are enabled. More...
 
virtual void init_data ()
 System initialization. Calls each physics implementation of init_variables() More...
 
void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Each Physics will register their postprocessed quantities with this call. More...
 
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. More...
 
virtual libMesh::AutoPtr< libMesh::DiffContext > build_context ()
 Override FEMSystem::build_context in order to use our own AssemblyContext. More...
 
virtual void init_context (libMesh::DiffContext &context)
 Context initialization. Calls each physics implementation of init_context() More...
 
virtual bool element_time_derivative (bool request_jacobian, libMesh::DiffContext &context)
 Element interior contributions to $F(u)$ which have time varying components. More...
 
virtual bool side_time_derivative (bool request_jacobian, libMesh::DiffContext &context)
 Boundary contributions to $F(u)$ which have time varying components. More...
 
virtual bool nonlocal_time_derivative (bool request_jacobian, libMesh::DiffContext &context)
 Contributions to $F(u)$ on SCALAR variables which have time varying components. More...
 
virtual bool element_constraint (bool request_jacobian, libMesh::DiffContext &context)
 
virtual bool side_constraint (bool request_jacobian, libMesh::DiffContext &context)
 Boundary contributions to $F(u)$ which do not have time varying components. More...
 
virtual bool nonlocal_constraint (bool request_jacobian, libMesh::DiffContext &context)
 Contributions to $F(u)$ on SCALAR variables which do not have time varying components. More...
 
virtual bool mass_residual (bool request_jacobian, libMesh::DiffContext &context)
 Contributions to $M(u)\dot{u}$. More...
 
virtual bool nonlocal_mass_residual (bool request_jacobian, libMesh::DiffContext &context)
 Contributions to $M(u)\dot{u}$ on SCALAR variables. More...
 
bool has_physics (const std::string physics_name) const
 Query to check if a particular physics has been enabled. More...
 
std::tr1::shared_ptr< GRINS::Physicsget_physics (const std::string physics_name)
 
std::tr1::shared_ptr< GRINS::Physicsget_physics (const std::string physics_name) const
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 

Private Types

typedef void(GRINS::Physics::* ResFuncType) (bool, AssemblyContext &, CachedValues &)
 
typedef void(GRINS::Physics::* CacheFuncType) (const AssemblyContext &, CachedValues &)
 

Private Member Functions

bool _general_residual (bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
 

Private Attributes

PhysicsList _physics_list
 Container of pointers to GRINS::Physics classes requested at runtime. More...
 
bool _use_numerical_jacobians_only
 

Detailed Description

Interface with libMesh for solving Multiphysics problems.

MultiphysicsSystem (through libMesh::FEMSystem) solves the following equation:

$M(u)\dot{u} = F(u)$

M = mass matrix u = solution vector F = time derivative

Note that for the nonlinear system that is solved for implicit time stepping is:

$M(u_{\theta})(u^n - u^{n+1}) + \Delta t F(u) = 0$

_time_derivative correspond to calculating terms for $F(u)$ _mass_residual correspond to calculating terms for $M(u)\dot{u}$

Definition at line 81 of file multiphysics_sys.h.

Member Typedef Documentation

typedef void(GRINS::Physics::* GRINS::MultiphysicsSystem::CacheFuncType) (const AssemblyContext &, CachedValues &)
private

Definition at line 191 of file multiphysics_sys.h.

typedef void(GRINS::Physics::* GRINS::MultiphysicsSystem::ResFuncType) (bool, AssemblyContext &, CachedValues &)
private

Definition at line 190 of file multiphysics_sys.h.

Constructor & Destructor Documentation

GRINS::MultiphysicsSystem::MultiphysicsSystem ( libMesh::EquationSystems &  es,
const std::string &  name,
const unsigned int  number 
)

Constructor. Will be called by libMesh only.

Definition at line 40 of file multiphysics_sys.C.

43  : FEMSystem(es, name, number),
45  {
46  return;
47  }
GRINS::MultiphysicsSystem::~MultiphysicsSystem ( )

Destructor. Clean up all physics allocations.

Definition at line 49 of file multiphysics_sys.C.

50  {
51  return;
52  }

Member Function Documentation

bool GRINS::MultiphysicsSystem::_general_residual ( bool  request_jacobian,
libMesh::DiffContext &  context,
ResFuncType  resfunc,
CacheFuncType  cachefunc 
)
private

Definition at line 228 of file multiphysics_sys.C.

References _physics_list, and _use_numerical_jacobians_only.

Referenced by element_constraint(), element_time_derivative(), mass_residual(), nonlocal_constraint(), nonlocal_mass_residual(), nonlocal_time_derivative(), side_constraint(), and side_time_derivative().

232  {
233  AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
234 
235  bool compute_jacobian = true;
236  if( !request_jacobian || _use_numerical_jacobians_only ) compute_jacobian = false;
237 
238  CachedValues cache;
239 
240  // Now compute cache for this element
241  for( PhysicsListIter physics_iter = _physics_list.begin();
242  physics_iter != _physics_list.end();
243  physics_iter++ )
244  {
245  // boost::shared_ptr gets confused by operator->*
246  ((*(physics_iter->second)).*cachefunc)( c, cache );
247  }
248 
249  // Loop over each physics and compute their contributions
250  for( PhysicsListIter physics_iter = _physics_list.begin();
251  physics_iter != _physics_list.end();
252  physics_iter++ )
253  {
254  if(c.has_elem())
255  {
256  if( (physics_iter->second)->enabled_on_elem( &c.get_elem() ) )
257  {
258  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
259  }
260  }
261  else
262  {
263  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
264  }
265  }
266 
267  // TODO: Need to think about the implications of this because there might be some
268  // TODO: jacobian terms we don't want to compute for efficiency reasons
269  return compute_jacobian;
270  }
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
void GRINS::MultiphysicsSystem::attach_physics_list ( PhysicsList  physics_list)

PhysicsList gets built by GRINS::PhysicsFactory and attached here.

Definition at line 54 of file multiphysics_sys.C.

References _physics_list.

Referenced by GRINS::Simulation::init_multiphysics_system().

55  {
56  _physics_list = physics_list;
57  return;
58  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
libMesh::AutoPtr< libMesh::DiffContext > GRINS::MultiphysicsSystem::build_context ( )
virtual

Override FEMSystem::build_context in order to use our own AssemblyContext.

Definition at line 159 of file multiphysics_sys.C.

160  {
161  AssemblyContext* context = new AssemblyContext(*this);
162 
163  libMesh::AutoPtr<libMesh::DiffContext> ap(context);
164 
165  libMesh::DifferentiablePhysics* phys = libMesh::FEMSystem::get_physics();
166 
167  libmesh_assert(phys);
168 
169  // If we are solving a moving mesh problem, tell that to the Context
170  context->set_mesh_system(phys->get_mesh_system());
171  context->set_mesh_x_var(phys->get_mesh_x_var());
172  context->set_mesh_y_var(phys->get_mesh_y_var());
173  context->set_mesh_z_var(phys->get_mesh_z_var());
174 
175  ap->set_deltat_pointer( &deltat );
176 
177  // If we are solving the adjoint problem, tell that to the Context
178  ap->is_adjoint() = this->get_time_solver().is_adjoint();
179 
180  return ap;
181  }
void GRINS::MultiphysicsSystem::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Definition at line 373 of file multiphysics_sys.C.

References _physics_list.

377  {
378  for( PhysicsListIter physics_iter = _physics_list.begin();
379  physics_iter != _physics_list.end();
380  physics_iter++ )
381  {
382  // Only compute if physics is active on current subdomain or globally
383  if( (physics_iter->second)->enabled_on_elem( &context.get_elem() ) )
384  {
385  (physics_iter->second)->compute_postprocessed_quantity( quantity_index, context, point, value );
386  }
387  }
388  return;
389  }
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
bool GRINS::MultiphysicsSystem::element_constraint ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Element interior contributions to $F(u)$ which do not have time varying components. Element interior contributions to $F(u)$ which do not have time varying components.

Definition at line 302 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_element_constraint_cache(), and GRINS::Physics::element_constraint().

304  {
305  return this->_general_residual
306  (request_jacobian,
307  context,
310  }
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_element_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:185
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
bool GRINS::MultiphysicsSystem::element_time_derivative ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Element interior contributions to $F(u)$ which have time varying components.

Definition at line 272 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_element_time_derivative_cache(), and GRINS::Physics::element_time_derivative().

274  {
275  return this->_general_residual
276  (request_jacobian,
277  context,
280  }
virtual void compute_element_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:167
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
std::tr1::shared_ptr< Physics > GRINS::MultiphysicsSystem::get_physics ( const std::string  physics_name)

Definition at line 352 of file multiphysics_sys.C.

References _physics_list.

Referenced by GRINS::Simulation::attach_dirichlet_bc_funcs(), and GRINS::Simulation::attach_neumann_bc_funcs().

353  {
354  if( _physics_list.find( physics_name ) == _physics_list.end() )
355  {
356  std::cerr << "Error: Could not find physics " << physics_name << std::endl;
357  libmesh_error();
358  }
359 
360  return _physics_list[physics_name];
361  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
std::tr1::shared_ptr< GRINS::Physics > GRINS::MultiphysicsSystem::get_physics ( const std::string  physics_name) const
inline

Definition at line 201 of file multiphysics_sys.h.

References _physics_list.

202  {
203  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
204 
205  return _physics_list.find(physics_name)->second;
206  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
bool GRINS::MultiphysicsSystem::has_physics ( const std::string  physics_name) const

Query to check if a particular physics has been enabled.

Definition at line 363 of file multiphysics_sys.C.

References _physics_list.

364  {
365  bool has_physics = false;
366 
367  if( _physics_list.find(physics_name) != _physics_list.end() )
368  has_physics = true;
369 
370  return has_physics;
371  }
bool has_physics(const std::string physics_name) const
Query to check if a particular physics has been enabled.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
void GRINS::MultiphysicsSystem::init_context ( libMesh::DiffContext &  context)
virtual

Context initialization. Calls each physics implementation of init_context()

Definition at line 212 of file multiphysics_sys.C.

References _physics_list.

213  {
214  AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
215 
216  //Loop over each physics to initialize relevant variable structures for assembling system
217  for( PhysicsListIter physics_iter = _physics_list.begin();
218  physics_iter != _physics_list.end();
219  physics_iter++ )
220  {
221  (physics_iter->second)->init_context( c );
222  }
223 
224  return;
225  }
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
virtual void init_context(libMesh::DiffContext &context)
Context initialization. Calls each physics implementation of init_context()
void GRINS::MultiphysicsSystem::init_data ( )
virtual

System initialization. Calls each physics implementation of init_variables()

Todo:
Figure out how to tell compilers not to fuse this loop when they want to be aggressive.

Definition at line 88 of file multiphysics_sys.C.

References _physics_list.

89  {
90  // Need this to be true because of our overloading of the
91  // mass_residual function.
92  // This is data in FEMSystem. MUST be set before FEMSystem::init_data.
93  use_fixed_solution = true;
94 
95  // Initalize all the variables. We pass this pointer for the system.
96  /* NOTE: We CANNOT fuse this loop with the others. This loop
97  MUST complete first. */
100  for( PhysicsListIter physics_iter = _physics_list.begin();
101  physics_iter != _physics_list.end();
102  physics_iter++ )
103  {
104  (physics_iter->second)->init_variables( this );
105  }
106 
107  // Now set time_evolving variables
108  for( PhysicsListIter physics_iter = _physics_list.begin();
109  physics_iter != _physics_list.end();
110  physics_iter++ )
111  {
112  (physics_iter->second)->set_time_evolving_vars( this );
113  }
114 
115  // Set whether the problem we're solving is steady or not
116  // Since the variable is static, just call one Physics class
117  {
118  (_physics_list.begin()->second)->set_is_steady((this->time_solver)->is_steady());
119  }
120 
121  for( PhysicsListIter physics_iter = _physics_list.begin();
122  physics_iter != _physics_list.end();
123  physics_iter++ )
124  {
125  // Initialize builtin BC's for each physics
126  (physics_iter->second)->init_bcs( this );
127  }
128 
129  // Next, call parent init_data function to intialize everything.
130  libMesh::FEMSystem::init_data();
131 
132  // After solution has been initialized we can project initial
133  // conditions to it
135  for( PhysicsListIter physics_iter = _physics_list.begin();
136  physics_iter != _physics_list.end();
137  physics_iter++ )
138  {
139  // Initialize builtin IC's for each physics
140  (physics_iter->second)->init_ics( this, ic_function );
141  }
142 
143  if (ic_function.n_subfunctions())
144  {
145  this->project_solution(&ic_function);
146  }
147 
148  // Now do any auxillary initialization required by each Physics
149  for( PhysicsListIter physics_iter = _physics_list.begin();
150  physics_iter != _physics_list.end();
151  physics_iter++ )
152  {
153  (physics_iter->second)->auxiliary_init( *this );
154  }
155 
156  return;
157  }
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
bool GRINS::MultiphysicsSystem::mass_residual ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Contributions to $M(u)\dot{u}$.

Definition at line 332 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_mass_residual_cache(), and GRINS::Physics::mass_residual().

334  {
335  return this->_general_residual
336  (request_jacobian,
337  context,
340  }
virtual void compute_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:203
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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
bool GRINS::MultiphysicsSystem::nonlocal_constraint ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Contributions to $F(u)$ on SCALAR variables which do not have time varying components.

Definition at line 322 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_nonlocal_constraint_cache(), and GRINS::Physics::nonlocal_constraint().

324  {
325  return this->_general_residual
326  (request_jacobian,
327  context,
330  }
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual void compute_nonlocal_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:197
virtual void nonlocal_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for scalar variables.
Definition: physics.C:250
bool GRINS::MultiphysicsSystem::nonlocal_mass_residual ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Contributions to $M(u)\dot{u}$ on SCALAR variables.

Definition at line 342 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_nonlocal_mass_residual_cache(), and GRINS::Physics::nonlocal_mass_residual().

344  {
345  return this->_general_residual
346  (request_jacobian,
347  context,
350  }
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for scalar variables.
Definition: physics.C:264
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual void compute_nonlocal_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:209
bool GRINS::MultiphysicsSystem::nonlocal_time_derivative ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Contributions to $F(u)$ on SCALAR variables which have time varying components.

Definition at line 292 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_nonlocal_time_derivative_cache(), and GRINS::Physics::nonlocal_time_derivative().

294  {
295  return this->_general_residual
296  (request_jacobian,
297  context,
300  }
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual void compute_nonlocal_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:179
void GRINS::MultiphysicsSystem::read_input_options ( const GetPot &  input)
virtual

Reads input options for this class and all physics that are enabled.

This function reads the input options for the MultiphysicsSystem class and then enables each of the requested physics in the system. Finally, the input options for each of the physics will be read.

Todo:
Remove old print_residual nomenclature

Definition at line 60 of file multiphysics_sys.C.

References _use_numerical_jacobians_only.

Referenced by GRINS::Simulation::init_multiphysics_system().

61  {
62  // Read options for MultiphysicsSystem first
63  this->verify_analytic_jacobians = input("linear-nonlinear-solver/verify_analytic_jacobians", 0.0 );
64  this->print_solution_norms = input("screen-options/print_solution_norms", false );
65  this->print_solutions = input("screen-options/print_solutions", false );
66  this->print_residual_norms = input("screen-options/print_residual_norms", false );
67 
68  // backwards compatibility with old config files.
70  this->print_residuals = input("screen-options/print_residual", false );
71  if (this->print_residuals)
72  libmesh_deprecated();
73 
74  this->print_residuals = input("screen-options/print_residuals", this->print_residuals );
75  this->print_jacobian_norms = input("screen-options/print_jacobian_norms", false );
76  this->print_jacobians = input("screen-options/print_jacobians", false );
77  this->print_element_solutions = input("screen-options/print_element_solutions", false );
78  this->print_element_residuals = input("screen-options/print_element_residuals", false );
79  this->print_element_jacobians = input("screen-options/print_element_jacobians", false );
80 
81  _use_numerical_jacobians_only = input("linear-nonlinear-solver/use_numerical_jacobians_only", false );
82 
83  numerical_jacobian_h =
84  input("linear-nonlinear-solver/numerical_jacobian_h",
85  numerical_jacobian_h);
86  }
void GRINS::MultiphysicsSystem::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
)

Each Physics will register its copy(s) of an independent variable.

Definition at line 197 of file multiphysics_sys.C.

Referenced by GRINS::ParameterManager::initialize().

199  {
200  //Loop over each physics to ask each for the requested parameter
201  for( PhysicsListIter physics_iter = _physics_list.begin();
202  physics_iter != _physics_list.end();
203  physics_iter++ )
204  {
205  (physics_iter->second)->register_parameter( param_name,
206  param_pointer );
207  }
208  }
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
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.
void GRINS::MultiphysicsSystem::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)

Each Physics will register their postprocessed quantities with this call.

Definition at line 183 of file multiphysics_sys.C.

References _physics_list.

Referenced by GRINS::Simulation::init_multiphysics_system().

185  {
186  for( PhysicsListIter physics_iter = _physics_list.begin();
187  physics_iter != _physics_list.end();
188  physics_iter++ )
189  {
190  (physics_iter->second)->register_postprocessing_vars( input, postprocessing );
191  }
192 
193  return;
194  }
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:60
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
bool GRINS::MultiphysicsSystem::side_constraint ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Boundary contributions to $F(u)$ which do not have time varying components.

Definition at line 312 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_side_constraint_cache(), and GRINS::Physics::side_constraint().

314  {
315  return this->_general_residual
316  (request_jacobian,
317  context,
320  }
virtual void compute_side_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:191
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
bool GRINS::MultiphysicsSystem::side_time_derivative ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

Boundary contributions to $F(u)$ which have time varying components.

Definition at line 282 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_side_time_derivative_cache(), and GRINS::Physics::side_time_derivative().

284  {
285  return this->_general_residual
286  (request_jacobian,
287  context,
290  }
virtual void compute_side_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:173
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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

Member Data Documentation

PhysicsList GRINS::MultiphysicsSystem::_physics_list
private

Container of pointers to GRINS::Physics classes requested at runtime.

Set using the attach_physics_list method as construction is taken care of by GRINS::PhysicsFactory.

Definition at line 181 of file multiphysics_sys.h.

Referenced by _general_residual(), attach_physics_list(), compute_postprocessed_quantity(), get_physics(), has_physics(), init_context(), init_data(), and register_postprocessing_vars().

bool GRINS::MultiphysicsSystem::_use_numerical_jacobians_only
private

Definition at line 183 of file multiphysics_sys.h.

Referenced by _general_residual(), and read_input_options().


The documentation for this class was generated from the following files:

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