MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tfe.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
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 
18 namespace mfem
19 {
20 
21 // Templated finite element classes, cf. fe.?pp
22 
23 template <typename real_t>
24 void CalcShapeMatrix(const FiniteElement &fe, const IntegrationRule &ir,
25  real_t *B, const Array<int> *dof_map = NULL)
26 {
27  // - B must be (nip x dof) with column major storage
28  // - The inverse of dof_map is applied to reorder the local dofs.
29  int nip = ir.GetNPoints();
30  int dof = fe.GetDof();
31  Vector shape(dof);
32 
33  for (int ip = 0; ip < nip; ip++)
34  {
35  fe.CalcShape(ir.IntPoint(ip), shape);
36  for (int id = 0; id < dof; id++)
37  {
38  int orig_id = dof_map ? (*dof_map)[id] : id;
39  B[ip+nip*id] = shape(orig_id);
40  }
41  }
42 }
43 
44 template <typename real_t>
45 void CalcGradTensor(const FiniteElement &fe, const IntegrationRule &ir,
46  real_t *G, const Array<int> *dof_map = NULL)
47 {
48  // - G must be (nip x dim x dof) with column major storage
49  // - The inverse of dof_map is applied to reorder the local dofs.
50  int dim = fe.GetDim();
51  int nip = ir.GetNPoints();
52  int dof = fe.GetDof();
53  DenseMatrix dshape(dof, dim);
54 
55  for (int ip = 0; ip < nip; ip++)
56  {
57  fe.CalcDShape(ir.IntPoint(ip), dshape);
58  for (int id = 0; id < dof; id++)
59  {
60  int orig_id = dof_map ? (*dof_map)[id] : id;
61  for (int d = 0; d < dim; d++)
62  {
63  G[ip+nip*(d+dim*id)] = dshape(orig_id, d);
64  }
65  }
66  }
67 }
68 
69 template <typename real_t>
70 void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir,
71  real_t *B, real_t *G, const Array<int> *dof_map)
72 {
73  if (B) { mfem::CalcShapeMatrix(fe, ir, B, dof_map); }
74  if (G) { mfem::CalcGradTensor(fe, ir, G, dof_map); }
75 }
76 
77 // H1 finite elements
78 
79 template <Geometry::Type G, int P>
81 
82 template <int P>
83 class H1_FiniteElement<Geometry::SEGMENT, P>
84 {
85 public:
87  static const int dim = 1;
88  static const int degree = P;
89  static const int dofs = P+1;
90 
91  static const bool tensor_prod = true;
92  static const int dofs_1d = P+1;
93 
94  // Type for run-time parameter for the constructor
96 
97 protected:
100  parameter_type type; // run-time specified basis type
101  void Init(const parameter_type type_)
102  {
103  type = type_;
104  if (type == H1_FECollection::GaussLobatto)
105  {
107  my_fe = fe;
108  my_dof_map = &fe->GetDofMap();
109  }
110  else if (type == H1_FECollection::Positive)
111  {
113  my_fe = fe;
114  my_dof_map = &fe->GetDofMap();
115  }
116  else
117  {
118  MFEM_ABORT("invalid basis type!");
119  }
120  }
121 
122 public:
124  {
125  Init(type_);
126  }
128  {
129  const H1_FECollection *h1_fec =
130  dynamic_cast<const H1_FECollection *>(&fec);
131  MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
132  Init(h1_fec->GetBasisType());
133  }
134  ~H1_FiniteElement() { delete my_fe; }
135 
136  template <typename real_t>
137  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
138  {
139  mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
140  }
141  template <typename real_t>
142  void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
143  {
144  CalcShapes(ir, B, G);
145  }
146  const Array<int> *GetDofMap() const { return my_dof_map; }
147 };
148 
149 template <int P>
150 class H1_FiniteElement<Geometry::TRIANGLE, P>
151 {
152 public:
154  static const int dim = 2;
155  static const int degree = P;
156  static const int dofs = ((P + 1)*(P + 2))/2;
157 
158  static const bool tensor_prod = false;
159 
160  // Type for run-time parameter for the constructor
162 
163 protected:
165  parameter_type type; // run-time specified basis type
166  void Init(const parameter_type type_)
167  {
168  type = type_;
169  if (type == H1_FECollection::GaussLobatto)
170  {
171  my_fe = new H1_TriangleElement(P);
172  }
173  else if (type == H1_FECollection::Positive)
174  {
175  my_fe = new H1Pos_TriangleElement(P);
176  }
177  else
178  {
179  MFEM_ABORT("invalid basis type!");
180  }
181  }
182 
183 public:
185  {
186  Init(type_);
187  }
189  {
190  const H1_FECollection *h1_fec =
191  dynamic_cast<const H1_FECollection *>(&fec);
192  MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
193  Init(h1_fec->GetBasisType());
194  }
195  ~H1_FiniteElement() { delete my_fe; }
196 
197  template <typename real_t>
198  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
199  {
200  mfem::CalcShapes(*my_fe, ir, B, G, NULL);
201  }
202  const Array<int> *GetDofMap() const { return NULL; }
203 };
204 
205 template <int P>
206 class H1_FiniteElement<Geometry::SQUARE, P>
207 {
208 public:
210  static const int dim = 2;
211  static const int degree = P;
212  static const int dofs = (P+1)*(P+1);
213 
214  static const bool tensor_prod = true;
215  static const int dofs_1d = P+1;
216 
217  // Type for run-time parameter for the constructor
219 
220 protected:
221  const FiniteElement *my_fe, *my_fe_1d;
223  parameter_type type; // run-time specified basis type
224  void Init(const parameter_type type_)
225  {
226  type = type_;
227  if (type == H1_FECollection::GaussLobatto)
228  {
230  my_fe = fe;
231  my_dof_map = &fe->GetDofMap();
232  my_fe_1d = new L2_SegmentElement(P, 1);
233  }
234  else if (type == H1_FECollection::Positive)
235  {
237  my_fe = fe;
238  my_dof_map = &fe->GetDofMap();
239  my_fe_1d = new L2Pos_SegmentElement(P);
240  }
241  else
242  {
243  MFEM_ABORT("invalid basis type!");
244  }
245  }
246 
247 public:
249  {
250  Init(type_);
251  }
253  {
254  const H1_FECollection *h1_fec =
255  dynamic_cast<const H1_FECollection *>(&fec);
256  MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
257  Init(h1_fec->GetBasisType());
258  }
259  ~H1_FiniteElement() { delete my_fe; delete my_fe_1d; }
260 
261  template <typename real_t>
262  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
263  {
264  mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
265  }
266  template <typename real_t>
267  void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
268  {
269  mfem::CalcShapes(*my_fe_1d, ir, B, G, NULL);
270  }
271  const Array<int> *GetDofMap() const { return my_dof_map; }
272 };
273 
274 template <int P>
275 class H1_FiniteElement<Geometry::TETRAHEDRON, P>
276 {
277 public:
279  static const int dim = 3;
280  static const int degree = P;
281  static const int dofs = ((P + 1)*(P + 2)*(P + 3))/6;
282 
283  static const bool tensor_prod = false;
284 
285  // Type for run-time parameter for the constructor
287 
288 protected:
290  parameter_type type; // run-time specified basis type
291  void Init(const parameter_type type_)
292  {
293  type = type_;
294  if (type == H1_FECollection::GaussLobatto)
295  {
296  my_fe = new H1_TetrahedronElement(P);
297  }
298  else if (type == H1_FECollection::Positive)
299  {
300  my_fe = new H1Pos_TetrahedronElement(P);
301  }
302  else
303  {
304  MFEM_ABORT("invalid basis type!");
305  }
306  }
307 
308 public:
310  {
311  Init(type_);
312  }
314  {
315  const H1_FECollection *h1_fec =
316  dynamic_cast<const H1_FECollection *>(&fec);
317  MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
318  Init(h1_fec->GetBasisType());
319  }
320  ~H1_FiniteElement() { delete my_fe; }
321 
322  template <typename real_t>
323  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
324  {
325  mfem::CalcShapes(*my_fe, ir, B, G, NULL);
326  }
327  const Array<int> *GetDofMap() const { return NULL; }
328 };
329 
330 template <int P>
331 class H1_FiniteElement<Geometry::CUBE, P>
332 {
333 public:
335  static const int dim = 3;
336  static const int degree = P;
337  static const int dofs = (P+1)*(P+1)*(P+1);
338 
339  static const bool tensor_prod = true;
340  static const int dofs_1d = P+1;
341 
342  // Type for run-time parameter for the constructor
344 
345 protected:
346  const FiniteElement *my_fe, *my_fe_1d;
348  parameter_type type; // run-time specified basis type
349 
350  void Init(const parameter_type type_)
351  {
352  type = type_;
353  if (type == H1_FECollection::GaussLobatto)
354  {
356  my_fe = fe;
357  my_dof_map = &fe->GetDofMap();
358  my_fe_1d = new L2_SegmentElement(P, 1);
359  }
360  else if (type == H1_FECollection::Positive)
361  {
363  my_fe = fe;
364  my_dof_map = &fe->GetDofMap();
365  my_fe_1d = new L2Pos_SegmentElement(P);
366  }
367  else
368  {
369  MFEM_ABORT("invalid basis type!");
370  }
371  }
372 
373 public:
375  {
376  Init(type_);
377  }
379  {
380  const H1_FECollection *h1_fec =
381  dynamic_cast<const H1_FECollection *>(&fec);
382  MFEM_ASSERT(h1_fec, "invalid FiniteElementCollection");
383  Init(h1_fec->GetBasisType());
384  }
385  ~H1_FiniteElement() { delete my_fe; delete my_fe_1d; }
386 
387  template <typename real_t>
388  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
389  {
390  mfem::CalcShapes(*my_fe, ir, B, G, my_dof_map);
391  }
392  template <typename real_t>
393  void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
394  {
395  mfem::CalcShapes(*my_fe_1d, ir, B, G, NULL);
396  }
397  const Array<int> *GetDofMap() const { return my_dof_map; }
398 };
399 
400 
401 // L2 finite elements
402 
403 template <Geometry::Type G, int P, typename L2_FE_type, typename L2Pos_FE_type,
404  int DOFS, bool TP>
406 {
407 public:
408  static const Geometry::Type geom = G;
410  static const int degree = P;
411  static const int dofs = DOFS;
412 
413  static const bool tensor_prod = TP;
414  static const int dofs_1d = P+1;
415 
416  // Type for run-time parameter for the constructor
418 
419 protected:
421  parameter_type type; // run-time specified basis type
422 
423  void Init(const parameter_type type_)
424  {
425  type = type_;
426  switch (type)
427  {
430  my_fe = new L2_FE_type(P, type);
431  my_fe_1d = (TP && dim != 1) ? new L2_SegmentElement(P, type) : NULL;
432  break;
434  my_fe = new L2Pos_FE_type(P);
435  my_fe_1d = (TP && dim != 1) ? new L2Pos_SegmentElement(P) : NULL;
436  break;
437  default:
438  MFEM_ABORT("invalid basis type");
439  }
440  }
441 
443  { Init(type); }
444 
446  {
447  const L2_FECollection *l2_fec =
448  dynamic_cast<const L2_FECollection *>(&fec);
449  MFEM_ASSERT(l2_fec, "invalid FiniteElementCollection");
450  Init(l2_fec->GetBasisType());
451  }
452 
453  ~L2_FiniteElement_base() { delete my_fe; delete my_fe_1d; }
454 
455 public:
456  template <typename real_t>
457  void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
458  {
459  mfem::CalcShapes(*my_fe, ir, B, Grad, NULL);
460  }
461  template <typename real_t>
462  void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
463  {
464  mfem::CalcShapes(dim == 1 ? *my_fe : *my_fe_1d, ir, B, Grad, NULL);
465  }
466  const Array<int> *GetDofMap() const { return NULL; }
467 };
468 
469 
470 template <Geometry::Type G, int P>
472 
473 
474 template <int P>
475 class L2_FiniteElement<Geometry::SEGMENT, P>
476  : public L2_FiniteElement_base<
477  Geometry::SEGMENT,P,L2_SegmentElement,L2Pos_SegmentElement,P+1,true>
478 {
479 protected:
482 public:
485  : base_class(type_) { }
487  : base_class(fec) { }
488 };
489 
490 
491 template <int P>
492 class L2_FiniteElement<Geometry::TRIANGLE, P>
493  : public L2_FiniteElement_base<Geometry::TRIANGLE,P,L2_TriangleElement,
494  L2Pos_TriangleElement,((P+1)*(P+2))/2,false>
495 {
496 protected:
498  L2Pos_TriangleElement,((P+1)*(P+2))/2,false> base_class;
499 public:
502  : base_class(type_) { }
504  : base_class(fec) { }
505 };
506 
507 
508 template <int P>
509 class L2_FiniteElement<Geometry::SQUARE, P>
510  : public L2_FiniteElement_base<Geometry::SQUARE,P,L2_QuadrilateralElement,
511  L2Pos_QuadrilateralElement,(P+1)*(P+1),true>
512 {
513 protected:
516 public:
519  : base_class(type_) { }
521  : base_class(fec) { }
522 };
523 
524 
525 template <int P>
526 class L2_FiniteElement<Geometry::TETRAHEDRON, P>
527  : public L2_FiniteElement_base<Geometry::TETRAHEDRON,P,L2_TetrahedronElement,
528  L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6,false>
529 {
530 protected:
532  L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6,false> base_class;
533 public:
536  : base_class(type_) { }
538  : base_class(fec) { }
539 };
540 
541 
542 template <int P>
543 class L2_FiniteElement<Geometry::CUBE, P>
544  : public L2_FiniteElement_base<Geometry::CUBE,P,L2_HexahedronElement,
545  L2Pos_HexahedronElement,(P+1)*(P+1)*(P+1),true>
546 {
547 protected:
549  L2Pos_HexahedronElement,(P+1)*(P+1)*(P+1),true> base_class;
550 public:
553  : base_class(type_) { }
555  : base_class(fec) { }
556 };
557 
558 } // namespace mfem
559 
560 #endif // MFEM_TEMPLATE_FINITE_ELEMENTS
const Array< int > * GetDofMap() const
Definition: tfe.hpp:327
void Init(const parameter_type type_)
Definition: tfe.hpp:291
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:228
Abstract class for Finite Elements.
Definition: fe.hpp:44
const Array< int > & GetDofMap() const
Definition: fe.hpp:1384
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:82
L2_FiniteElement(const parameter_type type_=L2_FECollection::GaussLegendre)
Definition: tfe.hpp:535
Class for integration rule.
Definition: intrules.hpp:83
L2_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:520
static const Geometry::Type geom
Definition: tfe.hpp:408
static const int dim
Definition: tfe.hpp:409
L2_FECollection::BasisType parameter_type
Definition: tfe.hpp:417
L2_FiniteElement_base< Geometry::SEGMENT, P, L2_SegmentElement, L2Pos_SegmentElement, P+1, true > base_class
Definition: tfe.hpp:481
const Array< int > * GetDofMap() const
Definition: tfe.hpp:271
L2_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:486
float real_t
Definition: mesh.cpp:4276
L2_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:554
H1_FECollection::BasisType parameter_type
Definition: tfe.hpp:95
H1_FiniteElement(const parameter_type type_=H1_FECollection::GaussLobatto)
Definition: tfe.hpp:184
const Geometry::Type geom
void Init(const parameter_type type_)
Definition: tfe.hpp:423
H1_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:313
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
void CalcShapeMatrix(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, const Array< int > *dof_map=NULL)
Definition: tfe.hpp:24
const Array< int > * GetDofMap() const
Definition: tfe.hpp:202
BasisType GetBasisType() const
Definition: fe_coll.hpp:170
const Array< int > & GetDofMap() const
Definition: fe.hpp:1362
L2_FiniteElement(const parameter_type type_=L2_FECollection::GaussLegendre)
Definition: tfe.hpp:484
H1_FECollection::BasisType parameter_type
Definition: tfe.hpp:286
H1_FECollection::BasisType parameter_type
Definition: tfe.hpp:161
H1_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:188
static const int degree
Definition: tfe.hpp:410
const Array< int > * GetDofMap() const
Definition: tfe.hpp:397
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:323
L2_FiniteElement(const parameter_type type_=L2_FECollection::GaussLegendre)
Definition: tfe.hpp:518
base_class::parameter_type parameter_type
Definition: tfe.hpp:500
L2_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:537
L2_FiniteElement(const parameter_type type_=L2_FECollection::GaussLegendre)
Definition: tfe.hpp:501
H1_FECollection::BasisType parameter_type
Definition: tfe.hpp:218
void Init(const parameter_type type_)
Definition: tfe.hpp:166
L2_FiniteElement_base(const FiniteElementCollection &fec)
Definition: tfe.hpp:445
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:231
int dim
Definition: ex3.cpp:47
const Array< int > & GetDofMap() const
Definition: fe.hpp:1343
void CalcGradTensor(const FiniteElement &fe, const IntegrationRule &ir, real_t *G, const Array< int > *dof_map=NULL)
Definition: tfe.hpp:45
H1_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:378
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
static const bool tensor_prod
Definition: tfe.hpp:413
L2_FiniteElement_base< Geometry::TETRAHEDRON, P, L2_TetrahedronElement, L2Pos_TetrahedronElement,((P+1)*(P+2)*(P+3))/6, false > base_class
Definition: tfe.hpp:532
H1_FiniteElement(const parameter_type type_=H1_FECollection::GaussLobatto)
Definition: tfe.hpp:309
base_class::parameter_type parameter_type
Definition: tfe.hpp:534
H1_FiniteElement(const parameter_type type_=H1_FECollection::GaussLobatto)
Definition: tfe.hpp:123
base_class::parameter_type parameter_type
Definition: tfe.hpp:517
base_class::parameter_type parameter_type
Definition: tfe.hpp:551
const Array< int > & GetDofMap() const
Definition: fe.hpp:1422
static const int dofs_1d
Definition: tfe.hpp:414
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:137
void Init(const parameter_type type_)
Definition: tfe.hpp:224
const Array< int > & GetDofMap() const
Definition: fe.hpp:1403
L2_FiniteElement_base(const parameter_type type)
Definition: tfe.hpp:442
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
Definition: tfe.hpp:457
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
H1_FiniteElement(const parameter_type type_=H1_FECollection::GaussLobatto)
Definition: tfe.hpp:374
parameter_type type
Definition: tfe.hpp:421
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:198
H1_FiniteElement(const parameter_type type_=H1_FECollection::GaussLobatto)
Definition: tfe.hpp:248
H1_FECollection::BasisType parameter_type
Definition: tfe.hpp:343
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:262
H1_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:252
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:267
H1_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:127
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:393
BasisType GetBasisType() const
Definition: fe_coll.hpp:105
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:142
const FiniteElement * my_fe_1d
Definition: tfe.hpp:420
Vector data type.
Definition: vector.hpp:33
const FiniteElement * my_fe
Definition: tfe.hpp:420
L2_FiniteElement(const parameter_type type_=L2_FECollection::GaussLegendre)
Definition: tfe.hpp:552
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:77
static const int dofs
Definition: tfe.hpp:411
base_class::parameter_type parameter_type
Definition: tfe.hpp:483
void CalcShapes(const IntegrationRule &ir, real_t *B, real_t *G) const
Definition: tfe.hpp:388
void Init(const parameter_type type_)
Definition: tfe.hpp:350
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
L2_FiniteElement_base< Geometry::SQUARE, P, L2_QuadrilateralElement, L2Pos_QuadrilateralElement,(P+1)*(P+1), true > base_class
Definition: tfe.hpp:515
void CalcShapes(const FiniteElement &fe, const IntegrationRule &ir, real_t *B, real_t *G, const Array< int > *dof_map)
Definition: tfe.hpp:70
void Calc1DShapes(const IntegrationRule &ir, real_t *B, real_t *Grad) const
Definition: tfe.hpp:462
const Array< int > * GetDofMap() const
Definition: tfe.hpp:146
void Init(const parameter_type type_)
Definition: tfe.hpp:101
L2_FiniteElement(const FiniteElementCollection &fec)
Definition: tfe.hpp:503
const Array< int > * GetDofMap() const
Definition: tfe.hpp:466
const Array< int > & GetDofMap() const
Definition: fe.hpp:1324
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:130