21 #ifndef __mast__heat_conduction_elem_base__ 22 #define __mast__heat_conduction_elem_base__ 30 class ElementPropertyCardBase;
31 class BoundaryConditionBase;
32 class FEMOperatorMatrix;
33 template <
typename ValType>
class FieldFunction;
102 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
111 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
126 const unsigned int s,
142 const unsigned int s,
151 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
159 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
167 const unsigned int s,
169 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
182 const unsigned int s,
201 const unsigned int s,
220 const unsigned int s,
232 const unsigned int s,
254 const unsigned int s,
275 const unsigned int s,
286 const unsigned int s,
307 const unsigned int s,
328 const unsigned int s,
358 const unsigned int s,
380 const unsigned int dim,
382 std::vector<MAST::FEMOperatorMatrix>& dBmat);
393 #endif // __mast__heat_conduction_elem_base__
virtual void surface_convection_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface convection.
void _initialize_mass_fem_operator(const unsigned int qp, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
When mass = false, initializes the FEM operator matrix to the shape functions as .
virtual void volume_heat_source_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
virtual void internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
virtual void surface_flux_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface flux on element side s...
virtual void surface_flux_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element side s.
virtual void volume_heat_source_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source on element volumetric domain...
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
HeatConductionElementBase(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
void _initialize_fem_gradient_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, std::vector< MAST::FEMOperatorMatrix > &dBmat)
For mass = true, the FEM operator matrix is initilized to the weak form of the Laplacian ...
void side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual
virtual void surface_radiation_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void surface_convection_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface convection.
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void surface_convection_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
virtual void volume_heat_source_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
void side_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the side external force contribution to system residual
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
virtual ~HeatConductionElementBase()
const MAST::GeomElem & elem() const
void volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual
void volume_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the volume external force contribution to system residual
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
boundary velocity contribution of volume external force.
virtual void surface_flux_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
This element implements the Galerkin discretization of the heat conduction problem with the flux pro...
virtual void velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
virtual void surface_radiation_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element side...
virtual void velocity_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the damping force contribution to system residual
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
virtual void surface_radiation_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface radiation flux on side s...
const MAST::ElementPropertyCardBase & _property
element property
virtual void velocity_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
This is the base class for elements that implement calculation of finite element quantities over the ...