MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
fe_pos.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_POS
13#define MFEM_FE_POS
14
15#include "fe_base.hpp"
16
17namespace mfem
18{
19
20/** @brief Class for finite elements utilizing the
21 always positive Bernstein basis. */
23{
24public:
25 /** @brief Construct PositiveFiniteElement with given
26 @param D Reference space dimension
27 @param G Geometry type (of type Geometry::Type)
28 @param Do Number of degrees of freedom in the FiniteElement
29 @param O Order/degree of the FiniteElement
30 @param F FunctionSpace type of the FiniteElement
31 */
32 PositiveFiniteElement(int D, Geometry::Type G, int Do, int O,
33 int F = FunctionSpace::Pk) :
34 ScalarFiniteElement(D, G, Do, O, F)
35 { }
36
38 DenseMatrix &I) const
39 { ScalarLocalInterpolation(Trans, I, *this); }
40
42 DenseMatrix &R) const
43 { ScalarLocalL2Restriction(Trans, R, *this); }
44
45 virtual void GetTransferMatrix(const FiniteElement &fe,
47 DenseMatrix &I) const
48 { CheckScalarFE(fe).ScalarLocalInterpolation(Trans, I, *this); }
49
51
52 // Low-order monotone "projection" (actually it is not a projection): the
53 // dofs are set to be the Coefficient values at the nodes.
54 virtual void Project(Coefficient &coeff,
55 ElementTransformation &Trans, Vector &dofs) const;
56
57 virtual void Project (VectorCoefficient &vc,
58 ElementTransformation &Trans, Vector &dofs) const;
59
60 virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
61 DenseMatrix &I) const;
62};
63
64
67{
68public:
69 PositiveTensorFiniteElement(const int dims, const int p,
70 const DofMapType dmtype);
71
73 DofToQuad::Mode mode) const override
74 {
75 return (mode == DofToQuad::FULL) ?
77 GetTensorDofToQuad(*this, ir, mode, basis1d, true, dof2quad_array);
78 }
79
80 void GetFaceMap(const int face_id, Array<int> &face_map) const override;
81};
82
83
84/// A 2D positive bi-quadratic element on a square utilizing the 2nd order
85/// Bernstein basis
87{
88public:
89 /// Construct the BiQuadPos2DFiniteElement
91 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
92 virtual void CalcDShape(const IntegrationPoint &ip,
93 DenseMatrix &dshape) const;
95 DenseMatrix &I) const;
97 virtual void Project(Coefficient &coeff, ElementTransformation &Trans,
98 Vector &dofs) const;
99 virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans,
100 Vector &dofs) const;
101 virtual void ProjectDelta(int vertex, Vector &dofs) const
102 { dofs = 0.; dofs(vertex) = 1.; }
103};
104
105
106/// A 1D quadratic positive element utilizing the 2nd order Bernstein basis
108{
109public:
110 /// Construct the QuadPos1DFiniteElement
112 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
113 virtual void CalcDShape(const IntegrationPoint &ip,
114 DenseMatrix &dshape) const;
115};
116
117
118/// Arbitrary order H1 elements in 1D utilizing the Bernstein basis
120{
121private:
122#ifndef MFEM_THREAD_SAFE
123 // This is to share scratch space between invocations, which helps speed
124 // things up, but with OpenMP, we need one copy per thread. Right now, we
125 // solve this by allocating this space within each function call every time
126 // we call it. Alternatively, we should do some sort thread private thing.
127 // Brunner, Jan 2014
128 mutable Vector shape_x, dshape_x;
129#endif
130
131public:
132 /// Construct the H1Pos_SegmentElement of order @a p
133 H1Pos_SegmentElement(const int p);
134 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
135 virtual void CalcDShape(const IntegrationPoint &ip,
136 DenseMatrix &dshape) const;
137 virtual void ProjectDelta(int vertex, Vector &dofs) const;
138};
139
140
141/// Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square
143{
144private:
145#ifndef MFEM_THREAD_SAFE
146 // See comment in H1Pos_SegmentElement
147 mutable Vector shape_x, shape_y, dshape_x, dshape_y;
148#endif
149
150public:
151 /// Construct the H1Pos_QuadrilateralElement of order @a p
152 H1Pos_QuadrilateralElement(const int p);
153 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
154 virtual void CalcDShape(const IntegrationPoint &ip,
155 DenseMatrix &dshape) const;
156 virtual void ProjectDelta(int vertex, Vector &dofs) const;
157};
158
159
160/// Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube
162{
163private:
164#ifndef MFEM_THREAD_SAFE
165 // See comment in H1Pos_SegmentElement.
166 mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
167#endif
168
169public:
170 /// Construct the H1Pos_HexahedronElement of order @a p
171 H1Pos_HexahedronElement(const int p);
172 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
173 virtual void CalcDShape(const IntegrationPoint &ip,
174 DenseMatrix &dshape) const;
175 virtual void ProjectDelta(int vertex, Vector &dofs) const;
176};
177
178
179/// Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle
181{
182protected:
183#ifndef MFEM_THREAD_SAFE
186#endif
188
189public:
190 /// Construct the H1Pos_TriangleElement of order @a p
191 H1Pos_TriangleElement(const int p);
192
193 // The size of shape is (p+1)(p+2)/2 (dof).
194 static void CalcShape(const int p, const real_t x, const real_t y,
195 real_t *shape);
196
197 // The size of dshape_1d is p+1; the size of dshape is (dof x dim).
198 static void CalcDShape(const int p, const real_t x, const real_t y,
199 real_t *dshape_1d, real_t *dshape);
200
201 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
202 virtual void CalcDShape(const IntegrationPoint &ip,
203 DenseMatrix &dshape) const;
204};
205
206
207/// Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a
208/// tetrahedron
210{
211protected:
212#ifndef MFEM_THREAD_SAFE
215#endif
217
218public:
219 /// Construct the H1Pos_TetrahedronElement of order @a p
220 H1Pos_TetrahedronElement(const int p);
221
222 // The size of shape is (p+1)(p+2)(p+3)/6 (dof).
223 static void CalcShape(const int p, const real_t x, const real_t y,
224 const real_t z, real_t *shape);
225
226 // The size of dshape_1d is p+1; the size of dshape is (dof x dim).
227 static void CalcDShape(const int p, const real_t x, const real_t y,
228 const real_t z, real_t *dshape_1d, real_t *dshape);
229
230 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
231 virtual void CalcDShape(const IntegrationPoint &ip,
232 DenseMatrix &dshape) const;
233};
234
235
236/// Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a wedge
238{
239protected:
240#ifndef MFEM_THREAD_SAFE
243#endif
245
248
249public:
250 /// Construct the H1Pos_WedgeElement of order @a p
251 H1Pos_WedgeElement(const int p);
252
253 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
254 virtual void CalcDShape(const IntegrationPoint &ip,
255 DenseMatrix &dshape) const;
256};
257
258
259/// Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment
261{
262private:
263#ifndef MFEM_THREAD_SAFE
264 mutable Vector shape_x, dshape_x;
265#endif
266
267public:
268 /// Construct the L2Pos_SegmentElement of order @a p
269 L2Pos_SegmentElement(const int p);
270 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
271 virtual void CalcDShape(const IntegrationPoint &ip,
272 DenseMatrix &dshape) const;
273 virtual void ProjectDelta(int vertex, Vector &dofs) const;
274};
275
276
277/// Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square
279{
280private:
281#ifndef MFEM_THREAD_SAFE
282 mutable Vector shape_x, shape_y, dshape_x, dshape_y;
283#endif
284
285public:
286 /// Construct the L2Pos_QuadrilateralElement of order @a p
287 L2Pos_QuadrilateralElement(const int p);
288 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
289 virtual void CalcDShape(const IntegrationPoint &ip,
290 DenseMatrix &dshape) const;
291 virtual void ProjectDelta(int vertex, Vector &dofs) const;
292};
293
294
295/// Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube
297{
298private:
299#ifndef MFEM_THREAD_SAFE
300 mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
301#endif
302
303public:
304 /// Construct the L2Pos_HexahedronElement of order @a p
305 L2Pos_HexahedronElement(const int p);
306 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
307 virtual void CalcDShape(const IntegrationPoint &ip,
308 DenseMatrix &dshape) const;
309 virtual void ProjectDelta(int vertex, Vector &dofs) const;
310};
311
312
313/// Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle
315{
316private:
317#ifndef MFEM_THREAD_SAFE
318 mutable Vector dshape_1d;
319#endif
320
321public:
322 /// Construct the L2Pos_TriangleElement of order @a p
323 L2Pos_TriangleElement(const int p);
324 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
325 virtual void CalcDShape(const IntegrationPoint &ip,
326 DenseMatrix &dshape) const;
327 virtual void ProjectDelta(int vertex, Vector &dofs) const;
328};
329
330
331/// Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a
332/// tetrahedron
334{
335private:
336#ifndef MFEM_THREAD_SAFE
337 mutable Vector dshape_1d;
338#endif
339
340public:
341 /// Construct the L2Pos_TetrahedronElement of order @a p
342 L2Pos_TetrahedronElement(const int p);
343 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
344 virtual void CalcDShape(const IntegrationPoint &ip,
345 DenseMatrix &dshape) const;
346 virtual void ProjectDelta(int vertex, Vector &dofs) const;
347};
348
349
350/// Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a wedge
352{
353protected:
354#ifndef MFEM_THREAD_SAFE
357#endif
359
362
363public:
364 /// Construct the L2Pos_WedgeElement of order @a p
365 L2Pos_WedgeElement(const int p);
366
367 virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
368 virtual void CalcDShape(const IntegrationPoint &ip,
369 DenseMatrix &dshape) const;
370};
371
372} // namespace mfem
373
374#endif
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_pos.cpp:220
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:118
BiQuadPos2DFiniteElement()
Construct the BiQuadPos2DFiniteElement.
Definition fe_pos.cpp:95
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_pos.hpp:101
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:142
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_pos.cpp:191
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
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
Abstract class for all finite elements.
Definition fe_base.hpp:239
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
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
Array< DofToQuad * > dof2quad_array
Container for all DofToQuad objects created by the FiniteElement.
Definition fe_base.hpp:258
@ Pk
Polynomials of order k.
Definition fe_base.hpp:226
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube.
Definition fe_pos.hpp:162
H1Pos_HexahedronElement(const int p)
Construct the H1Pos_HexahedronElement of order p.
Definition fe_pos.cpp:432
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:475
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:454
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_pos.cpp:499
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square.
Definition fe_pos.hpp:143
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:404
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:384
H1Pos_QuadrilateralElement(const int p)
Construct the H1Pos_QuadrilateralElement of order p.
Definition fe_pos.cpp:364
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_pos.cpp:425
Arbitrary order H1 elements in 1D utilizing the Bernstein basis.
Definition fe_pos.hpp:120
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_pos.cpp:357
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:317
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:337
H1Pos_SegmentElement(const int p)
Construct the H1Pos_SegmentElement of order p.
Definition fe_pos.cpp:299
H1Pos_TetrahedronElement(const int p)
Construct the H1Pos_TetrahedronElement of order p.
Definition fe_pos.cpp:651
static void CalcShape(const int p, const real_t x, const real_t y, const real_t z, real_t *shape)
Definition fe_pos.cpp:753
static void CalcDShape(const int p, const real_t x, const real_t y, const real_t z, real_t *dshape_1d, real_t *dshape)
Definition fe_pos.cpp:786
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle.
Definition fe_pos.hpp:181
H1Pos_TriangleElement(const int p)
Construct the H1Pos_TriangleElement of order p.
Definition fe_pos.cpp:506
static void CalcDShape(const int p, const real_t x, const real_t y, real_t *dshape_1d, real_t *dshape)
Definition fe_pos.cpp:587
static void CalcShape(const int p, const real_t x, const real_t y, real_t *shape)
Definition fe_pos.cpp:561
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a wedge.
Definition fe_pos.hpp:238
H1Pos_TriangleElement TriangleFE
Definition fe_pos.hpp:246
H1Pos_SegmentElement SegmentFE
Definition fe_pos.hpp:247
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1020
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1001
H1Pos_WedgeElement(const int p)
Construct the H1Pos_WedgeElement of order p.
Definition fe_pos.cpp:903
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube.
Definition fe_pos.hpp:297
L2Pos_HexahedronElement(const int p)
Construct the L2Pos_HexahedronElement of order p.
Definition fe_pos.cpp:1168
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1195
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1216
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_pos.cpp:1240
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square.
Definition fe_pos.hpp:279
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1114
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_pos.cpp:1153
L2Pos_QuadrilateralElement(const int p)
Construct the L2Pos_QuadrilateralElement of order p.
Definition fe_pos.cpp:1090
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1133
Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment.
Definition fe_pos.hpp:261
L2Pos_SegmentElement(const int p)
Construct the L2Pos_SegmentElement of order p.
Definition fe_pos.cpp:1045
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1066
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1072
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_pos.cpp:1083
L2Pos_TetrahedronElement(const int p)
Construct the L2Pos_TetrahedronElement of order p.
Definition fe_pos.cpp:1310
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1333
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1340
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_pos.cpp:1351
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle.
Definition fe_pos.hpp:315
L2Pos_TriangleElement(const int p)
Construct the L2Pos_TriangleElement of order p.
Definition fe_pos.cpp:1259
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_pos.cpp:1298
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1287
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1281
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a wedge.
Definition fe_pos.hpp:352
L2Pos_TriangleElement TriangleFE
Definition fe_pos.hpp:360
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:1426
L2Pos_WedgeElement(const int p)
Construct the L2Pos_WedgeElement of order p.
Definition fe_pos.cpp:1364
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:1407
L2Pos_SegmentElement SegmentFE
Definition fe_pos.hpp:361
Class for finite elements utilizing the always positive Bernstein basis.
Definition fe_pos.hpp:23
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_pos.hpp:41
PositiveFiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct PositiveFiniteElement with given.
Definition fe_pos.hpp:32
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_pos.cpp:25
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_pos.hpp:37
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_pos.hpp:45
PositiveTensorFiniteElement(const int dims, const int p, const DofMapType dmtype)
Definition fe_pos.cpp:81
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_pos.cpp:88
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_pos.hpp:72
A 1D quadratic positive element utilizing the 2nd order Bernstein basis.
Definition fe_pos.hpp:108
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_pos.cpp:278
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_pos.cpp:288
QuadPos1DFiniteElement()
Construct the QuadPos1DFiniteElement.
Definition fe_pos.cpp:270
Class for finite elements with basis functions that return scalar values.
Definition fe_base.hpp:656
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
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
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:80
float real_t
Definition config.hpp:43
real_t p(const Vector &x, real_t t)