GRINS-0.8.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 void assembly (bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false)
 Override FEMSystem::assembly. More...
 
virtual void reinit ()
 Override FEMSystem::reinit. 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 &)
 
typedef void(GRINS::Physics::* CacheFuncType) (AssemblyContext &)
 

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...
 
libMesh::UniquePtr< libMesh::System::Constraint > _constraint
 Constraint application object. 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 77 of file multiphysics_sys.h.

Member Typedef Documentation

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

Definition at line 227 of file multiphysics_sys.h.

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

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

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

Destructor. Clean up all physics allocations.

Definition at line 87 of file multiphysics_sys.h.

87 {};

Member Function Documentation

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

Definition at line 306 of file multiphysics_sys.C.

References _physics_list, _use_numerical_jacobians_only, GRINS::CachedValues::clear(), and GRINS::AssemblyContext::get_cached_values().

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().

310  {
311  AssemblyContext& c = libMesh::cast_ref<AssemblyContext&>(context);
312 
313  bool compute_jacobian = true;
314  if( !request_jacobian || _use_numerical_jacobians_only ) compute_jacobian = false;
315 
316  CachedValues & cache = c.get_cached_values();
317 
318  // Now compute cache for this element
319  for( PhysicsListIter physics_iter = _physics_list.begin();
320  physics_iter != _physics_list.end();
321  physics_iter++ )
322  {
323  // shared_ptr gets confused by operator->*
324  ((*(physics_iter->second)).*cachefunc)( c );
325  }
326 
327  // Loop over each physics and compute their contributions
328  for( PhysicsListIter physics_iter = _physics_list.begin();
329  physics_iter != _physics_list.end();
330  physics_iter++ )
331  {
332  if(c.has_elem())
333  {
334  if( (physics_iter->second)->enabled_on_elem( &c.get_elem() ) )
335  {
336  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c );
337  }
338  }
339  else
340  {
341  ((*(physics_iter->second)).*resfunc)( compute_jacobian, c );
342  }
343  }
344 
345  // We need to clear out the cache when we're done so we don't interfere
346  // with other residual functions
347  cache.clear();
348 
349  // TODO: Need to think about the implications of this because there might be some
350  // TODO: jacobian terms we don't want to compute for efficiency reasons
351  return compute_jacobian;
352  }
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 500 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().

502  {
503  AssemblyContext& assembly_context =
504  libMesh::cast_ref<AssemblyContext&>( context );
505 
506  std::vector<BoundaryID> ids;
507  assembly_context.side_boundary_ids(ids);
508 
509  bool compute_jacobian = request_jacobian;
510  if( !request_jacobian || _use_numerical_jacobians_only ) compute_jacobian = false;
511 
512  for( std::vector<BoundaryID>::const_iterator it = ids.begin();
513  it != ids.end(); it++ )
514  {
515  BoundaryID bc_id = *it;
516 
517  libmesh_assert_not_equal_to(bc_id, libMesh::BoundaryInfo::invalid_id);
518 
519  // Retreive the NeumannBCContainers that are active on the current bc_id
520  std::vector<SharedPtr<NeumannBCContainer> > active_neumann_bcs;
521  this->get_active_neumann_bcs( bc_id, _neumann_bcs, active_neumann_bcs );
522 
523  if( !active_neumann_bcs.empty() )
524  {
525  typedef std::vector<SharedPtr<NeumannBCContainer> >::iterator BCIt;
526 
527  for( BCIt container = active_neumann_bcs.begin(); container < active_neumann_bcs.end(); ++container )
528  {
529  SharedPtr<NeumannBCAbstract>& func = (*container)->get_func();
530 
531  const FEVariablesBase& var = (*container)->get_fe_var();
532 
533  func->eval_flux( compute_jacobian, assembly_context,
534  var.neumann_bc_sign(), Physics::is_axisymmetric() );
535  }
536  }
537  } // end loop over boundary ids
538 
539  return compute_jacobian;
540  }
static bool is_axisymmetric()
Definition: physics.h:132
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::assembly ( bool  get_residual,
bool  get_jacobian,
bool  apply_heterogeneous_constraints = false,
bool  apply_no_constraints = false 
)
virtual

Override FEMSystem::assembly.

This allows us to insert things like preassembly().

Definition at line 269 of file multiphysics_sys.C.

References _physics_list.

Referenced by GRINS::UnsteadyVisualization::output_residual().

