MFEM  v4.6.0
Finite element discretization library
fe_base.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_FE_BASE
13 #define MFEM_FE_BASE
14 
15 #include "../intrules.hpp"
16 #include "../geom.hpp"
17 #include "../doftrans.hpp"
18 
19 #include <map>
20 
21 namespace mfem
22 {
23 
24 /// Possible basis types. Note that not all elements can use all BasisType(s).
25 class BasisType
26 {
27 public:
28  enum
29  {
30  Invalid = -1,
31  GaussLegendre = 0, ///< Open type
32  GaussLobatto = 1, ///< Closed type
33  Positive = 2, ///< Bernstein polynomials
34  OpenUniform = 3, ///< Nodes: x_i = (i+1)/(n+1), i=0,...,n-1
35  ClosedUniform = 4, ///< Nodes: x_i = i/(n-1), i=0,...,n-1
36  OpenHalfUniform = 5, ///< Nodes: x_i = (i+1/2)/n, i=0,...,n-1
37  Serendipity = 6, ///< Serendipity basis (squares / cubes)
38  ClosedGL = 7, ///< Closed GaussLegendre
39  IntegratedGLL = 8, ///< Integrated GLL indicator functions
40  NumBasisTypes = 9 /**< Keep track of maximum types to prevent
41  hard-coding */
42  };
43  /** @brief If the input does not represents a valid BasisType, abort with an
44  error; otherwise return the input. */
45  static int Check(int b_type)
46  {
47  MFEM_VERIFY(0 <= b_type && b_type < NumBasisTypes,
48  "unknown BasisType: " << b_type);
49  return b_type;
50  }
51  /** @brief If the input does not represents a valid nodal BasisType, abort
52  with an error; otherwise return the input. */
53  static int CheckNodal(int b_type)
54  {
55  MFEM_VERIFY(Check(b_type) != Positive && b_type != IntegratedGLL,
56  "invalid nodal BasisType: " << Name(b_type));
57  return b_type;
58  }
59  /** @brief Get the corresponding Quadrature1D constant, when that makes
60  sense; otherwise return Quadrature1D::Invalid. */
61  static int GetQuadrature1D(int b_type)
62  {
63  switch (b_type)
64  {
67  case Positive: return Quadrature1D::ClosedUniform; // <-----
72  case ClosedGL: return Quadrature1D::ClosedGL;
74  }
75  return Quadrature1D::Invalid;
76  }
77  /// Return the nodal BasisType corresponding to the Quadrature1D type.
78  static int GetNodalBasis(int qpt_type)
79  {
80  switch (qpt_type)
81  {
87  case Quadrature1D::ClosedGL: return ClosedGL;
88  }
89  return Invalid;
90  }
91  /// Check and convert a BasisType constant to a string identifier.
92  static const char *Name(int b_type)
93  {
94  static const char *name[] =
95  {
96  "Gauss-Legendre", "Gauss-Lobatto", "Positive (Bernstein)",
97  "Open uniform", "Closed uniform", "Open half uniform",
98  "Serendipity", "Closed Gauss-Legendre",
99  "Integrated Gauss-Lobatto indicator"
100  };
101  return name[Check(b_type)];
102  }
103  /// Check and convert a BasisType constant to a char basis identifier.
104  static char GetChar(int b_type)
105  {
106  static const char ident[]
107  = { 'g', 'G', 'P', 'u', 'U', 'o', 'S', 'c', 'i' };
108  return ident[Check(b_type)];
109  }
110  /// Convert char basis identifier to a BasisType constant.
111  static int GetType(char b_ident)
112  {
113  switch (b_ident)
114  {
115  case 'g': return GaussLegendre;
116  case 'G': return GaussLobatto;
117  case 's': return GaussLobatto;
118  case 'P': return Positive;
119  case 'u': return OpenUniform;
120  case 'U': return ClosedUniform;
121  case 'o': return OpenHalfUniform;
122  case 'S': return Serendipity;
123  case 'c': return ClosedGL;
124  case 'i': return IntegratedGLL;
125  }
126  MFEM_ABORT("unknown BasisType identifier");
127  return -1;
128  }
129 };
130 
131 /** @brief Structure representing the matrices/tensors needed to evaluate (in
132  reference space) the values, gradients, divergences, or curls of a
133  FiniteElement at a the quadrature points of a given IntegrationRule. */
134 /** Object of this type are typically created and owned by the respective
135  FiniteElement object. */
137 {
138 public:
139  /// The FiniteElement that created and owns this object.
140  /** This pointer is not owned. */
141  const class FiniteElement *FE;
142 
143  /** @brief IntegrationRule that defines the quadrature points at which the
144  basis functions of the #FE are evaluated. */
145  /** This pointer is not owned. */
147 
148  /// Type of data stored in the arrays #B, #Bt, #G, and #Gt.
149  enum Mode
150  {
151  /** @brief Full multidimensional representation which does not use tensor
152  product structure. The ordering of the degrees of freedom is as
153  defined by #FE */
155 
156  /** @brief Tensor product representation using 1D matrices/tensors with
157  dimensions using 1D number of quadrature points and degrees of
158  freedom. */
159  /** When representing a vector-valued FiniteElement, two DofToQuad objects
160  are used to describe the "closed" and "open" 1D basis functions. */
162  };
163 
164  /// Describes the contents of the #B, #Bt, #G, and #Gt arrays, see #Mode.
166 
167  /** @brief Number of degrees of freedom = number of basis functions. When
168  #mode is TENSOR, this is the 1D number. */
169  int ndof;
170 
171  /** @brief Number of quadrature points. When #mode is TENSOR, this is the 1D
172  number. */
173  int nqpt;
174 
175  /// Basis functions evaluated at quadrature points.
176  /** The storage layout is column-major with dimensions:
177  - #nqpt x #ndof, for scalar elements, or
178  - #nqpt x dim x #ndof, for vector elements,
179 
180  where
181 
182  - dim = dimension of the finite element reference space when #mode is
183  FULL, and dim = 1 when #mode is TENSOR. */
185 
186  /// Transpose of #B.
187  /** The storage layout is column-major with dimensions:
188  - #ndof x #nqpt, for scalar elements, or
189  - #ndof x #nqpt x dim, for vector elements. */
191 
192  /** @brief Gradients/divergences/curls of basis functions evaluated at
193  quadrature points. */
194  /** The storage layout is column-major with dimensions:
195  - #nqpt x dim x #ndof, for scalar elements, or
196  - #nqpt x #ndof, for H(div) vector elements, or
197  - #nqpt x cdim x #ndof, for H(curl) vector elements,
198 
199  where
200 
201  - dim = dimension of the finite element reference space when #mode is
202  FULL, and 1 when #mode is TENSOR,
203  - cdim = 1/1/3 in 1D/2D/3D, respectively, when #mode is FULL, and cdim =
204  1 when #mode is TENSOR. */
206 
207  /// Transpose of #G.
208  /** The storage layout is column-major with dimensions:
209  - #ndof x #nqpt x dim, for scalar elements, or
210  - #ndof x #nqpt, for H(div) vector elements, or
211  - #ndof x #nqpt x cdim, for H(curl) vector elements. */
213 };
214 
215 /// Describes the function space on each element
217 {
218 public:
219  enum
220  {
221  Pk, ///< Polynomials of order k
222  Qk, ///< Tensor products of polynomials of order k
223  rQk ///< Refined tensor products of polynomials of order k
224  };
225 };
226 
227 class ElementTransformation;
228 class Coefficient;
229 class VectorCoefficient;
230 class MatrixCoefficient;
231 
232 /// Abstract class for all finite elements.
234 {
235 protected:
236  int dim; ///< Dimension of reference space
237  int vdim; ///< Vector dimension of vector-valued basis functions
238  int cdim; ///< Dimension of curl for vector-valued basis functions
239  Geometry::Type geom_type; ///< Geometry::Type of the reference element
242  mutable
243  int dof, ///< Number of degrees of freedom
244  order; ///< Order/degree of the shape functions
245  mutable int orders[Geometry::MaxDim]; ///< Anisotropic orders
247 #ifndef MFEM_THREAD_SAFE
248  mutable DenseMatrix vshape; // Dof x Dim
249 #endif
250  /// Container for all DofToQuad objects created by the FiniteElement.
251  /** Multiple DofToQuad objects may be needed when different quadrature rules
252  or different DofToQuad::Mode are used. */
254 
255 public:
256  /// Enumeration for range_type and deriv_range_type
258 
259  /** @brief Enumeration for MapType: defines how reference functions are
260  mapped to physical space.
261 
262  A reference function \f$ \hat u(\hat x) \f$ can be mapped to a function
263  \f$ u(x) \f$ on a general physical element in following ways:
264  - \f$ x = T(\hat x) \f$ is the image of the reference point \f$ \hat x \f$
265  - \f$ J = J(\hat x) \f$ is the Jacobian matrix of the transformation T
266  - \f$ w = w(\hat x) = det(J) \f$ is the transformation weight factor for square J
267  - \f$ w = w(\hat x) = det(J^t J)^{1/2} \f$ is the transformation weight factor in general
268  */
269  enum MapType
270  {
271  UNKNOWN_MAP_TYPE = -1, /**< Used to distinguish an unset MapType variable
272  from the known values below. */
273  VALUE, /**< For scalar fields; preserves point values
274  \f$ u(x) = \hat u(\hat x) \f$ */
275  INTEGRAL, /**< For scalar fields; preserves volume integrals
276  \f$ u(x) = (1/w) \hat u(\hat x) \f$ */
277  H_DIV, /**< For vector fields; preserves surface integrals of the
278  normal component \f$ u(x) = (J/w) \hat u(\hat x) \f$ */
279  H_CURL /**< For vector fields; preserves line integrals of the
280  tangential component
281  \f$ u(x) = J^{-t} \hat u(\hat x) \f$ (square J),
282  \f$ u(x) = J(J^t J)^{-1} \hat u(\hat x) \f$ (general J) */
283  };
284 
285  /** @brief Enumeration for DerivType: defines which derivative method
286  is implemented.
287 
288  Each FiniteElement class implements up to one type of derivative. The
289  value returned by GetDerivType() indicates which derivative method is
290  implemented.
291  */
293  {
294  NONE, ///< No derivatives implemented
295  GRAD, ///< Implements CalcDShape methods
296  DIV, ///< Implements CalcDivShape methods
297  CURL ///< Implements CalcCurlShape methods
298  };
299 
300  /** @brief Construct FiniteElement with given
301  @param D Reference space dimension
302  @param G Geometry type (of type Geometry::Type)
303  @param Do Number of degrees of freedom in the FiniteElement
304  @param O Order/degree of the FiniteElement
305  @param F FunctionSpace type of the FiniteElement
306  */
307  FiniteElement(int D, Geometry::Type G, int Do, int O,
308  int F = FunctionSpace::Pk);
309 
310  /// Returns the reference space dimension for the finite element
311  int GetDim() const { return dim; }
312 
313  /// Returns the vector dimension for vector-valued finite elements
314  int GetVDim() const { return vdim; }
315 
316  /// Returns the dimension of the curl for vector-valued finite elements
317  int GetCurlDim() const { return cdim; }
318 
319  /// Returns the Geometry::Type of the reference element
321 
322  /// Returns the number of degrees of freedom in the finite element
323  int GetDof() const { return dof; }
324 
325  /** @brief Returns the order of the finite element. In the case of
326  anisotropic orders, returns the maximum order. */
327  int GetOrder() const { return order; }
328 
329  /** @brief Returns true if the FiniteElement basis *may be using* different
330  orders/degrees in different spatial directions. */
331  bool HasAnisotropicOrders() const { return orders[0] != -1; }
332 
333  /// Returns an array containing the anisotropic orders/degrees.
334  const int *GetAnisotropicOrders() const { return orders; }
335 
336  /// Returns the type of FunctionSpace on the element.
337  int Space() const { return func_space; }
338 
339  /// Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
340  int GetRangeType() const { return range_type; }
341 
342  /** @brief Returns the FiniteElement::RangeType of the element derivative, either
343  SCALAR or VECTOR. */
344  int GetDerivRangeType() const { return deriv_range_type; }
345 
346  /** @brief Returns the FiniteElement::MapType of the element describing how reference
347  functions are mapped to physical space, one of {VALUE, INTEGRAL
348  H_DIV, H_CURL}. */
349  int GetMapType() const { return map_type; }
350 
351  /** @brief Returns the FiniteElement::DerivType of the element describing the
352  spatial derivative method implemented, one of {NONE, GRAD,
353  DIV, CURL}. */
354  int GetDerivType() const { return deriv_type; }
355 
356  /** @brief Returns the FiniteElement::DerivType of the element describing how
357  reference function derivatives are mapped to physical space, one of {VALUE,
358  INTEGRAL, H_DIV, H_CURL}. */
359  int GetDerivMapType() const { return deriv_map_type; }
360 
361  /** @brief Evaluate the values of all shape functions of a scalar finite
362  element in reference space at the given point @a ip. */
363  /** The size (#dof) of the result Vector @a shape must be set in advance. */
364  virtual void CalcShape(const IntegrationPoint &ip,
365  Vector &shape) const = 0;
366 
367  /** @brief Evaluate the values of all shape functions of a scalar finite
368  element in physical space at the point described by @a Trans. */
369  /** The size (#dof) of the result Vector @a shape must be set in advance. */
370  void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const;
371 
372  /** @brief Evaluate the gradients of all shape functions of a scalar finite
373  element in reference space at the given point @a ip. */
374  /** Each row of the result DenseMatrix @a dshape contains the derivatives of
375  one shape function. The size (#dof x #dim) of @a dshape must be set in
376  advance. */
377  virtual void CalcDShape(const IntegrationPoint &ip,
378  DenseMatrix &dshape) const = 0;
379 
380  /** @brief Evaluate the gradients of all shape functions of a scalar finite
381  element in physical space at the point described by @a Trans. */
382  /** Each row of the result DenseMatrix @a dshape contains the derivatives of
383  one shape function. The size (#dof x SDim) of @a dshape must be set in
384  advance, where SDim >= #dim is the physical space dimension as described
385  by @a Trans. */
386  void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const;
387 
388  /// Get a const reference to the nodes of the element
389  const IntegrationRule & GetNodes() const { return Nodes; }
390 
391  // virtual functions for finite elements on vector spaces
392 
393  /** @brief Evaluate the values of all shape functions of a *vector* finite
394  element in reference space at the given point @a ip. */
395  /** Each row of the result DenseMatrix @a shape contains the components of
396  one vector shape function. The size (#dof x #dim) of @a shape must be set
397  in advance. */
398  virtual void CalcVShape(const IntegrationPoint &ip,
399  DenseMatrix &shape) const;
400 
401  /** @brief Evaluate the values of all shape functions of a *vector* finite
402  element in physical space at the point described by @a Trans. */
403  /** Each row of the result DenseMatrix @a shape contains the components of
404  one vector shape function. The size (#dof x SDim) of @a shape must be set
405  in advance, where SDim >= #dim is the physical space dimension as
406  described by @a Trans. */
407  virtual void CalcVShape(ElementTransformation &Trans,
408  DenseMatrix &shape) const;
409 
410  /// Equivalent to the CalcVShape() method with the same arguments.
412  { CalcVShape(Trans, shape); }
413 
414  /** @brief Evaluate the divergence of all shape functions of a *vector*
415  finite element in reference space at the given point @a ip. */
416  /** The size (#dof) of the result Vector @a divshape must be set in advance.
417  */
418  virtual void CalcDivShape(const IntegrationPoint &ip,
419  Vector &divshape) const;
420 
421  /** @brief Evaluate the divergence of all shape functions of a *vector*
422  finite element in physical space at the point described by @a Trans. */
423  /** The size (#dof) of the result Vector @a divshape must be set in advance.
424  */
425  void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const;
426 
427  /** @brief Evaluate the curl of all shape functions of a *vector* finite
428  element in reference space at the given point @a ip. */
429  /** Each row of the result DenseMatrix @a curl_shape contains the components
430  of the curl of one vector shape function. The size (#dof x CDim) of
431  @a curl_shape must be set in advance, where CDim = 3 for #dim = 3 and
432  CDim = 1 for #dim = 2. */
433  virtual void CalcCurlShape(const IntegrationPoint &ip,
434  DenseMatrix &curl_shape) const;
435 
436  /** @brief Evaluate the curl of all shape functions of a *vector* finite
437  element in physical space at the point described by @a Trans. */
438  /** Each row of the result DenseMatrix @a curl_shape contains the components
439  of the curl of one vector shape function. The size (#dof x CDim) of
440  @a curl_shape must be set in advance, where CDim = 3 for #dim = 3 and
441  CDim = 1 for #dim = 2. */
442  virtual void CalcPhysCurlShape(ElementTransformation &Trans,
443  DenseMatrix &curl_shape) const;
444 
445  /** @brief Get the dofs associated with the given @a face.
446  @a *dofs is set to an internal array of the local dofc on the
447  face, while *ndofs is set to the number of dofs on that face.
448  */
449  virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const;
450 
451  /** @brief Evaluate the Hessians of all shape functions of a scalar finite
452  element in reference space at the given point @a ip. */
453  /** Each row of the result DenseMatrix @a Hessian contains upper triangular
454  part of the Hessian of one shape function.
455  The order in 2D is {u_xx, u_xy, u_yy}.
456  The size (#dof x (#dim (#dim+1)/2) of @a Hessian must be set in advance.*/
457  virtual void CalcHessian(const IntegrationPoint &ip,
458  DenseMatrix &Hessian) const;
459 
460  /** @brief Evaluate the Hessian of all shape functions of a scalar finite
461  element in reference space at the given point @a ip. */
462  /** The size (#dof, #dim*(#dim+1)/2) of @a Hessian must be set in advance. */
463  virtual void CalcPhysHessian(ElementTransformation &Trans,
464  DenseMatrix& Hessian) const;
465 
466  /** @brief Evaluate the Laplacian of all shape functions of a scalar finite
467  element in reference space at the given point @a ip. */
468  /** The size (#dof) of @a Laplacian must be set in advance. */
469  virtual void CalcPhysLaplacian(ElementTransformation &Trans,
470  Vector& Laplacian) const;
471 
472  virtual void CalcPhysLinLaplacian(ElementTransformation &Trans,
473  Vector& Laplacian) const;
474 
475  /** @brief Return the local interpolation matrix @a I (Dof x Dof) where the
476  fine element is the image of the base geometry under the given
477  transformation. */
478  virtual void GetLocalInterpolation(ElementTransformation &Trans,
479  DenseMatrix &I) const;
480 
481  /** @brief Return a local restriction matrix @a R (Dof x Dof) mapping fine
482  dofs to coarse dofs.
483 
484  The fine element is the image of the base geometry under the given
485  transformation, @a Trans.
486 
487  The assumption in this method is that a subset of the coarse dofs can be
488  expressed only in terms of the dofs of the given fine element.
489 
490  Rows in @a R corresponding to coarse dofs that cannot be expressed in
491  terms of the fine dofs will be marked as invalid by setting the first
492  entry (column 0) in the row to infinity().
493 
494  This method assumes that the dimensions of @a R are set before it is
495  called. */
496  virtual void GetLocalRestriction(ElementTransformation &Trans,
497  DenseMatrix &R) const;
498 
499  /** @brief Return interpolation matrix, @a I, which maps dofs from a coarse
500  element, @a fe, to the fine dofs on @a this finite element. */
501  /** @a Trans represents the mapping from the reference element of @a this
502  element into a subset of the reference space of the element @a fe, thus
503  allowing the "coarse" FiniteElement to be different from the "fine"
504  FiniteElement as when h-refinement is combined with p-refinement or
505  p-derefinement. It is assumed that both finite elements use the same
506  FiniteElement::MapType. */
507  virtual void GetTransferMatrix(const FiniteElement &fe,
508  ElementTransformation &Trans,
509  DenseMatrix &I) const;
510 
511  /** @brief Given a coefficient and a transformation, compute its projection
512  (approximation) in the local finite dimensional space in terms
513  of the degrees of freedom. */
514  /** The approximation used to project is usually local interpolation of
515  degrees of freedom. The derived class could use other methods not
516  implemented yet, e.g. local L2 projection. */
517  virtual void Project(Coefficient &coeff,
518  ElementTransformation &Trans, Vector &dofs) const;
519 
520  /** @brief Given a vector coefficient and a transformation, compute its
521  projection (approximation) in the local finite dimensional space
522  in terms of the degrees of freedom. (VectorFiniteElements) */
523  /** The approximation used to project is usually local interpolation of
524  degrees of freedom. The derived class could use other methods not
525  implemented yet, e.g. local L2 projection. */
526  virtual void Project(VectorCoefficient &vc,
527  ElementTransformation &Trans, Vector &dofs) const;
528 
529  /** @brief Given a vector of values at the finite element nodes and a
530  transformation, compute its projection (approximation) in the local
531  finite dimensional space in terms of the degrees of freedom. Valid for
532  VectorFiniteElements. */
533  virtual void ProjectFromNodes(Vector &vc, ElementTransformation &Trans,
534  Vector &dofs) const;
535 
536  /** @brief Given a matrix coefficient and a transformation, compute an
537  approximation ("projection") in the local finite dimensional space in
538  terms of the degrees of freedom. For VectorFiniteElements, the rows of
539  the coefficient are projected in the vector space. */
540  virtual void ProjectMatrixCoefficient(
541  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const;
542 
543  /** @brief Project a delta function centered on the given @a vertex in
544  the local finite dimensional space represented by the @a dofs. */
545  virtual void ProjectDelta(int vertex, Vector &dofs) const;
546 
547  /** @brief Compute the embedding/projection matrix from the given
548  FiniteElement onto 'this' FiniteElement. The ElementTransformation is
549  included to support cases when the projection depends on it. */
550  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
551  DenseMatrix &I) const;
552 
553  /** @brief Compute the discrete gradient matrix from the given FiniteElement
554  onto 'this' FiniteElement. The ElementTransformation is included to
555  support cases when the matrix depends on it. */
556  virtual void ProjectGrad(const FiniteElement &fe,
557  ElementTransformation &Trans,
558  DenseMatrix &grad) const;
559 
560  /** @brief Compute the discrete curl matrix from the given FiniteElement onto
561  'this' FiniteElement. The ElementTransformation is included to support
562  cases when the matrix depends on it. */
563  virtual void ProjectCurl(const FiniteElement &fe,
564  ElementTransformation &Trans,
565  DenseMatrix &curl) const;
566 
567  /** @brief Compute the discrete divergence matrix from the given
568  FiniteElement onto 'this' FiniteElement. The ElementTransformation is
569  included to support cases when the matrix depends on it. */
570  virtual void ProjectDiv(const FiniteElement &fe,
571  ElementTransformation &Trans,
572  DenseMatrix &div) const;
573 
574  /** @brief Return a DofToQuad structure corresponding to the given
575  IntegrationRule using the given DofToQuad::Mode. */
576  /** See the documentation for DofToQuad for more details. */
577  virtual const DofToQuad &GetDofToQuad(const IntegrationRule &ir,
578  DofToQuad::Mode mode) const;
579 
580 
581  /** @brief Return the mapping from lexicographic face DOFs to lexicographic
582  element DOFs for the given local face @a face_id. */
583  /** Given the @a ith DOF (lexicographically ordered) on the face referenced
584  by @a face_id, face_map[i] gives the corresponding index of the DOF in
585  the element (also lexicographically ordered).
586 
587  @note For L2 spaces, this is only well-defined for "closed" bases such as
588  the Gauss-Lobatto or Bernstein (positive) bases.
589 
590  @warning GetFaceMap() is currently only implemented for tensor-product
591  (quadrilateral and hexahedral) elements. Its functionality may change
592  when simplex elements are supported in the future. */
593  virtual void GetFaceMap(const int face_id, Array<int> &face_map) const;
594 
595  /** @brief Return a DoF transformation object for this particular type of
596  basis.
597  */
599  { return NULL; }
600 
601  /// Deconstruct the FiniteElement
602  virtual ~FiniteElement();
603 
604  /** @brief Return true if the BasisType of @a b_type is closed
605  (has Quadrature1D points on the boundary). */
606  static bool IsClosedType(int b_type)
607  {
608  const int q_type = BasisType::GetQuadrature1D(b_type);
609  return ((q_type != Quadrature1D::Invalid) &&
611  }
612 
613  /** @brief Return true if the BasisType of @a b_type is open
614  (doesn't have Quadrature1D points on the boundary). */
615  static bool IsOpenType(int b_type)
616  {
617  const int q_type = BasisType::GetQuadrature1D(b_type);
618  return ((q_type != Quadrature1D::Invalid) &&
620  }
621 
622  /** @brief Ensure that the BasisType of @a b_type is closed
623  (has Quadrature1D points on the boundary). */
624  static int VerifyClosed(int b_type)
625  {
626  MFEM_VERIFY(IsClosedType(b_type),
627  "invalid closed basis type: " << b_type);
628  return b_type;
629  }
630 
631  /** @brief Ensure that the BasisType of @a b_type is open
632  (doesn't have Quadrature1D points on the boundary). */
633  static int VerifyOpen(int b_type)
634  {
635  MFEM_VERIFY(IsOpenType(b_type), "invalid open basis type: " << b_type);
636  return b_type;
637  }
638 
639  /** @brief Ensure that the BasisType of @a b_type nodal
640  (satisfies the interpolation property). */
641  static int VerifyNodal(int b_type)
642  {
643  return BasisType::CheckNodal(b_type);
644  }
645 };
646 
647 /** @brief Class for finite elements with basis functions
648  that return scalar values. */
650 {
651 protected:
653  {
654  MFEM_VERIFY(fe.GetRangeType() == SCALAR,
655  "'fe' must be a ScalarFiniteElement");
656  return static_cast<const ScalarFiniteElement &>(fe);
657  }
658 
659 public:
660  /** @brief Construct ScalarFiniteElement with given
661  @param D Reference space dimension
662  @param G Geometry type (of type Geometry::Type)
663  @param Do Number of degrees of freedom in the FiniteElement
664  @param O Order/degree of the FiniteElement
665  @param F FunctionSpace type of the FiniteElement
666  */
667  ScalarFiniteElement(int D, Geometry::Type G, int Do, int O,
668  int F = FunctionSpace::Pk)
669  : FiniteElement(D, G, Do, O, F)
671 
672  /** @brief Set the FiniteElement::MapType of the element to either VALUE or
673  INTEGRAL. Also sets the FiniteElement::DerivType to GRAD if the
674  FiniteElement::MapType is VALUE. */
675  virtual void SetMapType(int M)
676  {
677  MFEM_VERIFY(M == VALUE || M == INTEGRAL, "unknown MapType");
678  map_type = M;
679  deriv_type = (M == VALUE) ? GRAD : NONE;
680  }
681 
682  /** @brief Get the matrix @a I that defines nodal interpolation
683  @a between this element and the refined element @a fine_fe. */
685  DenseMatrix &I,
686  const ScalarFiniteElement &fine_fe) const;
687 
688  /** @brief Get matrix @a I "Interpolation" defined through local
689  L2-projection in the space defined by the @a fine_fe. */
690  /** If the "fine" elements cannot represent all basis functions of the
691  "coarse" element, then boundary values from different sub-elements are
692  generally different. */
694  DenseMatrix &I,
695  const ScalarFiniteElement &fine_fe) const;
696 
697  /** @brief Get restriction matrix @a R defined through local L2-projection
698  in the space defined by the @a coarse_fe. */
699  /** If the "fine" elements cannot represent all basis functions of the
700  "coarse" element, then boundary values from different sub-elements are
701  generally different. */
703  DenseMatrix &R,
704  const ScalarFiniteElement &coarse_fe) const;
705 };
706 
707 /// Class for standard nodal finite elements.
709 {
710 protected:
712  void ProjectCurl_2D(const FiniteElement &fe,
713  ElementTransformation &Trans,
714  DenseMatrix &curl) const;
715 
716 public:
717  /** @brief Construct NodalFiniteElement with given
718  @param D Reference space dimension
719  @param G Geometry type (of type Geometry::Type)
720  @param Do Number of degrees of freedom in the FiniteElement
721  @param O Order/degree of the FiniteElement
722  @param F FunctionSpace type of the FiniteElement
723  */
724  NodalFiniteElement(int D, Geometry::Type G, int Do, int O,
725  int F = FunctionSpace::Pk)
726  : ScalarFiniteElement(D, G, Do, O, F) { }
727 
729  DenseMatrix &I) const override
730  { NodalLocalInterpolation(Trans, I, *this); }
731 
733  DenseMatrix &R) const override;
734 
736  ElementTransformation &Trans,
737  DenseMatrix &I) const override
738  { CheckScalarFE(fe).NodalLocalInterpolation(Trans, I, *this); }
739 
740  void Project(Coefficient &coeff,
741  ElementTransformation &Trans, Vector &dofs) const override;
742 
743  void Project(VectorCoefficient &vc,
744  ElementTransformation &Trans, Vector &dofs) const override;
745 
746  // (mc.height x mc.width) @ DOFs -> (Dof x mc.width x mc.height) in dofs
748  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const override;
749 
750  void Project(const FiniteElement &fe, ElementTransformation &Trans,
751  DenseMatrix &I) const override;
752 
753  void ProjectGrad(const FiniteElement &fe,
754  ElementTransformation &Trans,
755  DenseMatrix &grad) const override;
756 
757  void ProjectDiv(const FiniteElement &fe,
758  ElementTransformation &Trans,
759  DenseMatrix &div) const override;
760 
761  /** @brief Get an Array<int> that maps lexicographically ordered indices to
762  the indices of the respective nodes/dofs/basis functions.
763 
764  Lexicographic ordering of nodes is defined in terms of reference-space
765  coordinates (x,y,z). Lexicographically ordered nodes are listed first in
766  order of increasing x-coordinate, and then in order of increasing
767  y-coordinate, and finally in order of increasing z-coordinate.
768 
769  For example, the six nodes of a quadratic triangle are lexicographically
770  ordered as follows:
771 
772  5
773  |\
774  3 4
775  | \
776  0-1-2
777 
778  The resulting array may be empty if the DOFs are already ordered
779  lexicographically, or if the finite element does not support creating
780  this permutation. The array returned is the same as the array given by
781  TensorBasisElement::GetDofMap, but it is also available for non-tensor
782  elements. */
784 };
785 
786 /** @brief Intermediate class for finite elements whose basis functions return
787  vector values. */
789 {
790  // Hide the scalar functions CalcShape and CalcDShape.
791 private:
792  /// Overrides the scalar CalcShape function to print an error.
793  void CalcShape(const IntegrationPoint &ip,
794  Vector &shape) const override;
795 
796  /// Overrides the scalar CalcDShape function to print an error.
797  void CalcDShape(const IntegrationPoint &ip,
798  DenseMatrix &dshape) const override;
799 
800 protected:
801  bool is_nodal;
802 #ifndef MFEM_THREAD_SAFE
803  mutable DenseMatrix JtJ;
805 #endif
806  void SetDerivMembers();
807 
809  DenseMatrix &shape) const;
810 
812  DenseMatrix &shape) const;
813 
814  /** @brief Project a vector coefficient onto the RT basis functions
815  @param nk Face normal vectors for this element type
816  @param d2n Offset into nk for each degree of freedom
817  @param vc Vector coefficient to be projected
818  @param Trans Transformation from reference to physical coordinates
819  @param dofs Expansion coefficients for the approximation of vc
820  */
821  void Project_RT(const double *nk, const Array<int> &d2n,
823  Vector &dofs) const;
824 
825  /// Projects the vector of values given at FE nodes to RT space
826  /** Project vector values onto the RT basis functions
827  @param nk Face normal vectors for this element type
828  @param d2n Offset into nk for each degree of freedom
829  @param vc Vector values at each interpolation point
830  @param Trans Transformation from reference to physical coordinates
831  @param dofs Expansion coefficients for the approximation of vc
832  */
833  void Project_RT(const double *nk, const Array<int> &d2n,
834  Vector &vc, ElementTransformation &Trans,
835  Vector &dofs) const;
836 
837  /// Project the rows of the matrix coefficient in an RT space
839  const double *nk, const Array<int> &d2n,
840  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const;
841 
842  /** @brief Project vector-valued basis functions onto the RT basis functions
843  @param nk Face normal vectors for this element type
844  @param d2n Offset into nk for each degree of freedom
845  @param fe Vector-valued finite element basis
846  @param Trans Transformation from reference to physical coordinates
847  @param I Expansion coefficients for the approximation of each basis
848  function
849 
850  Note: If the FiniteElement, fe, is scalar-valued the projection will
851  assume that a FiniteElementSpace is being used to define a vector
852  field using the scalar basis functions for each component of the
853  vector field.
854  */
855  void Project_RT(const double *nk, const Array<int> &d2n,
856  const FiniteElement &fe, ElementTransformation &Trans,
857  DenseMatrix &I) const;
858 
859  // rotated gradient in 2D
860  void ProjectGrad_RT(const double *nk, const Array<int> &d2n,
861  const FiniteElement &fe, ElementTransformation &Trans,
862  DenseMatrix &grad) const;
863 
864  // Compute the curl as a discrete operator from ND FE (fe) to ND FE (this).
865  // The natural FE for the range is RT, so this is an approximation.
866  void ProjectCurl_ND(const double *tk, const Array<int> &d2t,
867  const FiniteElement &fe, ElementTransformation &Trans,
868  DenseMatrix &curl) const;
869 
870  void ProjectCurl_RT(const double *nk, const Array<int> &d2n,
871  const FiniteElement &fe, ElementTransformation &Trans,
872  DenseMatrix &curl) const;
873 
874  /** @brief Project a vector coefficient onto the ND basis functions
875  @param tk Edge tangent vectors for this element type
876  @param d2t Offset into tk for each degree of freedom
877  @param vc Vector coefficient to be projected
878  @param Trans Transformation from reference to physical coordinates
879  @param dofs Expansion coefficients for the approximation of vc
880  */
881  void Project_ND(const double *tk, const Array<int> &d2t,
883  Vector &dofs) const;
884 
885  /// Projects the vector of values given at FE nodes to ND space
886  /** Project vector values onto the ND basis functions
887  @param tk Edge tangent vectors for this element type
888  @param d2t Offset into tk for each degree of freedom
889  @param vc Vector values at each interpolation point
890  @param Trans Transformation from reference to physical coordinates
891  @param dofs Expansion coefficients for the approximation of vc
892  */
893  void Project_ND(const double *tk, const Array<int> &d2t,
894  Vector &vc, ElementTransformation &Trans,
895  Vector &dofs) const;
896 
897  /// Project the rows of the matrix coefficient in an ND space
899  const double *tk, const Array<int> &d2t,
900  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const;
901 
902  /** @brief Project vector-valued basis functions onto the ND basis functions
903  @param tk Edge tangent vectors for this element type
904  @param d2t Offset into tk for each degree of freedom
905  @param fe Vector-valued finite element basis
906  @param Trans Transformation from reference to physical coordinates
907  @param I Expansion coefficients for the approximation of each basis
908  function
909 
910  Note: If the FiniteElement, fe, is scalar-valued the projection will
911  assume that a FiniteElementSpace is being used to define a vector
912  field using the scalar basis functions for each component of the
913  vector field.
914  */
915  void Project_ND(const double *tk, const Array<int> &d2t,
916  const FiniteElement &fe, ElementTransformation &Trans,
917  DenseMatrix &I) const;
918 
919  void ProjectGrad_ND(const double *tk, const Array<int> &d2t,
920  const FiniteElement &fe, ElementTransformation &Trans,
921  DenseMatrix &grad) const;
922 
924  ElementTransformation &Trans,
925  DenseMatrix &I) const;
926 
928  const double *nk, const Array<int> &d2n,
929  ElementTransformation &Trans,
930  DenseMatrix &I) const;
931 
933  ElementTransformation &Trans,
934  DenseMatrix &I) const;
935 
937  const double *tk, const Array<int> &d2t,
938  ElementTransformation &Trans,
939  DenseMatrix &I) const;
940 
941  void LocalRestriction_RT(const double *nk, const Array<int> &d2n,
942  ElementTransformation &Trans,
943  DenseMatrix &R) const;
944 
945  void LocalRestriction_ND(const double *tk, const Array<int> &d2t,
946  ElementTransformation &Trans,
947  DenseMatrix &R) const;
948 
950  {
951  if (fe.GetRangeType() != VECTOR)
952  { mfem_error("'fe' must be a VectorFiniteElement"); }
953  return static_cast<const VectorFiniteElement &>(fe);
954  }
955 
956 public:
957  VectorFiniteElement(int D, Geometry::Type G, int Do, int O, int M,
958  int F = FunctionSpace::Pk);
959 };
960 
961 /// @brief Class for computing 1D special polynomials and their associated basis
962 /// functions
963 class Poly_1D
964 {
965 public:
966  /// One-dimensional basis evaluation type
967  enum EvalType
968  {
969  ChangeOfBasis = 0, ///< Use change of basis, O(p^2) Evals
970  Barycentric = 1, ///< Use barycentric Lagrangian interpolation, O(p) Evals
971  Positive = 2, ///< Fast evaluation of Bernstein polynomials
972  Integrated = 3, ///< Integrated indicator functions (cf. Gerritsma)
973  NumEvalTypes = 4 ///< Keep count of the number of eval types
974  };
975 
976  /// @brief Class for evaluating 1D nodal, positive (Bernstein), or integrated
977  /// (Gerritsma) bases.
978  class Basis
979  {
980  private:
981  EvalType etype; ///< Determines how the basis functions should be evaluated.
983  mutable Vector x, w;
984  /// The following data members are used for "integrated basis type", which
985  /// is defined in terms of nodal basis of one degree higher.
986  ///@{
987  mutable Vector u_aux, d_aux, d2_aux;
988  ///@}
989  /// @brief An auxiliary nodal basis used to evaluate the integrated basis.
990  /// This member variable is NULL whenever etype != Integrated.
991  Basis *auxiliary_basis;
992  /// Should the integrated basis functions be scaled? See ScaleIntegrated.
993  bool scale_integrated;
994 
995  public:
996  /// Create a nodal or positive (Bernstein) basis of degree @a p
997  Basis(const int p, const double *nodes, EvalType etype = Barycentric);
998  /// Evaluate the basis functions at point @a x in [0,1]
999  void Eval(const double x, Vector &u) const;
1000  /// @brief Evaluate the basis functions and their derivatives at point @a
1001  /// x in [0,1]
1002  void Eval(const double x, Vector &u, Vector &d) const;
1003  /// @brief Evaluate the basis functions and their first two derivatives at
1004  /// point @a x in [0,1]
1005  void Eval(const double x, Vector &u, Vector &d, Vector &d2) const;
1006  /// @brief Evaluate the "integrated" basis type using pre-computed closed
1007  /// basis derivatives.
1008  ///
1009  /// This basis is given by the negative partial sum of the corresponding
1010  /// closed basis derivatives. The closed basis derivatives are given by @a
1011  /// d, and the result is stored in @a i.
1012  void EvalIntegrated(const Vector &d, Vector &i) const;
1013  /// @brief Set whether the "integrated" basis should be scaled by the
1014  /// subcell sizes. Has no effect for non-integrated bases.
1015  ///
1016  /// Generally, this should be true for mfem::FiniteElement::MapType VALUE
1017  /// and false for all other map types. If this option is enabled, the
1018  /// basis functions will be scaled by the widths of the subintervals, so
1019  /// that the basis functions represent mean values. Otherwise, the basis
1020  /// functions represent integrated values.
1021  void ScaleIntegrated(bool scale_integrated_);
1022  /// Returns true if the basis is "integrated", false otherwise.
1023  bool IsIntegratedType() const { return etype == Integrated; }
1024  ~Basis();
1025  };
1026 
1027 private:
1028  typedef std::map< int, Array<double*>* > PointsMap;
1029  typedef std::map< int, Array<Basis*>* > BasisMap;
1030 
1031  MemoryType h_mt;
1032  PointsMap points_container;
1033  BasisMap bases_container;
1034 
1035  static Array2D<int> binom;
1036 
1037  static void CalcMono(const int p, const double x, double *u);
1038  static void CalcMono(const int p, const double x, double *u, double *d);
1039 
1040  static void CalcChebyshev(const int p, const double x, double *u);
1041  static void CalcChebyshev(const int p, const double x, double *u, double *d);
1042  static void CalcChebyshev(const int p, const double x, double *u, double *d,
1043  double *dd);
1044 
1045  QuadratureFunctions1D quad_func;
1046 
1047 public:
1048  Poly_1D(): h_mt(MemoryType::HOST) { }
1049 
1050  /** @brief Get a pointer to an array containing the binomial coefficients "p
1051  choose k" for k=0,...,p for the given p. */
1052  static const int *Binom(const int p);
1053 
1054  /** @brief Get the coordinates of the points of the given BasisType,
1055  @a btype.
1056 
1057  @param[in] p The polynomial degree; the number of points is `p+1`.
1058  @param[in] btype The BasisType.
1059 
1060  @return A pointer to an array containing the `p+1` coordinates of the
1061  points. Returns NULL if the BasisType has no associated set of
1062  points. */
1063  const double *GetPoints(const int p, const int btype);
1064 
1065  /// Get coordinates of an open (GaussLegendre) set of points if degree @a p
1066  const double *OpenPoints(const int p,
1067  const int btype = BasisType::GaussLegendre)
1068  { return GetPoints(p, btype); }
1069 
1070  /// Get coordinates of a closed (GaussLegendre) set of points if degree @a p
1071  const double *ClosedPoints(const int p,
1072  const int btype = BasisType::GaussLobatto)
1073  { return GetPoints(p, btype); }
1074 
1075  /** @brief Get a Poly_1D::Basis object of the given degree and BasisType,
1076  @a btype.
1077 
1078  @param[in] p The polynomial degree of the basis.
1079  @param[in] btype The BasisType.
1080 
1081  @return A reference to an object of type Poly_1D::Basis that represents
1082  the requested basis type. */
1083  Basis &GetBasis(const int p, const int btype);
1084 
1085  /** @brief Evaluate the values of a hierarchical 1D basis at point x
1086  hierarchical = k-th basis function is degree k polynomial */
1087  static void CalcBasis(const int p, const double x, double *u)
1088  // { CalcMono(p, x, u); }
1089  // Bernstein basis is not hierarchical --> does not work for triangles
1090  // and tetrahedra
1091  // { CalcBernstein(p, x, u); }
1092  // { CalcLegendre(p, x, u); }
1093  { CalcChebyshev(p, x, u); }
1094 
1095  /** @brief Evaluate the values of a hierarchical 1D basis at point x
1096  hierarchical = k-th basis function is degree k polynomial */
1097  static void CalcBasis(const int p, const double x, Vector &u)
1098  { CalcBasis(p, x, u.GetData()); }
1099 
1100  /// Evaluate the values and derivatives of a hierarchical 1D basis at point @a x
1101  static void CalcBasis(const int p, const double x, double *u, double *d)
1102  // { CalcMono(p, x, u, d); }
1103  // { CalcBernstein(p, x, u, d); }
1104  // { CalcLegendre(p, x, u, d); }
1105  { CalcChebyshev(p, x, u, d); }
1106 
1107  /** @brief Evaluate the values and derivatives of a hierarchical 1D basis at
1108  point @a x. */
1109  static void CalcBasis(const int p, const double x, Vector &u, Vector &d)
1110  { CalcBasis(p, x, u.GetData(), d.GetData()); }
1111 
1112  /// Evaluate the values, derivatives and second derivatives of a hierarchical 1D basis at point x
1113  static void CalcBasis(const int p, const double x, double *u, double *d,
1114  double *dd)
1115  // { CalcMono(p, x, u, d); }
1116  // { CalcBernstein(p, x, u, d); }
1117  // { CalcLegendre(p, x, u, d); }
1118  { CalcChebyshev(p, x, u, d, dd); }
1119 
1120  /** @brief Evaluate the values, derivatives and second derivatives of a
1121  hierarchical 1D basis at point @a x. */
1122  static void CalcBasis(const int p, const double x, Vector &u, Vector &d,
1123  Vector &dd)
1124  { CalcBasis(p, x, u.GetData(), d.GetData(), dd.GetData()); }
1125 
1126  /// Evaluate a representation of a Delta function at point x
1127  static double CalcDelta(const int p, const double x)
1128  { return pow(x, (double) p); }
1129 
1130  /** @brief Compute the points for the Chebyshev polynomials of order @a p
1131  and place them in the already allocated @a x array. */
1132  static void ChebyshevPoints(const int p, double *x);
1133 
1134  /** @brief Compute the @a p terms in the expansion of the binomial (x + y)^p
1135  and store them in the already allocated @a u array. */
1136  static void CalcBinomTerms(const int p, const double x, const double y,
1137  double *u);
1138  /** @brief Compute the terms in the expansion of the binomial (x + y)^p and
1139  their derivatives with respect to x assuming that dy/dx = -1. Store the
1140  results in the already allocated @a u and @a d arrays.*/
1141  static void CalcBinomTerms(const int p, const double x, const double y,
1142  double *u, double *d);
1143  /** @brief Compute the derivatives (w.r.t. x) of the terms in the expansion
1144  of the binomial (x + y)^p assuming that dy/dx = -1. Store the results
1145  in the already allocated @a d array.*/
1146  static void CalcDBinomTerms(const int p, const double x, const double y,
1147  double *d);
1148 
1149  /** @brief Compute the values of the Bernstein basis functions of order
1150  @a p at coordinate @a x and store the results in the already allocated
1151  @a u array. */
1152  static void CalcBernstein(const int p, const double x, double *u)
1153  { CalcBinomTerms(p, x, 1. - x, u); }
1154 
1155  /** @brief Compute the values of the Bernstein basis functions of order
1156  @a p at coordinate @a x and store the results in the already allocated
1157  @a u array. */
1158  static void CalcBernstein(const int p, const double x, Vector &u)
1159  { CalcBernstein(p, x, u.GetData()); }
1160 
1161  /** @brief Compute the values and derivatives of the Bernstein basis functions
1162  of order @a p at coordinate @a x and store the results in the already allocated
1163  @a u and @a d arrays. */
1164  static void CalcBernstein(const int p, const double x, double *u, double *d)
1165  { CalcBinomTerms(p, x, 1. - x, u, d); }
1166 
1167  /** @brief Compute the values and derivatives of the Bernstein basis
1168  functions of order @a p at coordinate @a x and store the results in the
1169  already allocated @a u and @a d arrays. */
1170  static void CalcBernstein(const int p, const double x, Vector &u, Vector &d)
1171  { CalcBernstein(p, x, u.GetData(), d.GetData()); }
1172 
1173  static void CalcLegendre(const int p, const double x, double *u);
1174  static void CalcLegendre(const int p, const double x, double *u, double *d);
1175 
1176  ~Poly_1D();
1177 };
1178 
1179 extern MFEM_EXPORT Poly_1D poly1d;
1180 
1181 /// An element defined as an ND tensor product of 1D elements on a segment,
1182 /// square, or cube
1184 {
1185 protected:
1186  int b_type;
1190 
1191 public:
1193  {
1196  Sr_DOF_MAP = 2, // Sr = Serendipity
1197  };
1198 
1199  TensorBasisElement(const int dims, const int p, const int btype,
1200  const DofMapType dmtype);
1201 
1202  int GetBasisType() const { return b_type; }
1203 
1204  const Poly_1D::Basis &GetBasis1D() const { return basis1d; }
1205 
1206  /** @brief Get an Array<int> that maps lexicographically ordered indices to
1207  the indices of the respective nodes/dofs/basis functions. If the dofs are
1208  ordered lexicographically, i.e. the mapping is identity, the returned
1209  Array will be empty. */
1210  const Array<int> &GetDofMap() const { return dof_map; }
1211 
1213  {
1214  switch (dim)
1215  {
1216  case 1: return Geometry::SEGMENT;
1217  case 2: return Geometry::SQUARE;
1218  case 3: return Geometry::CUBE;
1219  default:
1220  MFEM_ABORT("invalid dimension: " << dim);
1221  return Geometry::INVALID;
1222  }
1223  }
1224 
1225  /// Return @a base raised to the power @a dim.
1226  static int Pow(int base, int dim)
1227  {
1228  switch (dim)
1229  {
1230  case 1: return base;
1231  case 2: return base*base;
1232  case 3: return base*base*base;
1233  default: MFEM_ABORT("invalid dimension: " << dim); return -1;
1234  }
1235  }
1236 
1237  static const DofToQuad &GetTensorDofToQuad(
1238  const FiniteElement &fe, const IntegrationRule &ir,
1239  DofToQuad::Mode mode, const Poly_1D::Basis &basis, bool closed,
1240  Array<DofToQuad*> &dof2quad_array);
1241 };
1242 
1244  public TensorBasisElement
1245 {
1246 public:
1247  NodalTensorFiniteElement(const int dims, const int p, const int btype,
1248  const DofMapType dmtype);
1249 
1251  DofToQuad::Mode mode) const override
1252  {
1253  return (mode == DofToQuad::FULL) ?
1254  FiniteElement::GetDofToQuad(ir, mode) :
1255  GetTensorDofToQuad(*this, ir, mode, basis1d, true, dof2quad_array);
1256  }
1257 
1258  void SetMapType(const int map_type_) override;
1259 
1261  ElementTransformation &Trans,
1262  DenseMatrix &I) const override
1263  {
1264  if (basis1d.IsIntegratedType())
1265  {
1266  CheckScalarFE(fe).ScalarLocalInterpolation(Trans, I, *this);
1267  }
1268  else
1269  {
1271  }
1272  }
1273 
1274  void GetFaceMap(const int face_id, Array<int> &face_map) const override;
1275 };
1276 
1278  public TensorBasisElement
1279 {
1280 private:
1281  mutable Array<DofToQuad*> dof2quad_array_open;
1282 
1283 protected:
1285 
1286 public:
1287  VectorTensorFiniteElement(const int dims, const int d, const int p,
1288  const int cbtype, const int obtype,
1289  const int M, const DofMapType dmtype);
1290 
1291  // For 1D elements: there is only an "open basis", no "closed basis"
1292  VectorTensorFiniteElement(const int dims, const int d, const int p,
1293  const int obtype, const int M,
1294  const DofMapType dmtype);
1295 
1297  DofToQuad::Mode mode) const override
1298  {
1299  return (mode == DofToQuad::FULL) ?
1300  FiniteElement::GetDofToQuad(ir, mode) :
1301  GetTensorDofToQuad(*this, ir, mode, basis1d, true, dof2quad_array);
1302  }
1303 
1305  DofToQuad::Mode mode) const
1306  {
1307  MFEM_VERIFY(mode != DofToQuad::FULL, "invalid mode requested");
1308  return GetTensorDofToQuad(*this, ir, mode, obasis1d, false,
1309  dof2quad_array_open);
1310  }
1311 
1312  virtual ~VectorTensorFiniteElement();
1313 };
1314 
1315 void InvertLinearTrans(ElementTransformation &trans,
1316  const IntegrationPoint &pt, Vector &x);
1317 
1318 } // namespace mfem
1319 
1320 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:233
virtual void CalcPhysLinLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Definition: fe_base.cpp:244
Tensor products of polynomials of order k.
Definition: fe_base.hpp:222
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:412
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
Definition: fe_base.hpp:1023
ScalarFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct ScalarFiniteElement with given.
Definition: fe_base.hpp:667
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...
Definition: fe_base.hpp:1164
const class FiniteElement * FE
The FiniteElement that created and owns this object.
Definition: fe_base.hpp:141
Array< int > lex_ordering
Definition: fe_base.hpp:711
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
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...
Definition: fe_base.hpp:1170
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
Definition: fe_base.hpp:45
static double CalcDelta(const int p, const double x)
Evaluate a representation of a Delta function at point x.
Definition: fe_base.hpp:1127
EvalType
One-dimensional basis evaluation type.
Definition: fe_base.hpp:967
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
Definition: fe_base.cpp:1932
NodalTensorFiniteElement(const int dims, const int p, const int btype, const DofMapType dmtype)
Definition: fe_base.cpp:2500
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...
Definition: fe_base.cpp:138
Basis(const int p, const double *nodes, EvalType etype=Barycentric)
Create a nodal or positive (Bernstein) basis of degree p.
Definition: fe_base.cpp:1611
void LocalL2Projection_RT(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1353
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:161
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.
Definition: fe_base.hpp:78
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...
Definition: fe_base.hpp:1296
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...
Definition: fe_base.cpp:65
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.
Definition: fe_base.cpp:1017
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe_base.cpp:1120
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const override
Compute the discrete gradient matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:830
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...
Definition: fe_base.cpp:101
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...
Definition: fe_base.cpp:107
FiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct FiniteElement with given.
Definition: fe_base.cpp:23
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 ...
Definition: fe_base.hpp:1158
aka "open half" Newton-Cotes
Definition: intrules.hpp:400
int GetCurlDim() const
Returns the dimension of the curl for vector-valued finite elements.
Definition: fe_base.hpp:317
int dim
Dimension of reference space.
Definition: fe_base.hpp:236
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...
Definition: fe_base.hpp:1250
TensorBasisElement(const int dims, const int p, const int btype, const DofMapType dmtype)
Definition: fe_base.cpp:2259
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...
Definition: fe_base.cpp:203
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:173
Array< double > Gt
Transpose of G.
Definition: fe_base.hpp:212
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:169
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
Definition: fe_base.cpp:1675
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
Definition: fe_base.hpp:61
VectorFiniteElement(int D, Geometry::Type G, int Do, int O, int M, int F=FunctionSpace::Pk)
Definition: fe_base.cpp:891
bool HasAnisotropicOrders() const
Returns true if the FiniteElement basis may be using different orders/degrees in different spatial di...
Definition: fe_base.hpp:331
MapType
Enumeration for MapType: defines how reference functions are mapped to physical space.
Definition: fe_base.hpp:269
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
Definition: fe_base.hpp:411
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...
Definition: fe_base.hpp:1122
DenseMatrix vshape
Definition: fe_base.hpp:248
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe_base.cpp:957
static const DofToQuad & GetTensorDofToQuad(const FiniteElement &fe, const IntegrationRule &ir, DofToQuad::Mode mode, const Poly_1D::Basis &basis, bool closed, Array< DofToQuad *> &dof2quad_array)
Definition: fe_base.cpp:2459
Geometry::Type geom_type
Geometry::Type of the reference element.
Definition: fe_base.hpp:239
DerivType
Enumeration for DerivType: defines which derivative method is implemented.
Definition: fe_base.hpp:292
Intermediate class for finite elements whose basis functions return vector values.
Definition: fe_base.hpp:788
static void CalcLegendre(const int p, const double x, double *u)
Definition: fe_base.cpp:2085
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe_base.hpp:92
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe_base.cpp:1147
Use change of basis, O(p^2) Evals.
Definition: fe_base.hpp:969
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...
Definition: fe_base.hpp:728
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Definition: fe_base.hpp:111
const double * ClosedPoints(const int p, const int btype=BasisType::GaussLobatto)
Get coordinates of a closed (GaussLegendre) set of points if degree p.
Definition: fe_base.hpp:1071
static int CheckNodal(int b_type)
If the input does not represents a valid nodal BasisType, abort with an error; otherwise return the i...
Definition: fe_base.hpp:53
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe_base.hpp:25
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.
Definition: fe_base.cpp:1233
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe_base.hpp:1210
int cdim
Dimension of curl for vector-valued basis functions.
Definition: fe_base.hpp:238
int vdim
Vector dimension of vector-valued basis functions.
Definition: fe_base.hpp:237
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...
Definition: fe_base.hpp:1087
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...
Definition: fe_base.cpp:1991
const double * GetPoints(const int p, const int btype)
Get the coordinates of the points of the given BasisType, btype.
Definition: fe_base.cpp:2184
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition: fe_base.cpp:126
Polynomials of order k.
Definition: fe_base.hpp:221
Class for standard nodal finite elements.
Definition: fe_base.hpp:708
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
VectorTensorFiniteElement(const int dims, const int d, const int p, const int cbtype, const int obtype, const int M, const DofMapType dmtype)
Definition: fe_base.cpp:2530
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...
Definition: fe_base.cpp:365
Class for finite elements with basis functions that return scalar values.
Definition: fe_base.hpp:649
void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const override
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
Definition: fe_base.cpp:749
static const ScalarFiniteElement & CheckScalarFE(const FiniteElement &fe)
Definition: fe_base.hpp:652
Nodes: x_i = (i+1/2)/n, i=0,...,n-1.
Definition: fe_base.hpp:36
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.hpp:675
void LocalRestriction_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
Definition: fe_base.cpp:1569
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...
Definition: fe_base.cpp:96
const int * GetAnisotropicOrders() const
Returns an array containing the anisotropic orders/degrees.
Definition: fe_base.hpp:334
void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const override
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe_base.cpp:678
Poly_1D::Basis & basis1d
Definition: fe_base.hpp:1188
Fast evaluation of Bernstein polynomials.
Definition: fe_base.hpp:971
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
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...
Definition: fe_base.cpp:1959
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition: fe_base.cpp:150
Class for computing 1D special polynomials and their associated basis functions.
Definition: fe_base.hpp:963
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
Definition: fe_base.hpp:165
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
Definition: fe_base.cpp:638
Bernstein polynomials.
Definition: fe_base.hpp:33
const Poly_1D::Basis & GetBasis1D() const
Definition: fe_base.hpp:1204
void LocalL2Projection_ND(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1441
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...
Definition: fe_base.cpp:71
RangeType
Enumeration for range_type and deriv_range_type.
Definition: fe_base.hpp:257
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...
Definition: fe_base.hpp:1097
static int VerifyOpen(int b_type)
Ensure that the BasisType of b_type is open (doesn&#39;t have Quadrature1D points on the boundary)...
Definition: fe_base.hpp:633
int orders[Geometry::MaxDim]
Anisotropic orders.
Definition: fe_base.hpp:245
Closed GaussLegendre.
Definition: fe_base.hpp:38
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:161
virtual ~FiniteElement()
Deconstruct the FiniteElement.
Definition: fe_base.cpp:500
Use barycentric Lagrangian interpolation, O(p) Evals.
Definition: fe_base.hpp:970
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...
Definition: fe_base.cpp:119
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...
Definition: fe_base.cpp:58
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
Definition: fe_base.cpp:2208
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:168
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...
Definition: fe_base.cpp:192
const double * OpenPoints(const int p, const int btype=BasisType::GaussLegendre)
Get coordinates of an open (GaussLegendre) set of points if degree p.
Definition: fe_base.hpp:1066
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:175
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...
Definition: fe_base.hpp:1113
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
Definition: fe_base.cpp:1907
const Array< int > & GetLexicographicOrdering() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe_base.hpp:783
static const int MaxDim
Definition: geom.hpp:43
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...
Definition: fe_base.hpp:1260
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
IntegrationRule Nodes
Definition: fe_base.hpp:246
void LocalInterpolation_RT(const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1400
Array< double > Bt
Transpose of B.
Definition: fe_base.hpp:190
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
Definition: fe_base.cpp:144
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
Definition: intrules.cpp:911
Integrated indicator functions (cf. Gerritsma)
Definition: fe_base.hpp:972
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
Refined tensor products of polynomials of order k.
Definition: fe_base.hpp:223
int GetDerivType() const
Returns the FiniteElement::DerivType of the element describing the spatial derivative method implemen...
Definition: fe_base.hpp:354
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
Definition: fe_base.hpp:1304
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.
Definition: fe_base.hpp:1109
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
static bool IsClosedType(int b_type)
Return true if the BasisType of b_type is closed (has Quadrature1D points on the boundary).
Definition: fe_base.hpp:606
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.
Definition: fe_base.cpp:980
void LocalInterpolation_ND(const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1487
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...
Definition: fe_base.cpp:548
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.
Definition: fe_base.cpp:1942
Poly_1D poly1d
Definition: fe.cpp:28
Serendipity basis (squares / cubes)
Definition: fe_base.hpp:37
Array< int > inv_dof_map
Definition: fe_base.hpp:1189
Nodes: x_i = (i+1)/(n+1), i=0,...,n-1.
Definition: fe_base.hpp:34
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Definition: intrules.cpp:923
void SetMapType(const int map_type_) override
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.cpp:2511
Integrated GLL indicator functions.
Definition: fe_base.hpp:39
int GetDerivMapType() const
Returns the FiniteElement::DerivType of the element describing how reference function derivatives are...
Definition: fe_base.hpp:359
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:184
aka closed Gauss Legendre
Definition: intrules.hpp:401
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
Definition: fe_base.hpp:104
aka closed Newton-Cotes
Definition: intrules.hpp:399
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 ...
Definition: fe_base.cpp:182
Class for integration point with weight.
Definition: intrules.hpp:31
Host memory; using new[] and delete[].
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition: fe_base.hpp:149
Full multidimensional representation which does not use tensor product structure. The ordering of the...
Definition: fe_base.hpp:154
virtual StatelessDofTransformation * GetDofTransformation() const
Return a DoF transformation object for this particular type of basis.
Definition: fe_base.hpp:598
void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const override
Compute the discrete divergence matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:859
static bool IsOpenType(int b_type)
Return true if the BasisType of b_type is open (doesn&#39;t have Quadrature1D points on the boundary)...
Definition: fe_base.hpp:615
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe_base.cpp:1185
int dof
Number of degrees of freedom.
Definition: fe_base.hpp:243
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
Array< DofToQuad * > dof2quad_array
Container for all DofToQuad objects created by the FiniteElement.
Definition: fe_base.hpp:253
static int VerifyClosed(int b_type)
Ensure that the BasisType of b_type is closed (has Quadrature1D points on the boundary).
Definition: fe_base.hpp:624
Implements CalcDivShape methods.
Definition: fe_base.hpp:296
int dim
Definition: ex24.cpp:53
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...
Definition: fe_base.hpp:735
aka open Newton-Cotes
Definition: intrules.hpp:398
NodalFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct NodalFiniteElement with given.
Definition: fe_base.hpp:724
int GetDerivRangeType() const
Returns the FiniteElement::RangeType of the element derivative, either SCALAR or VECTOR.
Definition: fe_base.hpp:344
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:205
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...
Definition: fe_base.cpp:40
void GetFaceMap(const int face_id, Array< int > &face_map) const override
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
Definition: fe_base.cpp:2524
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
Definition: fe_base.hpp:978
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe_base.hpp:349
static int VerifyNodal(int b_type)
Ensure that the BasisType of b_type nodal (satisfies the interpolation property). ...
Definition: fe_base.hpp:641
static int Pow(int base, int dim)
Return base raised to the power dim.
Definition: fe_base.hpp:1226
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:35
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
Definition: fe_base.cpp:662
Vector data type.
Definition: vector.hpp:58
virtual void GetFaceMap(const int face_id, Array< int > &face_map) const
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
Definition: fe_base.cpp:494
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 ...
Definition: fe_base.hpp:1152
static Geometry::Type GetTensorProductGeometry(int dim)
Definition: fe_base.hpp:1212
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe_base.cpp:969
void LocalRestriction_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
Definition: fe_base.cpp:1526
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition: fe_base.hpp:340
Describes the function space on each element.
Definition: fe_base.hpp:216
static const VectorFiniteElement & CheckVectorFE(const FiniteElement &fe)
Definition: fe_base.hpp:949
void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const override
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition: fe_base.cpp:710
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
Implements CalcDShape methods.
Definition: fe_base.hpp:295
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...
Definition: fe_base.cpp:288
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe_base.cpp:113
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:389
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.
Definition: fe_base.cpp:1204
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.
Definition: fe_base.hpp:1101
A Class that defines 1-D numerical quadrature rules on [0,1].
Definition: intrules.hpp:366
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
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...
Definition: fe_base.cpp:589
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...
Definition: fe_base.cpp:2055
const IntegrationRule * IntRule
IntegrationRule that defines the quadrature points at which the basis functions of the FE are evaluat...
Definition: fe_base.hpp:146
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...
Definition: fe_base.cpp:52
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe_base.cpp:1332
int order
Order/degree of the shape functions.
Definition: fe_base.hpp:243
Keep count of the number of eval types.
Definition: fe_base.hpp:973
int GetBasisType() const
Definition: fe_base.hpp:1202
No derivatives implemented.
Definition: fe_base.hpp:294
Implements CalcCurlShape methods.
Definition: fe_base.hpp:297
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...
Definition: fe_base.cpp:509