GRINS-0.7.0
|
Interface with libMesh for solving Multiphysics problems. More...
#include <multiphysics_sys.h>
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 PhysicsList & | get_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 ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_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 which have time varying components. More... | |
virtual bool | side_time_derivative (bool request_jacobian, libMesh::DiffContext &context) |
Boundary contributions to which have time varying components. More... | |
virtual bool | nonlocal_time_derivative (bool request_jacobian, libMesh::DiffContext &context) |
Contributions to 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 which do not have time varying components. More... | |
virtual bool | nonlocal_constraint (bool request_jacobian, libMesh::DiffContext &context) |
Contributions to on SCALAR variables which do not have time varying components. More... | |
virtual bool | damping_residual (bool request_jacobian, libMesh::DiffContext &context) |
Contributions to . More... | |
virtual bool | mass_residual (bool request_jacobian, libMesh::DiffContext &context) |
Contributions to . More... | |
virtual bool | nonlocal_mass_residual (bool request_jacobian, libMesh::DiffContext &context) |
Contributions to 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::Physics > | get_physics (const std::string physics_name) |
SharedPtr< GRINS::Physics > | get_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... | |
Interface with libMesh for solving Multiphysics problems.
MultiphysicsSystem (through libMesh::FEMSystem) solves the following equation:
M = mass matrix u = solution vector F = time derivative
Note that for the nonlinear system that is solved for implicit time stepping is:
_time_derivative correspond to calculating terms for _mass_residual correspond to calculating terms for
Definition at line 82 of file multiphysics_sys.h.
|
private |
Definition at line 225 of file multiphysics_sys.h.
|
private |
Definition at line 224 of file multiphysics_sys.h.
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.
|
inline |
Destructor. Clean up all physics allocations.
Definition at line 92 of file multiphysics_sys.h.
|
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().
|
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().
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().
|
virtual |
Override FEMSystem::build_context in order to use our own AssemblyContext.
Definition at line 191 of file multiphysics_sys.C.
|
virtual |
Definition at line 421 of file multiphysics_sys.C.
References _physics_list.
|
virtual |
Contributions to .
Definition at line 370 of file multiphysics_sys.C.
References _general_residual(), GRINS::Physics::compute_damping_residual_cache(), and GRINS::Physics::damping_residual().
|
virtual |
Element interior contributions to which do not have time varying components. Element interior contributions to 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().
|
virtual |
Element interior contributions to 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().
|
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().
|
inline |
Definition at line 179 of file multiphysics_sys.h.
References _neumann_bcs.
|
inline |
Definition at line 182 of file multiphysics_sys.h.
References _neumann_bcs.
SharedPtr< Physics > GRINS::MultiphysicsSystem::get_physics | ( | const std::string | physics_name | ) |
Definition at line 400 of file multiphysics_sys.C.
References _physics_list.
|
inline |
Definition at line 244 of file multiphysics_sys.h.
References _physics_list.
|
inline |
Definition at line 97 of file multiphysics_sys.h.
References _physics_list.
Referenced by GRINS::OldStyleBCBuilder::build_bcs().
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.
|
virtual |
Context initialization. Calls each physics implementation of init_context()
Definition at line 244 of file multiphysics_sys.C.
References _physics_list.
|
virtual |
System initialization. Calls each physics implementation of init_variables()
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().
|
virtual |
Contributions to .
Definition at line 380 of file multiphysics_sys.C.
References _general_residual(), GRINS::Physics::compute_mass_residual_cache(), and GRINS::Physics::mass_residual().
|
virtual |
Contributions to 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().
|
virtual |
Contributions to 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().
|
virtual |
Contributions to 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().
|
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.
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().
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().
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().
|
virtual |
Boundary contributions to 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().
|
virtual |
Boundary contributions to 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().
|
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().
|
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().
|
private |
Definition at line 204 of file multiphysics_sys.h.
Referenced by init_data(), and read_input_options().
|
private |
Definition at line 201 of file multiphysics_sys.h.
Referenced by init_data(), and read_input_options().
|
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().
|
private |
Definition at line 197 of file multiphysics_sys.h.
Referenced by _general_residual(), apply_neumann_bcs(), and read_input_options().