273  {
274  // First do any preassembly that the Physics requires (which by default is none)
275  for( PhysicsListIter physics_iter = _physics_list.begin();
276  physics_iter != _physics_list.end();
277  physics_iter++ )
278  (physics_iter->second)->preassembly(*this);
279 
280  // Now do the assembly
281  libMesh::FEMSystem::assembly(get_residual,get_jacobian,
282  apply_heterogeneous_constraints,
283  apply_no_constraints);
284  }
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
void GRINS::MultiphysicsSystem::attach_physics_list ( PhysicsList  physics_list)

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

Definition at line 52 of file multiphysics_sys.C.

References _physics_list.

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

53  {
54  _physics_list = physics_list;
55  }
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 198 of file multiphysics_sys.C.

Referenced by GRINSTesting::SpectroscopicAbsorptionTest::elem_qoi_derivatives().

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

Definition at line 471 of file multiphysics_sys.C.

References _physics_list.

475  {
476  for( PhysicsListIter physics_iter = _physics_list.begin();
477  physics_iter != _physics_list.end();
478  physics_iter++ )
479  {
480  // Only compute if physics is active on current subdomain or globally
481  if( (physics_iter->second)->enabled_on_elem( &context.get_elem() ) )
482  {
483  (physics_iter->second)->compute_postprocessed_quantity( quantity_index, context, point, value );
484  }
485  }
486  return;
487  }
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 420 of file multiphysics_sys.C.

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

