21 #ifndef __mast_fe_base_h__ 22 #define __mast_fe_base_h__ 28 #include "libmesh/elem.h" 29 #include "libmesh/fe_base.h" 30 #include "libmesh/quadrature.h" 36 class SystemInitialization;
71 const std::vector<libMesh::Point>* pts =
nullptr);
83 bool if_calculate_dphi);
89 virtual const std::vector<Real>&
96 virtual const std::vector<libMesh::Point>&
102 virtual const std::vector<std::vector<Real> >&
105 virtual const std::vector<std::vector<libMesh::RealVectorValue> >&
108 virtual const std::vector<std::vector<libMesh::RealTensorValue>>&
111 virtual const std::vector<Real>&
114 virtual const std::vector<Real>&
117 virtual const std::vector<Real>&
120 virtual const std::vector<Real>&
123 virtual const std::vector<Real>&
126 virtual const std::vector<Real>&
129 virtual const std::vector<Real>&
132 virtual const std::vector<Real>&
135 virtual const std::vector<Real>&
138 virtual const std::vector<libMesh::RealVectorValue>&
141 virtual const std::vector<libMesh::RealVectorValue>&
144 virtual const std::vector<libMesh::RealVectorValue>&
147 virtual const std::vector<std::vector<Real> >&
150 virtual const std::vector<std::vector<Real> >&
153 virtual const std::vector<std::vector<Real> >&
160 virtual const std::vector<libMesh::Point>&
167 virtual const std::vector<libMesh::Point>&
171 virtual const std::vector<libMesh::Point>&
174 virtual const libMesh::QBase&
194 #endif // __mast_fe_base_h__ virtual const std::vector< Real > & get_detadz() const
virtual const std::vector< Real > & get_dzetadx() const
std::vector< libMesh::Point > _global_normals
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdzeta() const
virtual const std::vector< Real > & get_JxW() const
virtual const std::vector< Real > & get_dxidy() const
virtual const std::vector< Real > & get_dxidz() const
virtual const libMesh::QBase & get_qrule() const
FEBase(const MAST::SystemInitialization &sys)
virtual const std::vector< Real > & get_detadx() const
virtual const std::vector< libMesh::Point > & get_xyz() const
physical location of the quadrature point in the global coordinate system for the reference element ...
bool _init_second_order_derivatives
virtual const std::vector< std::vector< Real > > & get_dphidzeta() const
virtual const std::vector< libMesh::Point > & get_qpoints() const
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdeta() const
virtual const std::vector< Real > & get_dzetady() const
std::vector< libMesh::Point > _local_normals
std::vector< libMesh::Point > _qpoints
virtual void init_for_side(const MAST::GeomElem &elem, unsigned int s, bool if_calculate_dphi)
Initializes the quadrature and finite element for element side integration.
virtual const std::vector< Real > & get_dzetadz() const
const MAST::SystemInitialization & _sys
void set_evaluate_second_order_derivatives(bool f)
sets the flag for evaluation of second order derivative
virtual const std::vector< std::vector< Real > > & get_dphideta() const
void set_extra_quadrature_order(int n)
this is used, in addition to libMesh::System::extra_quadrature_order to set the quadrature rule...
virtual const std::vector< std::vector< Real > > & get_dphidxi() const
virtual unsigned int n_shape_functions() const
virtual const std::vector< libMesh::Point > & get_normals_for_reference_coordinate() const
normals defined in the global coordinate system for the reference element
std::vector< libMesh::Point > _global_xyz
virtual const std::vector< Real > & get_detady() const
virtual const std::vector< Real > & get_dxidx() const
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
virtual void init(const MAST::GeomElem &elem, bool init_grads, const std::vector< libMesh::Point > *pts=nullptr)
Initializes the quadrature and finite element for element volume integration.
virtual const std::vector< std::vector< libMesh::RealVectorValue > > & get_dphi() const
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdxi() const
virtual const std::vector< std::vector< libMesh::RealTensorValue > > & get_d2phi() const
virtual const std::vector< std::vector< Real > > & get_phi() const
libMesh::FEType get_fe_type() const
virtual const std::vector< libMesh::Point > & get_normals_for_local_coordinate() const
normals defined in the coordinate system for the local reference element.
const MAST::GeomElem * _elem
unsigned int _extra_quadrature_order