15 #include "../intrules.hpp"
16 #include "../geom.hpp"
47 "unknown BasisType: " << b_type);
55 "invalid nodal BasisType: " <<
Name(b_type));
91 static const char *
Name(
int b_type)
93 static const char *name[] =
95 "Gauss-Legendre",
"Gauss-Lobatto",
"Positive (Bernstein)",
96 "Open uniform",
"Closed uniform",
"Open half uniform",
97 "Serendipity",
"Closed Gauss-Legendre",
98 "Integrated Gauss-Lobatto indicator"
100 return name[
Check(b_type)];
105 static const char ident[]
106 = {
'g',
'G',
'P',
'u',
'U',
'o',
'S',
'c',
'i' };
107 return ident[
Check(b_type)];
125 MFEM_ABORT(
"unknown BasisType identifier");
229 class ElementTransformation;
231 class VectorCoefficient;
232 class MatrixCoefficient;
249 #ifndef MFEM_THREAD_SAFE
450 virtual void GetFaceDofs(
int face,
int **dofs,
int *ndofs)
const;
606 "invalid closed basis type: " << b_type);
614 MFEM_VERIFY(
IsOpenType(b_type),
"invalid open basis type: " << b_type);
632 #ifndef MFEM_THREAD_SAFE
639 "'fe' must be a ScalarFiniteElement");
657 #ifdef MFEM_THREAD_SAFE
800 #ifndef MFEM_THREAD_SAFE
950 {
mfem_error(
"'fe' must be a VectorFiniteElement"); }
986 mutable Vector u_aux, d_aux, d2_aux;
990 Basis *auxiliary_basis;
992 bool scale_integrated;
1027 typedef std::map< int, Array<double*>* > PointsMap;
1028 typedef std::map< int, Array<Basis*>* > BasisMap;
1031 PointsMap points_container;
1032 BasisMap bases_container;
1036 static void CalcMono(
const int p,
const double x,
double *
u);
1037 static void CalcMono(
const int p,
const double x,
double *
u,
double *d);
1039 static void CalcChebyshev(
const int p,
const double x,
double *
u);
1040 static void CalcChebyshev(
const int p,
const double x,
double *
u,
double *d);
1041 static void CalcChebyshev(
const int p,
const double x,
double *
u,
double *d,
1051 static const int *
Binom(
const int p);
1062 const double *
GetPoints(
const int p,
const int btype);
1082 Basis &
GetBasis(
const int p,
const int btype);
1092 { CalcChebyshev(p, x, u); }
1095 static void CalcBasis(
const int p,
const double x,
double *
u,
double *d)
1099 { CalcChebyshev(p, x, u, d); }
1102 static void CalcBasis(
const int p,
const double x,
double *
u,
double *d,
1107 { CalcChebyshev(p, x, u, d, dd); }
1111 {
return pow(x, (
double) p); }
1119 static void CalcBinomTerms(
const int p,
const double x,
const double y,
1124 static void CalcBinomTerms(
const int p,
const double x,
const double y,
1125 double *
u,
double *d);
1144 static void CalcLegendre(
const int p,
const double x,
double *
u);
1145 static void CalcLegendre(
const int p,
const double x,
double *
u,
double *d);
1192 MFEM_ABORT(
"invalid dimension: " << dim);
1202 case 1:
return base;
1203 case 2:
return base*base;
1204 case 3:
return base*base*base;
1205 default: MFEM_ABORT(
"invalid dimension: " << dim);
return -1;
1225 virtual void SetMapType(
const int map_type_);
1253 const int cbtype,
const int obtype,
1258 const int obtype,
const int M,
1269 const bool closed)
const;
Abstract class for all finite elements.
Tensor products of polynomials of order k.
int GetDim() const
Returns the reference space dimension for the finite element.
void trans(const Vector &u, Vector &x)
ScalarFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct ScalarFiniteElement with given.
static void CalcBernstein(const int p, const double x, double *u, double *d)
Compute the values and derivatives of the Bernstein basis functions of order p at coordinate x and st...
Array< int > lex_ordering
Class for an integration rule - an Array of IntegrationPoint.
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
virtual void CalcPhysLinLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
static double CalcDelta(const int p, const double x)
Evaluate a representation of a Delta function at point x.
void ProjectMatrixCoefficient_RT(const double *nk, const Array< int > &d2n, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an RT space.
EvalType
One-dimensional basis evaluation type.
void ProjectMatrixCoefficient_ND(const double *tk, const Array< int > &d2t, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an ND space.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
NodalTensorFiniteElement(const int dims, const int p, const int btype, const DofMapType dmtype)
void LocalInterpolation_RT(const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
Basis(const int p, const double *nodes, EvalType etype=Barycentric)
Create a nodal or positive (Bernstein) basis of degree p.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Base class for vector Coefficients that optionally depend on time and space.
static int GetNodalBasis(int qpt_type)
Return the nodal BasisType corresponding to the Quadrature1D type.
virtual void CalcPhysHessian(ElementTransformation &Trans, DenseMatrix &Hessian) const
Evaluate the Hessian of all shape functions of a scalar finite element in reference space at the give...
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in physical space at the po...
int GetDerivMapType() const
Returns the FiniteElement::DerivType of the element describing how reference function derivatives are...
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
virtual void CalcPhysLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the gi...
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
FiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct FiniteElement with given.
aka "open half" Newton-Cotes
int dim
Dimension of reference space.
TensorBasisElement(const int dims, const int p, const int btype, const DofMapType dmtype)
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Data type dense matrix using column-major storage.
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Array< double > Gt
Transpose of G.
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
~VectorTensorFiniteElement()
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
VectorFiniteElement(int D, Geometry::Type G, int Do, int O, int M, int F=FunctionSpace::Pk)
MapType
Enumeration for MapType: defines how reference functions are mapped to physical space.
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
int Space() const
Returns the type of FunctionSpace on the element.
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
virtual void SetMapType(const int map_type_)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Geometry::Type geom_type
Geometry::Type of the reference element.
DerivType
Enumeration for DerivType: defines which derivative method is implemented.
Intermediate class for finite elements whose basis functions return vector values.
static void CalcLegendre(const int p, const double x, double *u)
void NodalLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get the matrix I that defines nodal interpolation between this element and the refined element fine_f...
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Use change of basis, O(p^2) Evals.
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
const DofToQuad & GetTensorDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode, const bool closed) const
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
const double * ClosedPoints(const int p, const int btype=BasisType::GaussLobatto)
Get coordinates of a closed (GaussLegendre) set of points if degree p.
static int CheckNodal(int b_type)
If the input does not represents a valid nodal BasisType, abort with an error; otherwise return the i...
int GetCurlDim() const
Returns the dimension of the curl for vector-valued finite elements.
Possible basis types. Note that not all elements can use all BasisType(s).
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
int cdim
Dimension of curl for vector-valued basis functions.
int vdim
Vector dimension of vector-valued basis functions.
static void CalcBasis(const int p, const double x, double *u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
virtual void ProjectFromNodes(Vector &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector of values at the finite element nodes and a transformation, compute its projection (ap...
static void CalcBinomTerms(const int p, const double x, const double y, double *u)
Compute the p terms in the expansion of the binomial (x + y)^p and store them in the already allocate...
const double * GetPoints(const int p, const int btype)
Get the coordinates of the points of the given BasisType, btype.
Poly_1D::Basis & obasis1d
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Class for standard nodal finite elements.
void ScalarLocalRestriction(ElementTransformation &Trans, DenseMatrix &R, const ScalarFiniteElement &coarse_fe) const
Get restriction matrix R defined through local L2-projection in the space defined by the coarse_fe...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
VectorTensorFiniteElement(const int dims, const int d, const int p, const int cbtype, const int obtype, const int M, const DofMapType dmtype)
Class for finite elements with basis functions that return scalar values.
static const ScalarFiniteElement & CheckScalarFE(const FiniteElement &fe)
Nodes: x_i = (i+1/2)/n, i=0,...,n-1.
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Poly_1D::Basis & cbasis1d
Fast evaluation of Bernstein polynomials.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
static void ChebyshevPoints(const int p, double *x)
Compute the points for the Chebyshev polynomials of order p and place them in the already allocated x...
Class for computing 1D special polynomials and their associated basis functions.
const int * GetAnisotropicOrders() const
Returns an array containing the anisotropic orders/degrees.
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
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...
class FiniteElement * FE
The FiniteElement that created and owns this object.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
RangeType
Enumeration for range_type and deriv_range_type.
static int VerifyOpen(int b_type)
Ensure that the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary)...
int GetDerivRangeType() const
Returns the FiniteElement::RangeType of the element derivative, either SCALAR or VECTOR.
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
int orders[Geometry::MaxDim]
Anisotropic orders.
virtual ~FiniteElement()
Deconstruct the FiniteElement.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Use barycentric Lagrangian interpolation, O(p) Evals.
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
const double * OpenPoints(const int p, const int btype=BasisType::GaussLegendre)
Get coordinates of an open (GaussLegendre) set of points if degree p.
static void CalcBasis(const int p, const double x, double *u, double *d, double *dd)
Evaluate the values, derivatives and second derivatives of a hierarchical 1D basis at point x...
void ScalarLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get matrix I "Interpolation" defined through local L2-projection in the space defined by the fine_fe...
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
bool HasAnisotropicOrders() const
Returns true if the FiniteElement basis may be using different orders/degrees in different spatial di...
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Array< double > Bt
Transpose of B.
void LocalInterpolation_ND(const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
double p(const Vector &x, double t)
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
Integrated indicator functions (cf. Gerritsma)
Refined tensor products of polynomials of order k.
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
MemoryType
Memory types supported by MFEM.
static bool IsClosedType(int b_type)
Return true if the BasisType of b_type is closed (has Quadrature1D points on the boundary).
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &Hessian) const
Evaluate the Hessians of all shape functions of a scalar finite element in reference space at the giv...
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Base class for Matrix Coefficients that optionally depend on time and space.
static const int * Binom(const int p)
Get a pointer to an array containing the binomial coefficients "p choose k" for k=0,...,p for the given p.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Serendipity basis (squares / cubes)
Nodes: x_i = (i+1)/(n+1), i=0,...,n-1.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
const Array< int > & GetLexicographicOrdering() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Integrated GLL indicator functions.
Array< double > B
Basis functions evaluated at quadrature points.
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const
Get the dofs associated with the given face. *dofs is set to an internal array of the local dofc on t...
aka closed Gauss Legendre
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
void LocalRestriction_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
Class for integration point with weight.
Host memory; using new[] and delete[].
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
Full multidimensional representation which does not use tensor product structure. The ordering of the...
const Poly_1D::Basis & GetBasis1D() const
static bool IsOpenType(int b_type)
Return true if the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary)...
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
void Project_RT(const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the RT basis functions.
int dof
Number of degrees of freedom.
Array< DofToQuad * > dof2quad_array
Container for all DofToQuad objects created by the FiniteElement.
static int VerifyClosed(int b_type)
Ensure that the BasisType of b_type is closed (has Quadrature1D points on the boundary).
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Implements CalcDivShape methods.
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
NodalFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct NodalFiniteElement with given.
void ProjectCurl_2D(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
static int VerifyNodal(int b_type)
Ensure that the BasisType of b_type nodal (satisfies the interpolation property). ...
static int Pow(int base, int dim)
Return base raised to the power dim.
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
void LocalL2Projection_ND(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Nodes: x_i = i/(n-1), i=0,...,n-1.
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
static void CalcBernstein(const int p, const double x, double *u)
Compute the values of the Bernstein basis functions of order p at coordinate x and store the results ...
int GetDerivType() const
Returns the FiniteElement::DerivType of the element describing the spatial derivative method implemen...
static Geometry::Type GetTensorProductGeometry(int dim)
void LocalRestriction_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
Describes the function space on each element.
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
static const VectorFiniteElement & CheckVectorFE(const FiniteElement &fe)
double u(const Vector &xvec)
Implements CalcDShape methods.
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
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...
void Project_ND(const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the ND basis functions.
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
static void CalcBasis(const int p, const double x, double *u, double *d)
Evaluate the values and derivatives of a hierarchical 1D basis at point x.
A Class that defines 1-D numerical quadrature rules on [0,1].
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
static void CalcDBinomTerms(const int p, const double x, const double y, double *d)
Compute the derivatives (w.r.t. x) of the terms in the expansion of the binomial (x + y)^p assuming t...
const IntegrationRule * IntRule
IntegrationRule that defines the quadrature points at which the basis functions of the FE are evaluat...
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
const DofToQuad & GetTensorDofToQuad(const class TensorBasisElement &tb, const IntegrationRule &ir, DofToQuad::Mode mode) const
int order
Order/degree of the shape functions.
Keep count of the number of eval types.
void LocalL2Projection_RT(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
No derivatives implemented.
Implements CalcCurlShape methods.
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const