422  {
423  return this->_general_residual
424  (request_jacobian,
425  context,
428  }
virtual void compute_damping_residual_cache(AssemblyContext &)
Definition: physics.h:224
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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
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 390 of file multiphysics_sys.C.

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

392  {
393  return this->_general_residual
394  (request_jacobian,
395  context,
398  }
virtual void compute_element_constraint_cache(AssemblyContext &)
Definition: physics.h:218
virtual void element_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for element interiors.
Definition: physics.h:186
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 354 of file multiphysics_sys.C.

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

356  {
357  return this->_general_residual
358  (request_jacobian,
359  context,
362  }
virtual void compute_element_time_derivative_cache(AssemblyContext &)
Definition: physics.h:212
virtual void element_time_derivative(bool, AssemblyContext &)
Time dependent part(s) of physics for element interiors.
Definition: physics.h:174
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 489 of file multiphysics_sys.C.

Referenced by apply_neumann_bcs().

492  {
493  // Manually writing the loop since std::copy_if is C++11 only
494  for( std::vector<SharedPtr<NeumannBCContainer> >::const_iterator it = neumann_bcs.begin();
495  it < neumann_bcs.end(); ++it )
496  if( (*it)->has_bc_id( bc_id ) )
497  active_neumann_bcs.push_back( *it );
498  }
std::vector<SharedPtr<NeumannBCContainer> >& GRINS::MultiphysicsSystem::get_neumann_bcs ( )
inline

Definition at line 186 of file multiphysics_sys.h.

References _neumann_bcs.

187  { 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 189 of file multiphysics_sys.h.

References _neumann_bcs.

190  { 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 450 of file multiphysics_sys.C.

References _physics_list.

451  {
452  if( _physics_list.find( physics_name ) == _physics_list.end() )
453  {
454  std::cerr << "Error: Could not find physics " << physics_name << std::endl;
455  libmesh_error();
456  }
457 
458  return _physics_list[physics_name];
459  }
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 246 of file multiphysics_sys.h.

References _physics_list.

247  {
248  libmesh_assert(_physics_list.find( physics_name ) != _physics_list.end());
249 
250  return _physics_list.find(physics_name)->second;
251  }
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 92 of file multiphysics_sys.h.

References _physics_list.

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

93  { 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 461 of file multiphysics_sys.C.

References _physics_list.

462  {
463  bool has_physics = false;
464 
465  if( _physics_list.find(physics_name) != _physics_list.end() )
466  has_physics = true;
467 
468  return has_physics;
469  }
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 254 of file multiphysics_sys.C.

References _physics_list.

Referenced by GRINSTesting::SpectroscopicAbsorptionTest::elem_qoi_derivatives().

255  {
256  AssemblyContext& c = libMesh::cast_ref<AssemblyContext&>(context);
257 
258  //Loop over each physics to initialize relevant variable structures for assembling system
259  for( PhysicsListIter physics_iter = _physics_list.begin();
260  physics_iter != _physics_list.end();
261  physics_iter++ )
262  {
263  (physics_iter->second)->init_context( c );
264  }
265 
266  return;
267  }
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 117 of file multiphysics_sys.C.

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

118  {
119  // Need this to be true because of our overloading of the
120  // mass_residual function.
121  // This is data in FEMSystem. MUST be set before FEMSystem::init_data.
122  use_fixed_solution = true;
123 
124  // Initalize all the variables. We pass this pointer for the system.
125  /* NOTE: We CANNOT fuse this loop with the others. This loop
126  MUST complete first. */
129  for( PhysicsListIter physics_iter = _physics_list.begin();
130  physics_iter != _physics_list.end();
131  physics_iter++ )
132  {
133  (physics_iter->second)->init_variables( this );
134  }
135 
136  libmesh_assert(_input);
138 
139  this->_constraint =
141 
142  this->attach_constraint_object(*this->_constraint);
143 
144  // If any variables need custom numerical_jacobian_h, we can set those
145  // values now that variable names are all registered with the System
146  for (unsigned int i=0; i != _numerical_jacobian_h_values.size(); ++i)
147  {
148  unsigned int var_num =
149  this->variable_number(_numerical_jacobian_h_variables[i]);
150  this->set_numerical_jacobian_h_for_var
151  (var_num, _numerical_jacobian_h_values[i]);
152  }
153 
154  // Now set time_evolving variables
155  for( PhysicsListIter physics_iter = _physics_list.begin();
156  physics_iter != _physics_list.end();
157  physics_iter++ )
158  {
159  (physics_iter->second)->set_time_evolving_vars( this );
160  }
161 
162  // Set whether the problem we're solving is steady or not
163  // Since the variable is static, just call one Physics class
164  {
165  (_physics_list.begin()->second)->set_is_steady((this->time_solver)->is_steady());
166  }
167 
168  // Next, call parent init_data function to intialize everything.
169  libMesh::FEMSystem::init_data();
170 
171  // After solution has been initialized we can project initial
172  // conditions to it
174  for( PhysicsListIter physics_iter = _physics_list.begin();
175  physics_iter != _physics_list.end();
176  physics_iter++ )
177  {
178  // Initialize builtin IC's for each physics
179  (physics_iter->second)->init_ics( this, ic_function );
180  }
181 
182  if (ic_function.n_subfunctions())
183  {
184  this->project_solution(&ic_function);
185  }
186 
187  // Now do any auxillary initialization required by each Physics
188  for( PhysicsListIter physics_iter = _physics_list.begin();
189  physics_iter != _physics_list.end();
190  physics_iter++ )
191  {
192  (physics_iter->second)->auxiliary_init( *this );
193  }
194 
195  return;
196  }
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.
libMesh::UniquePtr< libMesh::System::Constraint > _constraint
Constraint application object.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
static libMesh::UniquePtr< libMesh::System::Constraint > build_constraint_object(const GetPot &input, MultiphysicsSystem &system)
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 430 of file multiphysics_sys.C.

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

432  {
433  return this->_general_residual
434  (request_jacobian,
435  context,
438  }
virtual void compute_mass_residual_cache(AssemblyContext &)
Definition: physics.h:226
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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 410 of file multiphysics_sys.C.

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

412  {
413  return this->_general_residual
414  (request_jacobian,
415  context,
418  }
virtual void nonlocal_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for scalar variables.
Definition: physics.h:194
virtual void compute_nonlocal_constraint_cache(AssemblyContext &)
Definition: physics.h:222
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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 440 of file multiphysics_sys.C.

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

442  {
443  return this->_general_residual
444  (request_jacobian,
445  context,
448  }
virtual void nonlocal_mass_residual(bool, AssemblyContext &)
Mass matrix part(s) for scalar variables.
Definition: physics.h:206
virtual void compute_nonlocal_mass_residual_cache(AssemblyContext &)
Definition: physics.h:228
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
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 380 of file multiphysics_sys.C.

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

382  {
383  return this->_general_residual
384  (request_jacobian,
385  context,
388  }
virtual void compute_nonlocal_time_derivative_cache(AssemblyContext &)
Definition: physics.h:216
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual void nonlocal_time_derivative(bool, AssemblyContext &)
Time dependent part(s) of physics for scalar variables.
Definition: physics.h:182
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 57 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(), and GRINSTesting::SystemHelper::setup_multiphysics_system().

58  {
59  // Cache this for building boundary condition later
60  _input = &input;
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 
87  const unsigned int n_numerical_jacobian_h_values =
88  input.vector_variable_size
89  ("linear-nonlinear-solver/numerical_jacobian_h_values");
90 
91  if (n_numerical_jacobian_h_values !=
92  input.vector_variable_size
93  ("linear-nonlinear-solver/numerical_jacobian_h_variables"))
94  {
95  std::cerr << "Error: found " << n_numerical_jacobian_h_values
96  << " numerical_jacobian_h_values" << std::endl;
97  std::cerr << " but "
98  << input.vector_variable_size
99  ("linear-nonlinear-solver/numerical_jacobian_h_variables")
100  << " numerical_jacobian_h_variables" << std::endl;
101  libmesh_error();
102  }
103 
104  _numerical_jacobian_h_variables.resize(n_numerical_jacobian_h_values);
105  _numerical_jacobian_h_values.resize(n_numerical_jacobian_h_values);
106  for (unsigned int i=0; i != n_numerical_jacobian_h_values; ++i)
107  {
109  input("linear-nonlinear-solver/numerical_jacobian_h_variables",
110  "", i);
112  input("linear-nonlinear-solver/numerical_jacobian_h_values",
113  libMesh::Real(0), i);
114  }
115  }
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 236 of file multiphysics_sys.C.

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

238  {
239  //Loop over each physics to ask each for the requested parameter
240  for( PhysicsListIter physics_iter = _physics_list.begin();
241  physics_iter != _physics_list.end();
242  physics_iter++ )
243  {
244  (physics_iter->second)->register_parameter( param_name,
245  param_pointer );
246  }
247 
248  for (unsigned int bc=0; bc<_neumann_bcs.size(); bc++)
249  _neumann_bcs[bc]->get_func()->register_parameter(param_name, param_pointer);
250  }
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
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 222 of file multiphysics_sys.C.

References _physics_list.

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

224  {
225  for( PhysicsListIter physics_iter = _physics_list.begin();
226  physics_iter != _physics_list.end();
227  physics_iter++ )
228  {
229  (physics_iter->second)->register_postprocessing_vars( input, postprocessing );
230  }
231 
232  return;
233  }
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
void GRINS::MultiphysicsSystem::reinit ( )
virtual

Override FEMSystem::reinit.

This will allow each Physics to reinit things internally that need it, such as point locators.

Definition at line 286 of file multiphysics_sys.C.

References _physics_list, and GRINS::CompositeQoI::reinit().

Referenced by GRINSTesting::IntegratedFunctionTest::reinit_through_system().

287  {
288  // First call Parent
289  FEMSystem::reinit();
290 
291  // Now do per Physics reinit (which by default is none)
292  for( PhysicsListIter physics_iter = _physics_list.begin();
293  physics_iter != _physics_list.end();
294  physics_iter++ )
295  (physics_iter->second)->reinit(*this);
296 
297  // And now reinit the QoI
298  if (this->qoi.size() > 0)
299  {
300  libMesh::DifferentiableQoI* diff_qoi = this->get_qoi();
301  CompositeQoI* qoi = libMesh::cast_ptr<CompositeQoI*>(diff_qoi);
302  qoi->reinit(*this);
303  }
304  }
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 400 of file multiphysics_sys.C.

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

402  {
403  return this->_general_residual
404  (request_jacobian,
405  context,
408  }
virtual void compute_side_constraint_cache(AssemblyContext &)
Definition: physics.h:220
virtual void side_constraint(bool, AssemblyContext &)
Constraint part(s) of physics for boundaries of elements on the domain boundary.
Definition: physics.h:190
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 364 of file multiphysics_sys.C.

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

366  {
367  bool jacobian_computed = this->apply_neumann_bcs(request_jacobian,
368  context);
369 
370  jacobian_computed = jacobian_computed &&
371  this->_general_residual
372  (request_jacobian,
373  context,
376 
377  return jacobian_computed;
378  }
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
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual void compute_side_time_derivative_cache(AssemblyContext &)
Definition: physics.h:214
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

libMesh::UniquePtr<libMesh::System::Constraint> GRINS::MultiphysicsSystem::_constraint
private

Constraint application object.

Definition at line 222 of file multiphysics_sys.h.

Referenced by init_data().

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 212 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 219 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 206 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 203 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 197 of file multiphysics_sys.h.

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

bool GRINS::MultiphysicsSystem::_use_numerical_jacobians_only
private

Definition at line 199 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 Tue Dec 19 2017 12:47:31 for GRINS-0.8.0 by  doxygen 1.8.9.1