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