GRINS-0.7.0
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
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...
 
const PhysicsListget_physics_list () const
 
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::ParameterMultiAccessor< libMesh::Number > &param_pointer)
 Each Physics will register its copy(s) of an independent variable. More...
 
virtual libMesh::UniquePtr< 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 damping_residual (bool request_jacobian, libMesh::DiffContext &context)
 Contributions to $C(u)\dot{u}$. 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...
 
SharedPtr< GRINS::Physicsget_physics (const std::string physics_name)
 
SharedPtr< 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)
 
std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs ()
 
const std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs () const
 

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)
 
void get_active_neumann_bcs (BoundaryID bc_id, const std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs, std::vector< SharedPtr< NeumannBCContainer > > &active_neumann_bcs)
 Extract the bcs from neumann_bcs that are active on bc_id and return them in active_neumann_bcs. More...
 
bool apply_neumann_bcs (bool request_jacobian, libMesh::DiffContext &context)
 Applies the subset of _neumann_bcs that are active on the current element side. More...
 

Private Attributes

PhysicsList _physics_list
 Container of pointers to GRINS::Physics classes requested at runtime. More...
 
bool _use_numerical_jacobians_only
 
std::vector< std::string > _numerical_jacobian_h_variables
 
std::vector< libMesh::Real > _numerical_jacobian_h_values
 
const GetPot * _input
 Cached for helping build boundary conditions. More...
 
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
 Neumann boundary conditions. More...
 

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 82 of file multiphysics_sys.h.

Member Typedef Documentation

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

Definition at line 225 of file multiphysics_sys.h.

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

Definition at line 224 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 43 of file multiphysics_sys.C.

46  : FEMSystem(es, name, number),
48  {}
GRINS::MultiphysicsSystem::~MultiphysicsSystem ( )
inline

Destructor. Clean up all physics allocations.

Definition at line 92 of file multiphysics_sys.h.

92 {};

Member Function Documentation

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

Definition at line 260 of file multiphysics_sys.C.

References _physics_list, and _use_numerical_jacobians_only.

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

