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");
226 class ElementTransformation;
228 class VectorCoefficient;
229 class MatrixCoefficient;
246 #ifndef MFEM_THREAD_SAFE 448 virtual void GetFaceDofs(
int face,
int **dofs,
int *ndofs)
const;
605 "invalid closed basis type: " << b_type);
613 MFEM_VERIFY(
IsOpenType(b_type),
"invalid open basis type: " << b_type);
633 "'fe' must be a ScalarFiniteElement");
772 Vector &shape)
const override;
780 #ifndef MFEM_THREAD_SAFE 930 {
mfem_error(
"'fe' must be a VectorFiniteElement"); }
965 mutable Vector u_aux, d_aux, d2_aux;
969 Basis *auxiliary_basis;
971 bool scale_integrated;
1006 typedef std::map< int, Array<double*>* > PointsMap;
1007 typedef std::map< int, Array<Basis*>* > BasisMap;
1010 PointsMap points_container;
1011 BasisMap bases_container;
1015 static void CalcMono(
const int p,
const double x,
double *
u);
1016 static void CalcMono(
const int p,
const double x,
double *
u,
double *d);
1018 static void CalcChebyshev(
const int p,
const double x,
double *
u);
1019 static void CalcChebyshev(
const int p,
const double x,
double *
u,
double *d);
1020 static void CalcChebyshev(
const int p,
const double x,
double *
u,
double *d,
1030 static const int *
Binom(
const int p);
1041 const double *
GetPoints(
const int p,
const int btype);
1061 Basis &
GetBasis(
const int p,
const int btype);
1071 { CalcChebyshev(
p, x,
u); }
1079 static void CalcBasis(
const int p,
const double x,
double *
u,
double *d)
1083 { CalcChebyshev(
p, x,
u, d); }
1091 static void CalcBasis(
const int p,
const double x,
double *
u,
double *d,
1096 { CalcChebyshev(
p, x,
u, d, dd); }
1106 {
return pow(x, (
double)
p); }
1114 static void CalcBinomTerms(
const int p,
const double x,
const double y,
1119 static void CalcBinomTerms(
const int p,
const double x,
const double y,
1120 double *
u,
double *d);
1151 static void CalcLegendre(
const int p,
const double x,
double *
u);
1152 static void CalcLegendre(
const int p,
const double x,
double *
u,
double *d);
1157 extern MFEM_EXPORT Poly_1D
poly1d;
1198 MFEM_ABORT(
"invalid dimension: " <<
dim);
1208 case 1:
return base;
1209 case 2:
return base*base;
1210 case 3:
return base*base*base;
1211 default: MFEM_ABORT(
"invalid dimension: " <<
dim);
return -1;
1236 void SetMapType(
const int map_type_)
override;
1264 const int cbtype,
const int obtype,
1269 const int obtype,
const int M,
1285 dof2quad_array_open);
1292 const IntegrationPoint &pt, Vector &x);
Abstract class for all finite elements.
virtual void CalcPhysLinLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Tensor products of polynomials of order k.
void trans(const Vector &u, Vector &x)
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
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...
const class FiniteElement * FE
The FiniteElement that created and owns this object.
Array< int > lex_ordering
Class for an integration rule - an Array of IntegrationPoint.
static void CalcBernstein(const int p, const double x, Vector &u, Vector &d)
Compute the values and derivatives of the Bernstein basis functions of order p at coordinate x and st...
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
static double CalcDelta(const int p, const double x)
Evaluate a representation of a Delta function at point x.
EvalType
One-dimensional basis evaluation type.
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
NodalTensorFiniteElement(const int dims, const int p, const int btype, const DofMapType dmtype)
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...
Basis(const int p, const double *nodes, EvalType etype=Barycentric)
Create a nodal or positive (Bernstein) basis of degree p.
void LocalL2Projection_RT(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
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.
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const override
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
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 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.
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
int Space() const
Returns the type of FunctionSpace on the element.
void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const override
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.
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...
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...
FiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct FiniteElement with given.
static void CalcBernstein(const int p, const double x, Vector &u)
Compute the values of the Bernstein basis functions of order p at coordinate x and store the results ...
aka "open half" Newton-Cotes
int GetCurlDim() const
Returns the dimension of the curl for vector-valued finite elements.
int dim
Dimension of reference space.
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const override
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
TensorBasisElement(const int dims, const int p, const int btype, const DofMapType dmtype)
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...
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.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
virtual ~VectorTensorFiniteElement()
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
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)
bool HasAnisotropicOrders() const
Returns true if the FiniteElement basis may be using different orders/degrees in different spatial di...
MapType
Enumeration for MapType: defines how reference functions are mapped to physical space.
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
static void CalcBasis(const int p, const double x, Vector &u, Vector &d, Vector &dd)
Evaluate the values, derivatives and second derivatives of a hierarchical 1D basis at point x...
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
Poly_1D::Basis & obasis1d
static const DofToQuad & GetTensorDofToQuad(const FiniteElement &fe, const IntegrationRule &ir, DofToQuad::Mode mode, const Poly_1D::Basis &basis, bool closed, Array< DofToQuad *> &dof2quad_array)
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)
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Use change of basis, O(p^2) Evals.
void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const override
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
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...
Possible basis types. Note that not all elements can use all BasisType(s).
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.
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
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...
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.
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Class for standard nodal finite elements.
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)
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...
Class for finite elements with basis functions that return scalar values.
void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const override
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
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...
void LocalRestriction_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
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...
const int * GetAnisotropicOrders() const
Returns an array containing the anisotropic orders/degrees.
void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const override
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
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...
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Class for computing 1D special polynomials and their associated basis functions.
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...
void ProjectCurl_2D(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
const Poly_1D::Basis & GetBasis1D() const
void LocalL2Projection_ND(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
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...
RangeType
Enumeration for range_type and deriv_range_type.
static void CalcBasis(const int p, const double x, Vector &u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
static int VerifyOpen(int b_type)
Ensure that the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary)...
int orders[Geometry::MaxDim]
Anisotropic orders.
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.
virtual ~FiniteElement()
Deconstruct the FiniteElement.
Use barycentric Lagrangian interpolation, O(p) 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...
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...
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
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.
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...
const double * OpenPoints(const int p, const int btype=BasisType::GaussLegendre)
Get coordinates of an open (GaussLegendre) set of points if degree p.
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.
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 EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
const Array< int > & GetLexicographicOrdering() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const override
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
double * GetData() const
Return a pointer to the beginning of the Vector data.
void LocalInterpolation_RT(const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
Array< double > Bt
Transpose of B.
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
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)
int GetDim() const
Returns the reference space dimension for the finite element.
Refined tensor products of polynomials of order k.
int GetDerivType() const
Returns the FiniteElement::DerivType of the element describing the spatial derivative method implemen...
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
static void CalcBasis(const int p, const double x, Vector &u, Vector &d)
Evaluate the values and derivatives of a hierarchical 1D basis at point x.
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).
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.
void LocalInterpolation_ND(const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
Base class for Matrix Coefficients that optionally depend on time and space.
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...
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.
Serendipity basis (squares / cubes)
Nodes: x_i = (i+1)/(n+1), i=0,...,n-1.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
void SetMapType(const int map_type_) override
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Integrated GLL indicator functions.
int GetDerivMapType() const
Returns the FiniteElement::DerivType of the element describing how reference function derivatives are...
Array< double > B
Basis functions evaluated at quadrature points.
aka closed Gauss Legendre
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
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 ...
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.
Full multidimensional representation which does not use tensor product structure. The ordering of the...
void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const override
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.
static bool IsOpenType(int b_type)
Return true if the BasisType of b_type is open (doesn't have Quadrature1D points on the boundary)...
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
int dof
Number of degrees of freedom.
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
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).
Implements CalcDivShape methods.
void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const override
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
NodalFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct NodalFiniteElement with given.
int GetDerivRangeType() const
Returns the FiniteElement::RangeType of the element derivative, either SCALAR or VECTOR.
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
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...
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
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.
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 ...
static Geometry::Type GetTensorProductGeometry(int dim)
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
void LocalRestriction_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Describes the function space on each element.
static const VectorFiniteElement & CheckVectorFE(const FiniteElement &fe)
void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const override
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
double u(const Vector &xvec)
Implements CalcDShape methods.
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 GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the 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...
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.
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].
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void ScalarLocalL2Restriction(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...
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 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...
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
int order
Order/degree of the shape functions.
Keep count of the number of eval types.
No derivatives implemented.
Implements CalcCurlShape methods.
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...