12 #ifndef MFEM_TEMPLATE_FINITE_ELEMENTS 13 #define MFEM_TEMPLATE_FINITE_ELEMENTS 15 #include "../config/tconfig.hpp" 35 template <
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);
73 template <
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);
98 template <
typename real_t>
100 real_t *B, real_t *G,
const Array<int> *dof_map)
108 template <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>
417 template <
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>
483 template <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>
573 #endif // MFEM_TEMPLATE_FINITE_ELEMENTS void Init(const parameter_type type_)
Abstract class for all finite elements.
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
const Array< int > * GetDofMap() const
const Array< int > * my_dof_map
const FiniteElement * my_fe_1d
int GetNPoints() const
Returns the number of the points in the integration rule.
const FiniteElement * my_fe
Class for an integration rule - an Array of IntegrationPoint.
L2_FiniteElement(const FiniteElementCollection &fec)
static const Geometry::Type geom
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle.
L2_FiniteElement_base< Geometry::SEGMENT, P, L2_SegmentElement, L2Pos_SegmentElement, P+1, true > base_class
L2_FiniteElement(const FiniteElementCollection &fec)
Arbitrary order L2 elements in 2D on a triangle.
Arbitrary order L2 elements in 3D on a cube.
L2_FiniteElement(const FiniteElementCollection &fec)
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment.
void Init(const parameter_type type_)
H1_FiniteElement(const FiniteElementCollection &fec)
Arbitrary order H1 elements in 2D on a square.
const Array< int > * my_dof_map
Data type dense matrix using column-major storage.
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
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 Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
H1_FiniteElement(const FiniteElementCollection &fec)
const Array< int > * my_dof_map
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
const FiniteElement * my_fe
base_class::parameter_type parameter_type
L2_FiniteElement(const FiniteElementCollection &fec)
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
void Init(const parameter_type type_)
L2_FiniteElement_base(const FiniteElementCollection &fec)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square.
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Arbitrary order L2 elements in 1D on a segment.
const FiniteElement * my_fe_1d
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...
H1_FiniteElement(const FiniteElementCollection &fec)
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...
static const bool tensor_prod
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
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Arbitrary order H1 elements in 3D on a tetrahedron.
base_class::parameter_type parameter_type
const FiniteElement * my_fe
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
base_class::parameter_type parameter_type
Arbitrary order L2 elements in 3D on a tetrahedron.
void Init(const parameter_type type_)
int GetDim() const
Returns the reference space dimension for the finite element.
L2_FiniteElement_base(const parameter_type type)
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Arbitrary order H1 elements in 3D on a cube.
const Array< int > * GetDofMap() const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
A class to initialize the size of a Tensor.
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
const Array< int > * GetDofMap() const
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
H1_FiniteElement(const FiniteElementCollection &fec)
H1_FiniteElement(const FiniteElementCollection &fec)
Arbitrary order H1 elements in 2D on a triangle.
Arbitrary order L2 elements in 2D on a square.
Arbitrary order H1 elements in 1D.
Arbitrary order H1 elements in 1D utilizing the Bernstein basis.
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square.
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube.
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
const FiniteElement * my_fe_1d
const Array< int > * GetDofMap() const
const FiniteElement * my_fe
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Arbitrary order H1-conforming (continuous) finite elements.
const Array< int > * GetDofMap() const
base_class::parameter_type parameter_type
const Array< int > * GetDofMap() const
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
void Init(const parameter_type type_)
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...
L2_FiniteElement_base< Geometry::SQUARE, P, L2_QuadrilateralElement, L2Pos_QuadrilateralElement,(P+1) *(P+1), true > base_class
void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, real_t *G, const Array< int > *dof_map)
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle.
void Init(const parameter_type type_)
L2_FiniteElement(const FiniteElementCollection &fec)
Arbitrary order "L2-conforming" discontinuous finite elements.