Go to the documentation of this file.
20 #ifndef LIBMESH_FE_BASE_H
21 #define LIBMESH_FE_BASE_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/compare_types.h"
26 #include "libmesh/fe_abstract.h"
27 #include "libmesh/fe_transformation_base.h"
28 #include "libmesh/point.h"
29 #include "libmesh/reference_counted_object.h"
30 #include "libmesh/tensor_tools.h"
31 #include "libmesh/type_n_tensor.h"
32 #include "libmesh/vector_value.h"
44 template <
typename T>
class DenseMatrix;
45 template <
typename T>
class DenseVector;
51 template <
typename T>
class NumericVector;
56 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
60 #ifdef LIBMESH_ENABLE_PERIODIC
65 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
66 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
83 template <
typename OutputType>
111 static std::unique_ptr<FEGenericBase>
build (
const unsigned int dim,
129 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
139 static std::unique_ptr<FEGenericBase>
build_InfFE (
const unsigned int dim,
144 #ifdef LIBMESH_ENABLE_AMR
154 const unsigned int variable_number,
166 const Elem * coarse_elem,
168 const unsigned int var,
169 const bool use_old_dof_indices =
false);
179 const Elem * coarse_elem,
181 const bool use_old_dof_indices =
false);
183 #endif // #ifdef LIBMESH_ENABLE_AMR
185 #ifdef LIBMESH_ENABLE_PERIODIC
197 const unsigned int variable_number,
200 #endif // LIBMESH_ENABLE_PERIODIC
206 const std::vector<std::vector<OutputShape>> &
get_phi()
const
214 const std::vector<std::vector<OutputGradient>> &
get_dphi()
const
230 const std::vector<std::vector<OutputDivergence>> &
get_div_phi()
const
238 const std::vector<std::vector<OutputShape>> &
get_dphidx()
const
246 const std::vector<std::vector<OutputShape>> &
get_dphidy()
const
254 const std::vector<std::vector<OutputShape>> &
get_dphidz()
const
282 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
288 const std::vector<std::vector<OutputTensor>> &
get_d2phi()
const
388 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
390 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
443 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
458 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
492 std::unique_ptr<FETransformationBase<OutputType>>
_fe_trans;
497 std::vector<std::vector<OutputShape>>
phi;
502 std::vector<std::vector<OutputGradient>>
dphi;
512 std::vector<std::vector<OutputDivergence>>
div_phi;
517 std::vector<std::vector<OutputShape>>
dphidxi;
532 std::vector<std::vector<OutputShape>>
dphidx;
537 std::vector<std::vector<OutputShape>>
dphidy;
542 std::vector<std::vector<OutputShape>>
dphidz;
545 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
550 std::vector<std::vector<OutputTensor>>
d2phi;
615 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
648 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
655 template <
unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
673 template <
typename OutputType>
689 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
704 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
714 template <
typename OutputType>
722 #endif // LIBMESH_FE_BASE_H
void print_d2phi(std::ostream &os) const
Prints the value of each shape function's second derivatives at each quadrature point.
The constraint matrix storage format.
std::vector< std::vector< OutputShape > > d2phidxdz
Shape function second derivatives in the x-z direction.
const std::vector< std::vector< OutputShape > > & get_d2phideta2() const
We're using a class instead of a typedef to allow forward declarations and future flexibility.
const std::vector< std::vector< OutputShape > > & get_d2phidzeta2() const
static void compute_periodic_constraints(DofConstraints &constraints, DofMap &dof_map, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for meshes with periodic boundary conditions) correspon...
std::unique_ptr< FETransformationBase< OutputType > > _fe_trans
Object that handles computing shape function values, gradients, etc in the physical domain.
bool calculate_dphi
Should we calculate shape function gradients?
const std::vector< std::vector< OutputShape > > & get_d2phidxidzeta() const
std::vector< std::vector< OutputShape > > d2phideta2
Shape function second derivatives in the eta direction.
std::vector< std::vector< OutputShape > > curl_phi
Shape function curl values.
FEGenericBase(const unsigned int dim, const FEType &fet)
Constructor.
OutputType OutputShape
Convenient typedefs for gradients of output, hessians of output, and potentially-complex-valued versi...
virtual ~FEGenericBase()
Destructor.
std::vector< std::vector< OutputShape > > dphidzeta
Shape function derivatives in the zeta direction.
std::vector< std::vector< OutputShape > > dphideta
Shape function derivatives in the eta direction.
const std::vector< std::vector< OutputShape > > & get_dphidxi() const
const std::vector< OutputGradient > & get_dphase() const
std::vector< std::vector< OutputShape > > dphidz
Shape function derivatives in the z direction.
A specific instantiation of the FEBase class.
const std::vector< std::vector< OutputShape > > & get_dphidx() const
std::vector< std::vector< OutputDivergence > > div_phi
Shape function divergence values.
std::vector< std::vector< OutputShape > > dphidy
Shape function derivatives in the y direction.
The libMesh namespace provides an interface to certain functionality in the library.
This class forms the foundation from which generic finite elements may be derived.
This class forms the foundation from which generic finite elements may be derived.
bool calculate_dphiref
Should we calculate reference shape function gradients?
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
bool calculate_phi
Should we calculate shape functions?
FEGenericBase< RealGradient > FEVectorBase
std::vector< std::vector< OutputShape > > d2phidxidzeta
Shape function second derivatives in the xi-zeta direction.
const std::vector< std::vector< OutputShape > > & get_d2phidy2() const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
bool calculate_curl_phi
Should we calculate shape function curls?
const std::vector< std::vector< OutputDivergence > > & get_div_phi() const
std::vector< std::vector< OutputShape > > phi
Shape function values.
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
std::vector< std::vector< OutputShape > > d2phidetadzeta
Shape function second derivatives in the eta-zeta direction.
std::vector< std::vector< OutputShape > > d2phidz2
Shape function second derivatives in the z direction.
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
const std::vector< std::vector< OutputShape > > & get_dphideta() const
std::vector< std::vector< OutputShape > > dphidxi
Shape function derivatives in the xi direction.
std::vector< std::vector< OutputShape > > d2phidzeta2
Shape function second derivatives in the zeta direction.
std::vector< RealGradient > dweight
Used for certain infinite element families: the global derivative of the additional radial weight ,...
const std::vector< std::vector< OutputShape > > & get_d2phidxideta() const
std::vector< std::vector< OutputShape > > d2phidxideta
Shape function second derivatives in the xi-eta direction.
This is the MeshBase class.
FEGenericBase< Real > FEBase
std::vector< std::vector< OutputTensor > > d2phi
Shape function second derivative values.
std::vector< std::vector< OutputShape > > d2phidydz
Shape function second derivatives in the y-z direction.
bool calculate_div_phi
Should we calculate shape function divergences?
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
const std::vector< std::vector< OutputShape > > & get_d2phidxdz() const
TensorTools::DecrementRank< OutputShape >::type OutputDivergence
std::vector< std::vector< OutputGradient > > dphi
Shape function derivative values.
std::vector< std::vector< OutputShape > > d2phidxi2
Shape function second derivatives in the xi direction.
void print_dphi(std::ostream &os) const
Prints the value of each shape function's derivative at each quadrature point.
const unsigned int dim
The dimensionality of the object.
static void coarsened_dof_values(const NumericVector< Number > &global_vector, const DofMap &dof_map, const Elem *coarse_elem, DenseVector< Number > &coarse_dofs, const unsigned int var, const bool use_old_dof_indices=false)
Creates a local projection on coarse_elem, based on the DoF values in global_vector for it's children...
TensorTools::MakeNumber< OutputShape >::type OutputNumber
const std::vector< std::vector< OutputShape > > & get_d2phidxi2() const
const std::vector< RealGradient > & get_Sobolev_dweight() const
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp)
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
TensorTools::DecrementRank< OutputNumber >::type OutputNumberDivergence
const std::vector< std::vector< OutputShape > > & get_d2phidydz() const
The Node constraint storage format.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
const std::vector< std::vector< OutputShape > > & get_dphidy() const
TensorTools::IncrementRank< OutputShape >::type OutputGradient
This class handles the numbering of degrees of freedom on a mesh.
const std::vector< std::vector< OutputShape > > & get_curl_phi() const
static void compute_proj_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
void print_phi(std::ostream &os) const
Prints the value of each shape function at each quadrature point.
std::vector< OutputGradient > dphase
Used for certain infinite element families: the first derivatives of the phase term in global coordin...
const std::vector< std::vector< OutputShape > > & get_d2phidxdy() const
This is the base class from which all geometric element types are derived.
std::vector< std::vector< OutputShape > > dphidx
Shape function derivatives in the x direction.
std::vector< std::vector< OutputShape > > d2phidxdy
Shape function second derivatives in the x-y direction.
void determine_calculations()
Determine which values are to be calculated, for both the FE itself and for the FEMap.
const std::vector< std::vector< OutputShape > > & get_dphidz() const
TensorTools::IncrementRank< OutputNumber >::type OutputNumberGradient
std::vector< std::vector< OutputShape > > d2phidy2
Shape function second derivatives in the y direction.
const std::vector< std::vector< OutputShape > > & get_dphidzeta() const
const std::vector< std::vector< OutputShape > > & get_phi() const
const std::vector< std::vector< OutputShape > > & get_d2phidx2() const
TensorTools::IncrementRank< OutputNumberGradient >::type OutputNumberTensor
const std::vector< std::vector< OutputShape > > & get_d2phidetadzeta() const
This is the base class for point locators.
std::vector< std::vector< OutputShape > > d2phidx2
Shape function second derivatives in the x direction.
virtual void init_base_shape_functions(const std::vector< Point > &qp, const Elem *e)=0
Initialize the data fields for the base of an an infinite element.
bool calculate_d2phi
Should we calculate shape function hessians?
std::vector< Real > weight
Used for certain infinite element families: the additional radial weight in local coordinates,...
TensorTools::IncrementRank< OutputGradient >::type OutputTensor
const std::vector< std::vector< OutputShape > > & get_d2phidz2() const
const std::vector< Real > & get_Sobolev_weight() const