12#ifndef MFEM_TEMPLATE_FINITE_ELEMENTS
13#define MFEM_TEMPLATE_FINITE_ELEMENTS
35template <
typename real_t>
45 for (
int ip = 0; ip < nip; ip++)
48 for (
int id = 0;
id < dof;
id++)
50 int orig_id = dof_map ? (*dof_map)[id] : id;
51 B[ip+nip*id] = shape(orig_id);
73template <
typename real_t>
84 for (
int ip = 0; ip < nip; ip++)
87 for (
int id = 0;
id < dof;
id++)
89 int orig_id = dof_map ? (*dof_map)[id] : id;
90 for (
int d = 0; d <
dim; d++)
92 G[ip+nip*(d+
dim*id)] = dshape(orig_id, d);
98template <
typename real_t>
108template <Geometry::Type G,
int P>
117 static const int degree = P;
118 static const int dofs = P+1;
120 static const bool tensor_prod =
true;
121 static const int dofs_1d = P+1;
157 MFEM_ASSERT(h1_fec,
"invalid FiniteElementCollection");
162 template <
typename real_t>
167 template <
typename real_t>
181 static const int degree = P;
182 static const int dofs = ((P + 1)*(P + 2))/2;
184 static const bool tensor_prod =
false;
215 MFEM_ASSERT(h1_fec,
"invalid FiniteElementCollection");
220 template <
typename real_t>
234 static const int degree = P;
235 static const int dofs = (P+1)*(P+1);
237 static const bool tensor_prod =
true;
238 static const int dofs_1d = P+1;
276 MFEM_ASSERT(h1_fec,
"invalid FiniteElementCollection");
281 template <
typename real_t>
286 template <
typename real_t>
300 static const int degree = P;
301 static const int dofs = ((P + 1)*(P + 2)*(P + 3))/6;
303 static const bool tensor_prod =
false;
334 MFEM_ASSERT(h1_fec,
"invalid FiniteElementCollection");
339 template <
typename real_t>
353 static const int degree = P;
354 static const int dofs = (P+1)*(P+1)*(P+1);
356 static const bool tensor_prod =
true;
357 static const int dofs_1d = P+1;
396 MFEM_ASSERT(h1_fec,
"invalid FiniteElementCollection");
401 template <
typename real_t>
406 template <
typename real_t>
417template <
Geometry::Type G,
int P,
typename L2_FE_type,
typename L2Pos_FE_type,
443 my_fe =
new L2Pos_FE_type(P);
449 my_fe =
new L2_FE_type(P, pt_type);
462 MFEM_ASSERT(l2_fec,
"invalid FiniteElementCollection");
469 template <
typename real_t>
474 template <
typename real_t>
483template <Geometry::Type G,
int P>
490 Geometry::SEGMENT,P,L2_SegmentElement,L2Pos_SegmentElement,P+1,true>
507 L2Pos_TriangleElement,((P+1)*(P+2))/2,false>
524 L2Pos_QuadrilateralElement,(P+1)*(P+1),true>
541 L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6,false>
558 L2Pos_HexahedronElement,(P+1)*(P+1)*(P+1),true>
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
@ GaussLobatto
Closed type.
@ GaussLegendre
Open type.
@ Positive
Bernstein polynomials.
Data type dense matrix using column-major storage.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Abstract class for all finite elements.
int GetDim() const
Returns the reference space dimension for the finite element.
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube.
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square.
Arbitrary order H1 elements in 1D utilizing the Bernstein basis.
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle.
Arbitrary order H1-conforming (continuous) finite elements.
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
const Array< int > * GetDofMap() const
const Array< int > * my_dof_map
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
void Init(const parameter_type type_)
const FiniteElement * my_fe
H1_FiniteElement(const FiniteElementCollection &fec)
const FiniteElement * my_fe
H1_FiniteElement(const FiniteElementCollection &fec)
void Init(const parameter_type type_)
const Array< int > * GetDofMap() const
const Array< int > * my_dof_map
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
const FiniteElement * my_fe
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
void Init(const parameter_type type_)
const Array< int > * GetDofMap() const
H1_FiniteElement(const FiniteElementCollection &fec)
const Array< int > * my_dof_map
const Array< int > * GetDofMap() const
H1_FiniteElement(const FiniteElementCollection &fec)
void Init(const parameter_type type_)
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
const FiniteElement * my_fe
H1_FiniteElement(const FiniteElementCollection &fec)
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
void Init(const parameter_type type_)
const Array< int > * GetDofMap() const
const FiniteElement * my_fe
Arbitrary order H1 elements in 3D on a cube.
Arbitrary order H1 elements in 2D on a square.
Arbitrary order H1 elements in 1D.
Arbitrary order H1 elements in 3D on a tetrahedron.
Arbitrary order H1 elements in 2D on a triangle.
A class to initialize the size of a Tensor.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube.
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square.
Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment.
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle.
Arbitrary order "L2-conforming" discontinuous finite elements.
L2_FiniteElement(const FiniteElementCollection &fec)
base_class::parameter_type parameter_type
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
L2_FiniteElement_base< Geometry::SEGMENT, P, L2_SegmentElement, L2Pos_SegmentElement, P+1, true > base_class
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
L2_FiniteElement(const FiniteElementCollection &fec)
base_class::parameter_type parameter_type
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
L2_FiniteElement(const FiniteElementCollection &fec)
L2_FiniteElement_base< Geometry::SQUARE, P, L2_QuadrilateralElement, L2Pos_QuadrilateralElement,(P+1) *(P+1), true > base_class
base_class::parameter_type parameter_type
base_class::parameter_type parameter_type
L2_FiniteElement(const FiniteElementCollection &fec)
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
L2_FiniteElement_base< Geometry::TETRAHEDRON, P, L2_TetrahedronElement, L2Pos_TetrahedronElement,((P+1) *(P+2) *(P+3))/6, false > base_class
base_class::parameter_type parameter_type
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
L2_FiniteElement(const FiniteElementCollection &fec)
const FiniteElement * my_fe
static const bool tensor_prod
L2_FiniteElement_base(const FiniteElementCollection &fec)
const Array< int > * GetDofMap() const
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
L2_FiniteElement_base(const parameter_type type)
const FiniteElement * my_fe_1d
static const Geometry::Type geom
void Init(const parameter_type type_)
Arbitrary order L2 elements in 3D on a cube.
Arbitrary order L2 elements in 2D on a square.
Arbitrary order L2 elements in 1D on a segment.
Arbitrary order L2 elements in 3D on a tetrahedron.
Arbitrary order L2 elements in 2D on a triangle.
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
void CalcGradTensor(const FiniteElement &fe, const IntegrationRule &ir, real_t *G, const Array< int > *dof_map=NULL)
store gradient matrix G for each integration point on the reference element. For tensor product evalu...
void CalcShapeMatrix(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, const Array< int > *dof_map=NULL)
Store mass-like matrix B for each integration point on the reference element. For tensor product eval...
void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, real_t *G, const Array< int > *dof_map)