264  {
265  AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
266 
267  bool compute_jacobian = true;
268  if( !request_jacobian || _use_numerical_jacobians_only ) compute_jacobian = false;
269 
270  CachedValues cache;
271 
272  // Now compute cache for this element
273  for( PhysicsListIter physics_iter = _physics_list.begin();
274  physics_iter != _physics_list.end();
275  physics_iter++ )
276  {
277  // shared_ptr gets confused by operator->*
278  ((*(physics_iter->second)).*cachefunc)( c, cache );
279  }
280 
281  // Loop over each physics and compute their contributions
282  for( PhysicsListIter physics_iter = _physics_list.begin();
283  physics_iter != _physics_list.end();
284  physics_iter++ )
285  {
286  if(c.has_elem())
287  {
288  if( (physics_iter->second)->enabled_on_elem( &c.get_elem() ) )
289  {
290  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
291  }
292  }
293  else
294  {
295  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
296  }
297  }
298 
299  // TODO: Need to think about the implications of this because there might be some
300  // TODO: jacobian terms we don't want to compute for efficiency reasons
301  return compute_jacobian;
302  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
bool GRINS::MultiphysicsSystem::apply_neumann_bcs ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
private

Applies the subset of _neumann_bcs that are active on the current element side.

Definition at line 450 of file multiphysics_sys.C.

References _neumann_bcs, _use_numerical_jacobians_only, get_active_neumann_bcs(), GRINS::Physics::is_axisymmetric(), and GRINS::FEVariablesBase::neumann_bc_sign().

Referenced by side_time_derivative().

452  {
453  AssemblyContext& assembly_context =
454  libMesh::libmesh_cast_ref<AssemblyContext&>( context );
455 
456  std::vector<BoundaryID> ids = assembly_context.side_boundary_ids();
457 
458  bool compute_jacobian = request_jacobian;
459  if( !request_jacobian || _use_numerical_jacobians_only ) compute_jacobian = false;
460 
461  for( std::vector<BoundaryID>::const_iterator it = ids.begin();
462  it != ids.end(); it++ )
463  {
464  BoundaryID bc_id = *it;
465 
466  libmesh_assert_not_equal_to(bc_id, libMesh::BoundaryInfo::invalid_id);
467 
468  // Retreive the NeumannBCContainers that are active on the current bc_id
469  std::vector<SharedPtr<NeumannBCContainer> > active_neumann_bcs;
470  this->get_active_neumann_bcs( bc_id, _neumann_bcs, active_neumann_bcs );
471 
472  if( !active_neumann_bcs.empty() )
473  {
474  typedef std::vector<SharedPtr<NeumannBCContainer> >::iterator BCIt;
475 
476  for( BCIt container = active_neumann_bcs.begin(); container < active_neumann_bcs.end(); ++container )
477  {
478  SharedPtr<NeumannBCAbstract>& func = (*container)->get_func();
479 
480  const FEVariablesBase& var = (*container)->get_fe_var();
481 
482  func->eval_flux( compute_jacobian, assembly_context,
483  var.neumann_bc_sign(), Physics::is_axisymmetric() );
484  }
485  }
486  } // end loop over boundary ids
487 
488  return compute_jacobian;
489  }
static bool is_axisymmetric()
Definition: physics.h:133
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
libMesh::boundary_id_type BoundaryID
More descriptive name of the type used for boundary ids.
Definition: var_typedefs.h:56
void get_active_neumann_bcs(BoundaryID bc_id, const std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs, std::vector< SharedPtr< NeumannBCContainer > > &active_neumann_bcs)
Extract the bcs from neumann_bcs that are active on bc_id and return them in active_neumann_bcs.
void GRINS::MultiphysicsSystem::attach_physics_list ( PhysicsList  physics_list)

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

Definition at line 50 of file multiphysics_sys.C.

References _physics_list.

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

51  {
52  _physics_list = physics_list;
53  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
libMesh::UniquePtr< libMesh::DiffContext > GRINS::MultiphysicsSystem::build_context ( )
virtual

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

Definition at line 191 of file multiphysics_sys.C.

192  {
193  AssemblyContext* context = new AssemblyContext(*this);
194 
195  libMesh::UniquePtr<libMesh::DiffContext> ap(context);
196 
197  libMesh::DifferentiablePhysics* phys = libMesh::FEMSystem::get_physics();
198 
199  libmesh_assert(phys);
200 
201  // If we are solving a moving mesh problem, tell that to the Context
202  context->set_mesh_system(phys->get_mesh_system());
203  context->set_mesh_x_var(phys->get_mesh_x_var());
204  context->set_mesh_y_var(phys->get_mesh_y_var());
205  context->set_mesh_z_var(phys->get_mesh_z_var());
206 
207  ap->set_deltat_pointer( &deltat );
208 
209  // If we are solving the adjoint problem, tell that to the Context
210  ap->is_adjoint() = this->get_time_solver().is_adjoint();
211 
212  return ap;
213  }
void GRINS::MultiphysicsSystem::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Definition at line 421 of file multiphysics_sys.C.

References _physics_list.

425  {
426  for( PhysicsListIter physics_iter = _physics_list.begin();
427  physics_iter != _physics_list.end();
428  physics_iter++ )
429  {
430  // Only compute if physics is active on current subdomain or globally
431  if( (physics_iter->second)->enabled_on_elem( &context.get_elem() ) )
432  {
433  (physics_iter->second)->compute_postprocessed_quantity( quantity_index, context, point, value );
434  }
435  }
436  return;
437  }
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)
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
bool GRINS::MultiphysicsSystem::damping_residual ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

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

Definition at line 370 of file multiphysics_sys.C.

References _general_residual(), GRINS::Physics::compute_damping_residual_cache(), and GRINS::Physics::damping_residual().

372  {
373  return this->_general_residual
374  (request_jacobian,
375  context,
378  }
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
virtual void compute_damping_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:170
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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 340 of file multiphysics_sys.C.

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

342  {
343  return this->_general_residual
344  (request_jacobian,
345  context,
348  }
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_element_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:152
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 304 of file multiphysics_sys.C.

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

306  {
307  return this->_general_residual
308  (request_jacobian,
309  context,
312  }
virtual void compute_element_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:134
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
void GRINS::MultiphysicsSystem::get_active_neumann_bcs ( BoundaryID  bc_id,
const std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs,
std::vector< SharedPtr< NeumannBCContainer > > &  active_neumann_bcs 
)
private

Extract the bcs from neumann_bcs that are active on bc_id and return them in active_neumann_bcs.

Definition at line 439 of file multiphysics_sys.C.

Referenced by apply_neumann_bcs().

442  {
443  // Manually writing the loop since std::copy_if is C++11 only
444  for( std::vector<SharedPtr<NeumannBCContainer> >::const_iterator it = neumann_bcs.begin();
445  it < neumann_bcs.end(); ++it )
446  if( (*it)->has_bc_id( bc_id ) )
447  active_neumann_bcs.push_back( *it );
448  }
std::vector<SharedPtr<NeumannBCContainer> >& GRINS::MultiphysicsSystem::get_neumann_bcs ( )
inline

Definition at line 179 of file multiphysics_sys.h.

References _neumann_bcs.

180  { return _neumann_bcs; }
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
const std::vector<SharedPtr<NeumannBCContainer> >& GRINS::MultiphysicsSystem::get_neumann_bcs ( ) const
inline

Definition at line 182 of file multiphysics_sys.h.

References _neumann_bcs.

183  { return _neumann_bcs; }
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
SharedPtr< Physics > GRINS::MultiphysicsSystem::get_physics ( const std::string  physics_name)

Definition at line 400 of file multiphysics_sys.C.

References _physics_list.

401  {
402  if( _physics_list.find( physics_name ) == _physics_list.end() )
403  {
404  std::cerr << "Error: Could not find physics " << physics_name << std::endl;
405  libmesh_error();
406  }
407 
408  return _physics_list[physics_name];
409  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
SharedPtr< GRINS::Physics > GRINS::MultiphysicsSystem::get_physics ( const std::string  physics_name) const
inline

Definition at line 244 of file multiphysics_sys.h.

References _physics_list.

245  {
246  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
247 
248  return _physics_list.find(physics_name)->second;
249  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
const PhysicsList& GRINS::MultiphysicsSystem::get_physics_list ( ) const
inline

Definition at line 97 of file multiphysics_sys.h.

References _physics_list.

Referenced by GRINS::OldStyleBCBuilder::build_bcs().

98  { return _physics_list; }
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 411 of file multiphysics_sys.C.

References _physics_list.

412  {
413  bool has_physics = false;
414 
415  if( _physics_list.find(physics_name) != _physics_list.end() )
416  has_physics = true;
417 
418  return has_physics;
419  }
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 244 of file multiphysics_sys.C.

References _physics_list.

245  {
246  AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
247 
248  //Loop over each physics to initialize relevant variable structures for assembling system
249  for( PhysicsListIter physics_iter = _physics_list.begin();
250  physics_iter != _physics_list.end();
251  physics_iter++ )
252  {
253  (physics_iter->second)->init_context( c );
254  }
255 
256  return;
257  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
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 115 of file multiphysics_sys.C.

References _input, _neumann_bcs, _numerical_jacobian_h_values, _numerical_jacobian_h_variables, _physics_list, and GRINS::BCBuilder::build_boundary_conditions().

116  {
117  // Need this to be true because of our overloading of the
118  // mass_residual function.
119  // This is data in FEMSystem. MUST be set before FEMSystem::init_data.
120  use_fixed_solution = true;
121 
122  // Initalize all the variables. We pass this pointer for the system.
123  /* NOTE: We CANNOT fuse this loop with the others. This loop
124  MUST complete first. */
127  for( PhysicsListIter physics_iter = _physics_list.begin();
128  physics_iter != _physics_list.end();
129  physics_iter++ )
130  {
131  (physics_iter->second)->init_variables( this );
132  }
133 
134  libmesh_assert(_input);
136 
137  // If any variables need custom numerical_jacobian_h, we can set those
138  // values now that variable names are all registered with the System
139  for (unsigned int i=0; i != _numerical_jacobian_h_values.size(); ++i)
140  {
141  unsigned int var_num =
142  this->variable_number(_numerical_jacobian_h_variables[i]);
143  this->set_numerical_jacobian_h_for_var
144  (var_num, _numerical_jacobian_h_values[i]);
145  }
146 
147  // Now set time_evolving variables
148  for( PhysicsListIter physics_iter = _physics_list.begin();
149  physics_iter != _physics_list.end();
150  physics_iter++ )
151  {
152  (physics_iter->second)->set_time_evolving_vars( this );
153  }
154 
155  // Set whether the problem we're solving is steady or not
156  // Since the variable is static, just call one Physics class
157  {
158  (_physics_list.begin()->second)->set_is_steady((this->time_solver)->is_steady());
159  }
160 
161  // Next, call parent init_data function to intialize everything.
162  libMesh::FEMSystem::init_data();
163 
164  // After solution has been initialized we can project initial
165  // conditions to it
167  for( PhysicsListIter physics_iter = _physics_list.begin();
168  physics_iter != _physics_list.end();
169  physics_iter++ )
170  {
171  // Initialize builtin IC's for each physics
172  (physics_iter->second)->init_ics( this, ic_function );
173  }
174 
175  if (ic_function.n_subfunctions())
176  {
177  this->project_solution(&ic_function);
178  }
179 
180  // Now do any auxillary initialization required by each Physics
181  for( PhysicsListIter physics_iter = _physics_list.begin();
182  physics_iter != _physics_list.end();
183  physics_iter++ )
184  {
185  (physics_iter->second)->auxiliary_init( *this );
186  }
187 
188  return;
189  }
static void build_boundary_conditions(const GetPot &input, MultiphysicsSystem &system, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Definition: bc_builder.C:41
std::vector< std::string > _numerical_jacobian_h_variables
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
std::vector< libMesh::Real > _numerical_jacobian_h_values
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
const GetPot * _input
Cached for helping build boundary conditions.
bool GRINS::MultiphysicsSystem::mass_residual ( bool  request_jacobian,
libMesh::DiffContext &  context 
)
virtual

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

Definition at line 380 of file multiphysics_sys.C.

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

382  {
383  return this->_general_residual
384  (request_jacobian,
385  context,
388  }
virtual void compute_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:176
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:237
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 360 of file multiphysics_sys.C.

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

362  {
363  return this->_general_residual
364  (request_jacobian,
365  context,
368  }
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:164
virtual void nonlocal_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for scalar variables.
Definition: physics.C:223
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 390 of file multiphysics_sys.C.

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

392  {
393  return this->_general_residual
394  (request_jacobian,
395  context,
398  }
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for scalar variables.
Definition: physics.C:244
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:182
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 330 of file multiphysics_sys.C.

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

332  {
333  return this->_general_residual
334  (request_jacobian,
335  context,
338  }
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
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:146
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 55 of file multiphysics_sys.C.

References _input, _numerical_jacobian_h_values, _numerical_jacobian_h_variables, and _use_numerical_jacobians_only.

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

56  {
57  // Cache this for building boundary condition later
58  _input = &input;
59 
60  // Read options for MultiphysicsSystem first
61  this->verify_analytic_jacobians = input("linear-nonlinear-solver/verify_analytic_jacobians", 0.0 );
62  this->print_solution_norms = input("screen-options/print_solution_norms", false );
63  this->print_solutions = input("screen-options/print_solutions", false );
64  this->print_residual_norms = input("screen-options/print_residual_norms", false );
65 
66  // backwards compatibility with old config files.
68  this->print_residuals = input("screen-options/print_residual", false );
69  if (this->print_residuals)
70  libmesh_deprecated();
71 
72  this->print_residuals = input("screen-options/print_residuals", this->print_residuals );
73  this->print_jacobian_norms = input("screen-options/print_jacobian_norms", false );
74  this->print_jacobians = input("screen-options/print_jacobians", false );
75  this->print_element_solutions = input("screen-options/print_element_solutions", false );
76  this->print_element_residuals = input("screen-options/print_element_residuals", false );
77  this->print_element_jacobians = input("screen-options/print_element_jacobians", false );
78 
79  _use_numerical_jacobians_only = input("linear-nonlinear-solver/use_numerical_jacobians_only", false );
80 
81  numerical_jacobian_h =
82  input("linear-nonlinear-solver/numerical_jacobian_h",
83  numerical_jacobian_h);
84 
85  const unsigned int n_numerical_jacobian_h_values =
86  input.vector_variable_size
87  ("linear-nonlinear-solver/numerical_jacobian_h_values");
88 
89  if (n_numerical_jacobian_h_values !=
90  input.vector_variable_size
91  ("linear-nonlinear-solver/numerical_jacobian_h_variables"))
92  {
93  std::cerr << "Error: found " << n_numerical_jacobian_h_values
94  << " numerical_jacobian_h_values" << std::endl;
95  std::cerr << " but "
96  << input.vector_variable_size
97  ("linear-nonlinear-solver/numerical_jacobian_h_variables")
98  << " numerical_jacobian_h_variables" << std::endl;
99  libmesh_error();
100  }
101 
102  _numerical_jacobian_h_variables.resize(n_numerical_jacobian_h_values);
103  _numerical_jacobian_h_values.resize(n_numerical_jacobian_h_values);
104  for (unsigned int i=0; i != n_numerical_jacobian_h_values; ++i)
105  {
107  input("linear-nonlinear-solver/numerical_jacobian_h_variables",
108  "", i);
110  input("linear-nonlinear-solver/numerical_jacobian_h_values",
111  libMesh::Real(0), i);
112  }
113  }
std::vector< std::string > _numerical_jacobian_h_variables
std::vector< libMesh::Real > _numerical_jacobian_h_values
const GetPot * _input
Cached for helping build boundary conditions.
void GRINS::MultiphysicsSystem::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiAccessor< libMesh::Number > &  param_pointer 
)

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

Definition at line 229 of file multiphysics_sys.C.

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

231  {
232  //Loop over each physics to ask each for the requested parameter
233  for( PhysicsListIter physics_iter = _physics_list.begin();
234  physics_iter != _physics_list.end();
235  physics_iter++ )
236  {
237  (physics_iter->second)->register_parameter( param_name,
238  param_pointer );
239  }
240  }
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer)
Each Physics will register its copy(s) of an independent variable.
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
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 215 of file multiphysics_sys.C.

References _physics_list.

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

217  {
218  for( PhysicsListIter physics_iter = _physics_list.begin();
219  physics_iter != _physics_list.end();
220  physics_iter++ )
221  {
222  (physics_iter->second)->register_postprocessing_vars( input, postprocessing );
223  }
224 
225  return;
226  }
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList.
Definition: var_typedefs.h:62
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 350 of file multiphysics_sys.C.

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

352  {
353  return this->_general_residual
354  (request_jacobian,
355  context,
358  }
virtual void compute_side_constraint_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:158
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
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 314 of file multiphysics_sys.C.

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

316  {
317  bool jacobian_computed = this->apply_neumann_bcs(request_jacobian,
318  context);
319 
320  jacobian_computed = jacobian_computed &&
321  this->_general_residual
322  (request_jacobian,
323  context,
326 
327  return jacobian_computed;
328  }
virtual void compute_side_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
Definition: physics.C:140
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:195
bool apply_neumann_bcs(bool request_jacobian, libMesh::DiffContext &context)
Applies the subset of _neumann_bcs that are active on the current element side.

Member Data Documentation

const GetPot* GRINS::MultiphysicsSystem::_input
private

Cached for helping build boundary conditions.

We can't make a copy because it will muck up the UFO detection amongst other things. So, we keep a raw pointer. We don't own this so we MUST not delete.

Definition at line 210 of file multiphysics_sys.h.

Referenced by init_data(), and read_input_options().

std::vector<SharedPtr<NeumannBCContainer> > GRINS::MultiphysicsSystem::_neumann_bcs
private

Neumann boundary conditions.

Store each NeumannBCContainer for each set of BoundaryIDs and Variables, as specified in the input file. The container knows what BoundaryIDs and Variables it applies to. We use SharedPtr here because libMesh::UniquePtr may still actually be an AutoPtr.

Definition at line 217 of file multiphysics_sys.h.

Referenced by apply_neumann_bcs(), get_neumann_bcs(), and init_data().

std::vector<libMesh::Real> GRINS::MultiphysicsSystem::_numerical_jacobian_h_values
private

Definition at line 204 of file multiphysics_sys.h.

Referenced by init_data(), and read_input_options().

std::vector<std::string> GRINS::MultiphysicsSystem::_numerical_jacobian_h_variables
private

Definition at line 201 of file multiphysics_sys.h.

Referenced by init_data(), and read_input_options().

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 195 of file multiphysics_sys.h.

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

bool GRINS::MultiphysicsSystem::_use_numerical_jacobians_only
private

Definition at line 197 of file multiphysics_sys.h.

Referenced by _general_residual(), apply_neumann_bcs(), and read_input_options().


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

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