21 #ifndef __mast__level_set_elem_base__ 22 #define __mast__level_set_elem_base__ 30 class ElementPropertyCardBase;
31 class BoundaryConditionBase;
32 class FEMOperatorMatrix;
33 template <
typename ValType>
class FieldFunction;
101 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
110 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
134 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
143 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
191 const libMesh::Point& p,
194 const std::vector<MAST::FEMOperatorMatrix>& dBmat,
204 const std::vector<MAST::FEMOperatorMatrix>& dBmat,
213 const unsigned int qp,
214 const std::vector<MAST::FEMOperatorMatrix>& dB_mat,
220 const unsigned int qp,
238 std::vector<MAST::FEMOperatorMatrix>& dBmat);
261 #endif // __mast__level_set_elem_base__ virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
void _calculate_dxidX(const MAST::FEBase &fe, const unsigned int qp, RealMatrixX &dxi_dX)
const RealVectorX & sol(bool if_sens=false) const
const MAST::FieldFunction< Real > * _phi_vel
element property
virtual bool velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
void _tau(const MAST::FEBase &fe, unsigned int qp, const MAST::FEMOperatorMatrix &Bmat, const std::vector< MAST::FEMOperatorMatrix > &dBmat, const RealVectorX &vel, RealMatrixX &tau)
initializes the tau operator
const MAST::GeomElem & elem() const
Real perimeter()
Approximates the integral of the Dirac delta function to approximate the perimeter.
Real perimeter_sensitivity()
computes the partial derivative of the integral of the Dirac delta function using the solution and se...
RealVectorX _ref_sol
reference solution for reinitialization of the level set
Real volume_boundary_velocity_on_side(unsigned int s)
MAST::AssemblyBase & assembly()
void set_reference_solution_for_initialization(const RealVectorX &sol)
For reinitialization to , the solution before initialization is used to calculate the source and velo...
bool volume_external_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the volume external force contribution to system residual
virtual ~LevelSetElementBase()
void set_propagation_mode(bool f)
This can operate in one of two modes: propagation of level set given Vn, or reinitialization of level...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
bool side_external_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the side external force contribution to system residual
virtual bool velocity_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
Matrix< Real, Dynamic, 1 > RealVectorX
void set_velocity_function(const MAST::FieldFunction< Real > &vel)
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
bool _if_propagation
this can operate in one of two modes: propagation of level set given Vn, or reinitialization of level...
void _initialize_fem_operators(const unsigned int qp, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat, std::vector< MAST::FEMOperatorMatrix > &dBmat)
When mass = false, initializes the FEM operator matrix to the shape functions as ...
bool 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 _velocity_and_source(const unsigned int qp, const libMesh::Point &p, const Real t, const MAST::FEMOperatorMatrix &Bmat, const std::vector< MAST::FEMOperatorMatrix > &dBmat, RealVectorX &vel, Real &source)
calculates the velocity at the quadrature point
LevelSetElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem)
Constructor.
virtual bool internal_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the internal force contribution to system residual
bool 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
void _dc_operator(const MAST::FEBase &fe, const unsigned int qp, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const RealVectorX &vel, Real &dc)
This is the base class for elements that implement calculation of finite element quantities over the ...