26 #ifndef GRINS_MULTIPHYSICS_SYS_H
27 #define GRINS_MULTIPHYSICS_SYS_H
33 #include "grins_config.h"
38 #include "libmesh/fem_system.h"
45 class EquationSystems;
48 template <
typename Scalar>
55 template <
typename Scalar>
83 const std::string& name,
84 const unsigned int number );
113 (
const std::string & param_name,
117 virtual libMesh::UniquePtr<libMesh::DiffContext>
build_context();
120 virtual void init_context( libMesh::DiffContext &context );
124 virtual void assembly(
bool get_residual,
126 bool apply_heterogeneous_constraints =
false,
127 bool apply_no_constraints =
false );
139 libMesh::DiffContext& context );
143 libMesh::DiffContext& context );
147 libMesh::DiffContext& context );
152 libMesh::DiffContext& context );
156 libMesh::DiffContext& context );
160 libMesh::DiffContext& context );
164 libMesh::DiffContext& context );
168 libMesh::DiffContext& context );
172 libMesh::DiffContext& context );
175 bool has_physics(
const std::string physics_name )
const;
177 SharedPtr<GRINS::Physics>
get_physics(
const std::string physics_name );
179 SharedPtr<GRINS::Physics>
get_physics(
const std::string physics_name )
const;
183 const libMesh::Point& point,
184 libMesh::Real& value );
231 libMesh::DiffContext & context,
237 const std::vector<SharedPtr<NeumannBCContainer> >& neumann_bcs,
238 std::vector<SharedPtr<NeumannBCContainer> >& active_neumann_bcs );
242 libMesh::DiffContext& context );
255 #endif // GRINS_MULTIPHYSICS_SYS_H
bool _use_numerical_jacobians_only
bool has_physics(const std::string physics_name) const
Query to check if a particular physics has been enabled.
const std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs() const
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
std::vector< std::string > _numerical_jacobian_h_variables
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions.
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
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.
libMesh::boundary_id_type BoundaryID
More descriptive name of the type used for boundary ids.
virtual bool damping_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to .
SharedPtr< GRINS::Physics > get_physics(const std::string physics_name)
libMesh::UniquePtr< libMesh::System::Constraint > _constraint
Constraint application object.
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.
std::vector< SharedPtr< NeumannBCContainer > > & get_neumann_bcs()
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.
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.
virtual bool nonlocal_mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to on SCALAR variables.
void(GRINS::Physics::* ResFuncType)(bool, AssemblyContext &)
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Element interior contributions to which have time varying components.
virtual void reinit()
Override FEMSystem::reinit.
virtual libMesh::UniquePtr< 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.
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.
Interface with libMesh for solving Multiphysics problems.
std::map< std::string, SharedPtr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
std::vector< libMesh::Real > _numerical_jacobian_h_values
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false)
Override FEMSystem::assembly.
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 .
bool apply_neumann_bcs(bool request_jacobian, libMesh::DiffContext &context)
Applies the subset of _neumann_bcs that are active on the current element side.
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()
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.
void(GRINS::Physics::* CacheFuncType)(AssemblyContext &)
const GetPot * _input
Cached for helping build boundary conditions.