MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
tfe.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_TEMPLATE_FINITE_ELEMENTS
13#define MFEM_TEMPLATE_FINITE_ELEMENTS
14
15#include "../config/tconfig.hpp"
16#include "fe_coll.hpp"
17
18namespace mfem
19{
20
21// Templated finite element classes, cf. fe.?pp
22
23/** @brief Store mass-like matrix B for each integration point on the reference
24 element.
25 For tensor product evaluation, this is only called on the 1D reference
26 element, and higher dimensions are put together from that.
27 The element mass matrix can be written $ M_E = B^T D_E B $ where the B
28 built here is the B, and is unchanging across the mesh. The diagonal matrix
29 $ D_E $ then contains all the element-specific geometry and physics data.
30 @param fe the element we are calculating on
31 @param ir the integration rule to calculate the shape matrix on
32 @param B must be (nip x dof) with column major storage
33 @param dof_map the inverse of dof_map is applied to reorder local dofs.
34*/
35template <typename real_t>
37 real_t *B, const Array<int> *dof_map = NULL)
38{
39 // - B must be (nip x dof) with column major storage
40 // - The inverse of dof_map is applied to reorder the local dofs.
41 int nip = ir.GetNPoints();
42 int dof = fe.GetDof();
43 Vector shape(dof);
44
45 for (int ip = 0; ip < nip; ip++)
46 {
47 fe.CalcShape(ir.IntPoint(ip), shape);
48 for (int id = 0; id < dof; id++)
49 {
50 int orig_id = dof_map ? (*dof_map)[id] : id;
51 B[ip+nip*id] = shape(orig_id);
52 }
53 }
54}
55
56/** @brief store gradient matrix G for each integration point on the reference
57 element.
58 For tensor product evaluation, this is only called on the 1D reference
59 element, and higher dimensions are put together from that.
60 The element stiffness matrix can be written
61 $$
62 S_E = \sum_{k=1}^{nq} G_{k,i}^T (D_E^G)_{k,k} G_{k,j}
63 $$
64 where $ nq $ is the number of quadrature points, $ D_E^G $ contains
65 all the information about the element geometry and coefficients (Jacobians
66 etc.), and $ G $ is the matrix built in this routine, which is the same
67 for all elements in a mesh.
68 @param fe the element we are calculating on
69 @param ir the integration rule to calculate the gradients on
70 @param[out] G must be (nip x dim x dof) with column major storage
71 @param[in] dof_map the inverse of dof_map is applied to reorder local dofs.
72*/
73template <typename real_t>
75 real_t *G, const Array<int> *dof_map = NULL)
76{
77 // - G must be (nip x dim x dof) with column major storage
78 // - The inverse of dof_map is applied to reorder the local dofs.
79 int dim = fe.GetDim();
80 int nip = ir.GetNPoints();
81 int dof = fe.GetDof();
82 DenseMatrix dshape(dof, dim);
83
84 for (int ip = 0; ip < nip; ip++)
85 {
86 fe.CalcDShape(ir.IntPoint(ip), dshape);
87 for (int id = 0; id < dof; id++)
88 {
89 int orig_id = dof_map ? (*dof_map)[id] : id;
90 for (int d = 0; d < dim; d++)
91 {
92 G[ip+nip*(d+dim*id)] = dshape(orig_id, d);
93 }
94 }
95 }
96}
97
98template <typename real_t>
99void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir,
100 real_t *B, real_t *G, const Array<int> *dof_map)
101{
102 if (B) { mfem::CalcShapeMatrix(fe, ir, B, dof_map); }
103 if (G) { mfem::CalcGradTensor(fe, ir, G, dof_map); }
104}
105
106// H1 finite elements
107
108template <Geometry::Type G, int P>
110
111template <int P>
112class H1_FiniteElement<Geometry::SEGMENT, P>
113{
114public:
116 static const int dim = 1;
117 static const int degree = P;
118 static const int dofs = P+1;
119
120 static const bool tensor_prod = true;
121 static const int dofs_1d = P+1;
122
123 // Type for run-time parameter for the constructor
124 typedef int parameter_type;
125
126protected:
129 parameter_type type; // run-time specified basis type
130 void Init(const parameter_type type_)
131 {
132 type = type_;
133 if (type == BasisType::Positive)
134 {
136 my_fe = fe;
137 my_dof_map = &fe->GetDofMap();
138 }
139 else
140 {
141 int pt_type = BasisType::GetQuadrature1D(type);
142 H1_SegmentElement *fe = new H1_SegmentElement(P, pt_type);
143 my_fe = fe;
144 my_dof_map = &fe->GetDofMap();
145 }
146 }
147
148public:
150 {
151 Init(type_);
152 }
154 {
155 const H1_FECollection *h1_fec =
156 dynamic_cast<const H1_FECollection *>(&fec);
157 MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
158 Init(h1_fec->GetBasisType());
159 }
160 ~H1_FiniteElement() { delete my_fe; }
161
162 template <typename real_t>
163 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
164 {
165 mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
166 }
167 template <typename real_t>
168 void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
169 {
170 CalcShapes(ir, B, G);
171 }
172 const Array<int> *GetDofMap() const { return my_dof_map; }
173};
174
175template <int P>
176class H1_FiniteElement<Geometry::TRIANGLE, P>
177{
178public:
180 static const int dim = 2;
181 static const int degree = P;
182 static const int dofs = ((P + 1)*(P + 2))/2;
183
184 static const bool tensor_prod = false;
185
186 // Type for run-time parameter for the constructor
187 typedef int parameter_type;
188
189protected:
191 parameter_type type; // run-time specified basis type
192 void Init(const parameter_type type_)
193 {
194 type = type_;
195 if (type == BasisType::Positive)
196 {
197 my_fe = new H1Pos_TriangleElement(P);
198 }
199 else
200 {
201 int pt_type = BasisType::GetQuadrature1D(type);
202 my_fe = new H1_TriangleElement(P, pt_type);
203 }
204 }
205
206public:
208 {
209 Init(type_);
210 }
212 {
213 const H1_FECollection *h1_fec =
214 dynamic_cast<const H1_FECollection *>(&fec);
215 MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
216 Init(h1_fec->GetBasisType());
217 }
218 ~H1_FiniteElement() { delete my_fe; }
219
220 template <typename real_t>
221 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
222 {
223 mfem::CalcShapes(*my_fe, ir, B, G, NULL);
224 }
225 const Array<int> *GetDofMap() const { return NULL; }
226};
227
228template <int P>
229class H1_FiniteElement<Geometry::SQUARE, P>
230{
231public:
232 static const Geometry::Type geom = Geometry::SQUARE;
233 static const int dim = 2;
234 static const int degree = P;
235 static const int dofs = (P+1)*(P+1);
236
237 static const bool tensor_prod = true;
238 static const int dofs_1d = P+1;
239
240 // Type for run-time parameter for the constructor
241 typedef int parameter_type;
242
243protected:
244 const FiniteElement *my_fe, *my_fe_1d;
246 parameter_type type; // run-time specified basis type
247 void Init(const parameter_type type_)
248 {
249 type = type_;
250 if (type == BasisType::Positive)
251 {
253 my_fe = fe;
254 my_dof_map = &fe->GetDofMap();
255 my_fe_1d = new L2Pos_SegmentElement(P);
256 }
257 else
258 {
259 int pt_type = BasisType::GetQuadrature1D(type);
261 my_fe = fe;
262 my_dof_map = &fe->GetDofMap();
263 my_fe_1d = new L2_SegmentElement(P, pt_type);
264 }
265 }
266
267public:
269 {
270 Init(type_);
271 }
273 {
274 const H1_FECollection *h1_fec =
275 dynamic_cast<const H1_FECollection *>(&fec);
276 MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
277 Init(h1_fec->GetBasisType());
278 }
279 ~H1_FiniteElement() { delete my_fe; delete my_fe_1d; }
280
281 template <typename real_t>
282 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
283 {
284 mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
285 }
286 template <typename real_t>
287 void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
288 {
289 mfem::CalcShapes(*my_fe_1d, ir, B, G, NULL);
290 }
291 const Array<int> *GetDofMap() const { return my_dof_map; }
292};
293
294template <int P>
295class H1_FiniteElement<Geometry::TETRAHEDRON, P>
296{
297public:
299 static const int dim = 3;
300 static const int degree = P;
301 static const int dofs = ((P + 1)*(P + 2)*(P + 3))/6;
302
303 static const bool tensor_prod = false;
304
305 // Type for run-time parameter for the constructor
306 typedef int parameter_type;
307
308protected:
310 parameter_type type; // run-time specified basis type
311 void Init(const parameter_type type_)
312 {
313 type = type_;
314 if (type == BasisType::Positive)
315 {
316 my_fe = new H1Pos_TetrahedronElement(P);
317 }
318 else
319 {
320 int pt_type = BasisType::GetQuadrature1D(type);
321 my_fe = new H1_TetrahedronElement(P, pt_type);
322 }
323 }
324
325public:
327 {
328 Init(type_);
329 }
331 {
332 const H1_FECollection *h1_fec =
333 dynamic_cast<const H1_FECollection *>(&fec);
334 MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
335 Init(h1_fec->GetBasisType());
336 }
337 ~H1_FiniteElement() { delete my_fe; }
338
339 template <typename real_t>
340 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
341 {
342 mfem::CalcShapes(*my_fe, ir, B, G, NULL);
343 }
344 const Array<int> *GetDofMap() const { return NULL; }
345};
346
347template <int P>
349{
350public:
351 static const Geometry::Type geom = Geometry::CUBE;
352 static const int dim = 3;
353 static const int degree = P;
354 static const int dofs = (P+1)*(P+1)*(P+1);
355
356 static const bool tensor_prod = true;
357 static const int dofs_1d = P+1;
358
359 // Type for run-time parameter for the constructor
360 typedef int parameter_type;
361
362protected:
363 const FiniteElement *my_fe, *my_fe_1d;
365 parameter_type type; // run-time specified basis type
366
367 void Init(const parameter_type type_)
368 {
369 type = type_;
370 if (type == BasisType::Positive)
371 {
373 my_fe = fe;
374 my_dof_map = &fe->GetDofMap();
375 my_fe_1d = new L2Pos_SegmentElement(P);
376 }
377 else
378 {
379 int pt_type = BasisType::GetQuadrature1D(type);
380 H1_HexahedronElement *fe = new H1_HexahedronElement(P, pt_type);
381 my_fe = fe;
382 my_dof_map = &fe->GetDofMap();
383 my_fe_1d = new L2_SegmentElement(P, pt_type);
384 }
385 }
386
387public:
389 {
390 Init(type_);
391 }
393 {
394 const H1_FECollection *h1_fec =
395 dynamic_cast<const H1_FECollection *>(&fec);
396 MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
397 Init(h1_fec->GetBasisType());
398 }
399 ~H1_FiniteElement() { delete my_fe; delete my_fe_1d; }
400
401 template <typename real_t>
402 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
403 {
404 mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
405 }
406 template <typename real_t>
407 void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
408 {
409 mfem::CalcShapes(*my_fe_1d, ir, B, G, NULL);
410 }
411 const Array<int> *GetDofMap() const { return my_dof_map; }
412};
413
414
415// L2 finite elements
416
417template <Geometry::Type G, int P, typename L2_FE_type, typename L2Pos_FE_type,
418 int DOFS, bool TP>
420{
421public:
422 static const Geometry::Type geom = G;
424 static const int degree = P;
425 static const int dofs = DOFS;
426
427 static const bool tensor_prod = TP;
428 static const int dofs_1d = P+1;
429
430 // Type for run-time parameter for the constructor
431 typedef int parameter_type;
432
433protected:
435 parameter_type type; // run-time specified basis type
436
437 void Init(const parameter_type type_)
438 {
439 type = type_;
440 switch (type)
441 {
443 my_fe = new L2Pos_FE_type(P);
444 my_fe_1d = (TP && dim != 1) ? new L2Pos_SegmentElement(P) : NULL;
445 break;
446
447 default:
448 int pt_type = BasisType::GetQuadrature1D(type);
449 my_fe = new L2_FE_type(P, pt_type);
450 my_fe_1d = (TP && dim != 1) ? new L2_SegmentElement(P, pt_type) :
451 NULL;
452 }
453 }
454
457
459 {
460 const L2_FECollection *l2_fec =
461 dynamic_cast<const L2_FECollection *>(&fec);
462 MFEM_ASSERT(l2_fec, "invalid FiniteElementCollection");
463 Init(l2_fec->GetBasisType());
464 }
465
467
468public:
469 template <typename real_t>
470 void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
471 {
472 mfem::CalcShapes(*my_fe, ir, B, Grad, NULL);
473 }
474 template <typename real_t>
475 void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
476 {
477 mfem::CalcShapes(dim == 1 ? *my_fe : *my_fe_1d, ir, B, Grad, NULL);
478 }
479 const Array<int> *GetDofMap() const { return NULL; }
480};
481
482
483template <Geometry::Type G, int P>
485
486
487template <int P>
488class L2_FiniteElement<Geometry::SEGMENT, P>
489 : public L2_FiniteElement_base<
490 Geometry::SEGMENT,P,L2_SegmentElement,L2Pos_SegmentElement,P+1,true>
491{
492protected:
495public:
501};
502
503
504template <int P>
505class L2_FiniteElement<Geometry::TRIANGLE, P>
506 : public L2_FiniteElement_base<Geometry::TRIANGLE,P,L2_TriangleElement,
507 L2Pos_TriangleElement,((P+1)*(P+2))/2,false>
508{
509protected:
511 L2Pos_TriangleElement,((P+1)*(P+2))/2,false> base_class;
512public:
518};
519
520
521template <int P>
522class L2_FiniteElement<Geometry::SQUARE, P>
523 : public L2_FiniteElement_base<Geometry::SQUARE,P,L2_QuadrilateralElement,
524 L2Pos_QuadrilateralElement,(P+1)*(P+1),true>
525{
526protected:
529public:
535};
536
537
538template <int P>
539class L2_FiniteElement<Geometry::TETRAHEDRON, P>
540 : public L2_FiniteElement_base<Geometry::TETRAHEDRON,P,L2_TetrahedronElement,
541 L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6,false>
542{
543protected:
545 L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6,false> base_class;
546public:
552};
553
554
555template <int P>
557 : public L2_FiniteElement_base<Geometry::CUBE,P,L2_HexahedronElement,
558 L2Pos_HexahedronElement,(P+1)*(P+1)*(P+1),true>
559{
560protected:
562 L2Pos_HexahedronElement,(P+1)*(P+1)*(P+1),true> base_class;
563public:
569};
570
571} // namespace mfem
572
573#endif // MFEM_TEMPLATE_FINITE_ELEMENTS
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
Definition fe_base.hpp:61
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ GaussLegendre
Open type.
Definition fe_base.hpp:31
@ Positive
Bernstein polynomials.
Definition fe_base.hpp:33
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:316
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...
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...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube.
Definition fe_pos.hpp:162
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square.
Definition fe_pos.hpp:143
Arbitrary order H1 elements in 1D utilizing the Bernstein basis.
Definition fe_pos.hpp:120
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle.
Definition fe_pos.hpp:181
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
int GetBasisType() const
Definition fe_coll.hpp:285
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Definition tfe.hpp:388
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:407
const Array< int > * GetDofMap() const
Definition tfe.hpp:411
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:402
void Init(const parameter_type type_)
Definition tfe.hpp:367
H1_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:392
H1_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:153
void Init(const parameter_type type_)
Definition tfe.hpp:130
const Array< int > * GetDofMap() const
Definition tfe.hpp:172
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:168
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:163
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Definition tfe.hpp:149
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Definition tfe.hpp:268
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:287
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:282
void Init(const parameter_type type_)
Definition tfe.hpp:247
const Array< int > * GetDofMap() const
Definition tfe.hpp:291
H1_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:272
const Array< int > * GetDofMap() const
Definition tfe.hpp:344
H1_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:330
void Init(const parameter_type type_)
Definition tfe.hpp:311
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:340
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Definition tfe.hpp:326
H1_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:211
H1_FiniteElement(const parameter_type type_=BasisType::GaussLobatto)
Definition tfe.hpp:207
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition tfe.hpp:221
void Init(const parameter_type type_)
Definition tfe.hpp:192
const Array< int > * GetDofMap() const
Definition tfe.hpp:225
Arbitrary order H1 elements in 3D on a cube.
Definition fe_h1.hpp:63
Arbitrary order H1 elements in 2D on a square.
Definition fe_h1.hpp:42
Arbitrary order H1 elements in 1D.
Definition fe_h1.hpp:22
Arbitrary order H1 elements in 3D on a tetrahedron.
Definition fe_h1.hpp:106
Arbitrary order H1 elements in 2D on a triangle.
Definition fe_h1.hpp:84
A class to initialize the size of a Tensor.
Definition dtensor.hpp:55
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube.
Definition fe_pos.hpp:297
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square.
Definition fe_pos.hpp:279
Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment.
Definition fe_pos.hpp:261
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle.
Definition fe_pos.hpp:315
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
int GetBasisType() const
Definition fe_coll.hpp:373
L2_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:567
base_class::parameter_type parameter_type
Definition tfe.hpp:564
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
Definition tfe.hpp:565
L2_FiniteElement_base< Geometry::SEGMENT, P, L2_SegmentElement, L2Pos_SegmentElement, P+1, true > base_class
Definition tfe.hpp:494
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
Definition tfe.hpp:497
L2_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:499
base_class::parameter_type parameter_type
Definition tfe.hpp:496
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
Definition tfe.hpp:531
L2_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:533
L2_FiniteElement_base< Geometry::SQUARE, P, L2_QuadrilateralElement, L2Pos_QuadrilateralElement,(P+1) *(P+1), true > base_class
Definition tfe.hpp:528
base_class::parameter_type parameter_type
Definition tfe.hpp:530
base_class::parameter_type parameter_type
Definition tfe.hpp:547
L2_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:550
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
Definition tfe.hpp:548
L2_FiniteElement_base< Geometry::TETRAHEDRON, P, L2_TetrahedronElement, L2Pos_TetrahedronElement,((P+1) *(P+2) *(P+3))/6, false > base_class
Definition tfe.hpp:545
base_class::parameter_type parameter_type
Definition tfe.hpp:513
L2_FiniteElement(const parameter_type type_=BasisType::GaussLegendre)
Definition tfe.hpp:514
L2_FiniteElement(const FiniteElementCollection &fec)
Definition tfe.hpp:516
const FiniteElement * my_fe
Definition tfe.hpp:434
static const bool tensor_prod
Definition tfe.hpp:427
parameter_type type
Definition tfe.hpp:435
static const int dim
Definition tfe.hpp:423
static const int dofs
Definition tfe.hpp:425
L2_FiniteElement_base(const FiniteElementCollection &fec)
Definition tfe.hpp:458
static const int degree
Definition tfe.hpp:424
const Array< int > * GetDofMap() const
Definition tfe.hpp:479
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
Definition tfe.hpp:475
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
Definition tfe.hpp:470
L2_FiniteElement_base(const parameter_type type)
Definition tfe.hpp:455
const FiniteElement * my_fe_1d
Definition tfe.hpp:434
static const Geometry::Type geom
Definition tfe.hpp:422
void Init(const parameter_type type_)
Definition tfe.hpp:437
static const int dofs_1d
Definition tfe.hpp:428
Arbitrary order L2 elements in 3D on a cube.
Definition fe_l2.hpp:79
Arbitrary order L2 elements in 2D on a square.
Definition fe_l2.hpp:45
Arbitrary order L2 elements in 1D on a segment.
Definition fe_l2.hpp:22
Arbitrary order L2 elements in 3D on a tetrahedron.
Definition fe_l2.hpp:139
Arbitrary order L2 elements in 2D on a triangle.
Definition fe_l2.hpp:109
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
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
void CalcGradTensor(const FiniteElement &fe, const IntegrationRule &ir, real_t *G, const Array< int > *dof_map=NULL)
store gradient matrix G for each integration point on the reference element. For tensor product evalu...
Definition tfe.hpp:74
void CalcShapeMatrix(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, const Array< int > *dof_map=NULL)
Store mass-like matrix B for each integration point on the reference element. For tensor product eval...
Definition tfe.hpp:36
void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, real_t *G, const Array< int > *dof_map)
Definition tfe.hpp:99
float real_t
Definition config.hpp:43