GRINS-0.6.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... | |
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::ParameterMultiPointer< libMesh::Number > ¶m_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 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 | 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... | |
std::tr1::shared_ptr< GRINS::Physics > | get_physics (const std::string physics_name) |
std::tr1::shared_ptr< 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) |
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 |
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 81 of file multiphysics_sys.h.
|
private |
Definition at line 191 of file multiphysics_sys.h.
|
private |
Definition at line 190 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 40 of file multiphysics_sys.C.
GRINS::MultiphysicsSystem::~MultiphysicsSystem | ( | ) |
Destructor. Clean up all physics allocations.
Definition at line 49 of file multiphysics_sys.C.
|
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().
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().
|
virtual |
Override FEMSystem::build_context in order to use our own AssemblyContext.
Definition at line 159 of file multiphysics_sys.C.
|
virtual |
Definition at line 373 of file multiphysics_sys.C.
References _physics_list.
|
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 302 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 272 of file multiphysics_sys.C.
References _general_residual(), GRINS::Physics::compute_element_time_derivative_cache(), and GRINS::Physics::element_time_derivative().
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().
|
inline |
Definition at line 201 of file multiphysics_sys.h.
References _physics_list.
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.
|
virtual |
Context initialization. Calls each physics implementation of init_context()
Definition at line 212 of file multiphysics_sys.C.
References _physics_list.
|
virtual |
System initialization. Calls each physics implementation of init_variables()
Definition at line 88 of file multiphysics_sys.C.
References _physics_list.
|
virtual |
Contributions to .
Definition at line 332 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 322 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 342 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 292 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 60 of file multiphysics_sys.C.
References _use_numerical_jacobians_only.
Referenced by GRINS::Simulation::init_multiphysics_system().
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().
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().
|
virtual |
Boundary contributions to 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().
|
virtual |
Boundary contributions to 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().
|
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().
|
private |
Definition at line 183 of file multiphysics_sys.h.
Referenced by _general_residual(), and read_input_options().