26 #ifndef GRINS_MULTIPHYSICS_SYS_H
27 #define GRINS_MULTIPHYSICS_SYS_H
33 #include "grins_config.h"
37 #include "libmesh/fem_system.h"
39 #ifdef GRINS_HAVE_GRVY
49 class EquationSystems;
52 template <
typename Scalar>
59 template <
typename Scalar>
87 const std::string& name,
88 const unsigned int number );
114 (
const std::string & param_name,
118 virtual libMesh::AutoPtr<libMesh::DiffContext>
build_context();
121 virtual void init_context( libMesh::DiffContext &context );
128 libMesh::DiffContext& context );
132 libMesh::DiffContext& context );
136 libMesh::DiffContext& context );
141 libMesh::DiffContext& context );
145 libMesh::DiffContext& context );
149 libMesh::DiffContext& context );
153 libMesh::DiffContext& context );
157 libMesh::DiffContext& context );
160 bool has_physics(
const std::string physics_name )
const;
162 std::tr1::shared_ptr<GRINS::Physics>
get_physics(
const std::string physics_name );
164 std::tr1::shared_ptr<GRINS::Physics>
get_physics(
const std::string physics_name )
const;
168 const libMesh::Point& point,
169 libMesh::Real& value );
171 #ifdef GRINS_USE_GRVY_TIMERS
172 void attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer );
185 #ifdef GRINS_USE_GRVY_TIMERS
186 GRVY::GRVY_Timer_Class* _timer;
195 libMesh::DiffContext& context,
210 #endif // GRINS_MULTIPHYSICS_SYS_H
void(GRINS::Physics::* ResFuncType)(bool, AssemblyContext &, CachedValues &)
bool _use_numerical_jacobians_only
bool has_physics(const std::string physics_name) const
Query to check if a particular physics has been enabled.
void(GRINS::Physics::* CacheFuncType)(const AssemblyContext &, CachedValues &)
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
std::tr1::shared_ptr< GRINS::Physics > get_physics(const std::string physics_name)
virtual void init_data()
System initialization. Calls each physics implementation of init_variables()
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
virtual bool side_constraint(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to which do not have time varying components.
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.
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime.
virtual bool nonlocal_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables which have time varying components.
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
~MultiphysicsSystem()
Destructor. Clean up all physics allocations.
virtual bool nonlocal_mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables.
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Element interior contributions to which have time varying components.
virtual libMesh::AutoPtr< libMesh::DiffContext > build_context()
Override FEMSystem::build_context in order to use our own AssemblyContext.
virtual bool nonlocal_constraint(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables which do not have time varying components.
Interface with libMesh for solving Multiphysics problems.
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
virtual bool element_constraint(bool request_jacobian, libMesh::DiffContext &context)
virtual bool mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to .
MultiphysicsSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Constructor. Will be called by libMesh only.
virtual bool side_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to which have time varying components.
virtual void init_context(libMesh::DiffContext &context)
Context initialization. Calls each physics implementation of init_context()
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled.