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