MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_coll.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_FE_COLLECTION
13 #define MFEM_FE_COLLECTION
14 
15 #include "../config/config.hpp"
16 #include "geom.hpp"
17 #include "fe.hpp"
18 
19 namespace mfem
20 {
21 
22 /// Possible basis types. Note that not all elements can use all BasisType(s).
23 class BasisType
24 {
25 public:
26  enum
27  {
28  GaussLegendre = 0, ///< Open type
29  GaussLobatto = 1, ///< Closed type
30  Positive = 2, ///< Bernstein polynomials
31  OpenUniform = 3, ///< Nodes: x_i = (i+1)/(n+1), i=0,...,n-1
32  ClosedUniform = 4, ///< Nodes: x_i = i/(n-1), i=0,...,n-1
33  OpenHalfUniform = 5 ///< Nodes: x_i = (i+1/2)/n, i=0,...,n-1
34  };
35  /** @brief If the input does not represents a valid BasisType, abort with an
36  error; otherwise return the input. */
37  static int Check(int b_type)
38  {
39  MFEM_VERIFY(0 <= b_type && b_type < 6, "unknown BasisType: " << b_type);
40  return b_type;
41  }
42  /** @brief Get the corresponding Quadrature1D constant, when that makes
43  sense; otherwise return Quadrature1D::Invalid. */
44  static int GetQuadrature1D(int b_type)
45  {
46  switch (b_type)
47  {
53  }
54  return Quadrature1D::Invalid;
55  }
56  /// Check and convert a BasisType constant to a string identifier.
57  static const char *Name(int b_type)
58  {
59  static const char *name[] =
60  {
61  "Gauss-Legendre", "Gauss-Lobatto", "Positive (Bernstein)",
62  "Open uniform", "Closed uniform", "Open half uniform"
63  };
64  return name[Check(b_type)];
65  }
66  /// Check and convert a BasisType constant to a char basis identifier.
67  static char GetChar(int b_type)
68  {
69  static const char ident[] = { 'g', 'G', 'P', 'u', 'U', 'o' };
70  return ident[Check(b_type)];
71  }
72  /// Convert char basis identifier to a BasisType constant.
73  static int GetType(char b_ident)
74  {
75  switch (b_ident)
76  {
77  case 'g': return GaussLegendre;
78  case 'G': return GaussLobatto;
79  case 'P': return Positive;
80  case 'u': return OpenUniform;
81  case 'U': return ClosedUniform;
82  case 'o': return OpenHalfUniform;
83  }
84  MFEM_ABORT("unknown BasisType identifier");
85  return -1;
86  }
87 };
88 
89 /** Collection of finite elements from the same family in multiple dimensions.
90  This class is used to match the degrees of freedom of a FiniteElementSpace
91  between elements, and to provide the finite element restriction from an
92  element to its boundary. */
94 {
95 protected:
96  template <Geometry::Type geom>
97  static inline void GetNVE(int &nv, int &ne);
98 
99  template <Geometry::Type geom, typename v_t>
100  static inline void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo,
101  const int edge_info);
102 
103  template <Geometry::Type geom, Geometry::Type f_geom,
104  typename v_t, typename e_t, typename eo_t>
105  static inline void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo,
106  int &nf, int &f, int &fg, int &fo,
107  const int face_info);
108 
109 public:
110  virtual const FiniteElement *
111  FiniteElementForGeometry(int GeomType) const = 0;
112 
113  virtual int DofForGeometry(int GeomType) const = 0;
114 
115  virtual int * DofOrderForOrientation(int GeomType, int Or) const = 0;
116 
117  virtual const char * Name() const { return "Undefined"; }
118 
119  int HasFaceDofs(int GeomType) const;
120 
122  int GeomType) const
123  {
124  return FiniteElementForGeometry(GeomType);
125  }
126 
128 
130 
131  /** @brief Factory method: return a newly allocated FiniteElementCollection
132  according to the given name. */
133  static FiniteElementCollection *New(const char *name);
134 
135  /** @brief Get the local dofs for a given sub-manifold.
136 
137  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex,
138  1D - edge, 2D - face) including those on its boundary. The local index of
139  the sub-manifold (inside Geom) and its orientation are given by the
140  parameter Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed
141  that 0 <= SDim <= Dim(Geom). */
142  void SubDofOrder(int Geom, int SDim, int Info, Array<int> &dofs) const;
143 };
144 
145 /// Arbitrary order H1-conforming (continuous) finite elements.
147 {
148 
149 protected:
150  int m_type;
151  char h1_name[32];
154  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
155 
156 public:
157  explicit H1_FECollection(const int p, const int dim = 3,
158  const int type = BasisType::GaussLobatto);
159 
160  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
161  { return H1_Elements[GeomType]; }
162  virtual int DofForGeometry(int GeomType) const
163  { return H1_dof[GeomType]; }
164  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
165  virtual const char *Name() const { return h1_name; }
167 
168  int GetBasisType() const { return m_type; }
169  /// Get the Cartesian to local H1 dof map
170  const int *GetDofMap(int GeomType) const;
171 
172  virtual ~H1_FECollection();
173 };
174 
175 /** Arbitrary order H1-conforming (continuous) finite elements with positive
176  basis functions. */
178 {
179 public:
180  explicit H1Pos_FECollection(const int p, const int dim = 3)
181  : H1_FECollection(p, dim, BasisType::Positive) { }
182 };
183 
184 /** Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the
185  interface between mesh elements (faces,edges,vertices); these are the trace
186  FEs of the H1-conforming FEs. */
188 {
189 public:
190  H1_Trace_FECollection(const int p, const int dim,
191  const int type = BasisType::GaussLobatto);
192 };
193 
194 /// Arbitrary order "L2-conforming" discontinuous finite elements.
196 {
197 private:
198  int m_type; // BasisType
199  char d_name[32];
202  int *SegDofOrd[2]; // for rotating segment dofs in 1D
203  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
204  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
205 
206 public:
207  L2_FECollection(const int p, const int dim,
208  const int type = BasisType::GaussLegendre,
209  const int map_type = FiniteElement::VALUE);
210 
211  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
212  { return L2_Elements[GeomType]; }
213  virtual int DofForGeometry(int GeomType) const
214  {
215  if (L2_Elements[GeomType])
216  {
217  return L2_Elements[GeomType]->GetDof();
218  }
219  return 0;
220  }
221  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
222  virtual const char *Name() const { return d_name; }
223 
225  int GeomType) const
226  {
227  return Tr_Elements[GeomType];
228  }
229 
230  int GetBasisType() const { return m_type; }
231 
232  virtual ~L2_FECollection();
233 };
234 
235 // Declare an alternative name for L2_FECollection = DG_FECollection
237 
238 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
240 {
241 protected:
242  int ob_type; // open BasisType
243  char rt_name[32];
246  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
247 
248  // Initialize only the face elements
249  void InitFaces(const int p, const int dim, const int map_type,
250  const bool signs);
251 
252  // Constructor used by the constructor of the RT_Trace_FECollection and
253  // DG_Interface_FECollection classes
254  RT_FECollection(const int p, const int dim, const int map_type,
255  const bool signs,
256  const int ob_type = BasisType::GaussLegendre);
257 
258 public:
259  RT_FECollection(const int p, const int dim,
260  const int cb_type = BasisType::GaussLobatto,
261  const int ob_type = BasisType::GaussLegendre);
262 
263  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
264  { return RT_Elements[GeomType]; }
265  virtual int DofForGeometry(int GeomType) const
266  { return RT_dof[GeomType]; }
267  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
268  virtual const char *Name() const { return rt_name; }
270 
271  virtual ~RT_FECollection();
272 };
273 
274 /** Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the
275  interface between mesh elements (faces); these are the normal trace FEs of
276  the H(div)-conforming FEs. */
278 {
279 public:
280  RT_Trace_FECollection(const int p, const int dim,
281  const int map_type = FiniteElement::INTEGRAL,
282  const int ob_type = BasisType::GaussLegendre);
283 };
284 
285 /** Arbitrary order discontinuous finite elements defined on the interface
286  between mesh elements (faces). The functions in this space are single-valued
287  on each face and are discontinuous across its boundary. */
289 {
290 public:
291  DG_Interface_FECollection(const int p, const int dim,
292  const int map_type = FiniteElement::VALUE,
293  const int ob_type = BasisType::GaussLegendre);
294 };
295 
296 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
298 {
299 protected:
300  char nd_name[32];
303  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
304 
305 public:
306  ND_FECollection(const int p, const int dim,
307  const int cb_type = BasisType::GaussLobatto,
308  const int ob_type = BasisType::GaussLegendre);
309 
310  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
311  { return ND_Elements[GeomType]; }
312  virtual int DofForGeometry(int GeomType) const
313  { return ND_dof[GeomType]; }
314  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
315  virtual const char *Name() const { return nd_name; }
317 
318  virtual ~ND_FECollection();
319 };
320 
321 /** Arbitrary order H(curl)-trace finite elements defined on the interface
322  between mesh elements (faces,edges); these are the tangential trace FEs of
323  the H(curl)-conforming FEs. */
325 {
326 public:
327  ND_Trace_FECollection(const int p, const int dim,
328  const int cb_type = BasisType::GaussLobatto,
329  const int ob_type = BasisType::GaussLegendre);
330 };
331 
332 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
334 {
335 private:
336  NURBS1DFiniteElement *SegmentFE;
337  NURBS2DFiniteElement *QuadrilateralFE;
338  NURBS3DFiniteElement *ParallelepipedFE;
339 
340  char name[16];
341 
342  void Allocate(int Order);
343  void Deallocate();
344 
345 public:
346  explicit NURBSFECollection(int Order) { Allocate(Order); }
347 
348  int GetOrder() const { return SegmentFE->GetOrder(); }
349 
350  /// Change the order of the collection
351  void UpdateOrder(int Order) { Deallocate(); Allocate(Order); }
352 
353  void Reset() const
354  {
355  SegmentFE->Reset();
356  QuadrilateralFE->Reset();
357  ParallelepipedFE->Reset();
358  }
359 
360  virtual const FiniteElement *
361  FiniteElementForGeometry(int GeomType) const;
362 
363  virtual int DofForGeometry(int GeomType) const;
364 
365  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
366 
367  virtual const char *Name() const { return name; }
368 
370 
371  virtual ~NURBSFECollection() { Deallocate(); }
372 };
373 
374 
375 /// Piecewise-(bi)linear continuous finite elements.
377 {
378 private:
379  const PointFiniteElement PointFE;
380  const Linear1DFiniteElement SegmentFE;
381  const Linear2DFiniteElement TriangleFE;
382  const BiLinear2DFiniteElement QuadrilateralFE;
383  const Linear3DFiniteElement TetrahedronFE;
384  const TriLinear3DFiniteElement ParallelepipedFE;
385 public:
387 
388  virtual const FiniteElement *
389  FiniteElementForGeometry(int GeomType) const;
390 
391  virtual int DofForGeometry(int GeomType) const;
392 
393  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
394 
395  virtual const char * Name() const { return "Linear"; }
396 };
397 
398 /// Piecewise-(bi)quadratic continuous finite elements.
400 {
401 private:
402  const PointFiniteElement PointFE;
403  const Quad1DFiniteElement SegmentFE;
404  const Quad2DFiniteElement TriangleFE;
405  const BiQuad2DFiniteElement QuadrilateralFE;
406  const Quadratic3DFiniteElement TetrahedronFE;
407  const LagrangeHexFiniteElement ParallelepipedFE;
408 
409 public:
410  QuadraticFECollection() : ParallelepipedFE(2) { }
411 
412  virtual const FiniteElement *
413  FiniteElementForGeometry(int GeomType) const;
414 
415  virtual int DofForGeometry(int GeomType) const;
416 
417  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
418 
419  virtual const char * Name() const { return "Quadratic"; }
420 };
421 
422 /// Version of QuadraticFECollection with positive basis functions.
424 {
425 private:
426  const QuadPos1DFiniteElement SegmentFE;
427  const BiQuadPos2DFiniteElement QuadrilateralFE;
428 
429 public:
431 
432  virtual const FiniteElement *
433  FiniteElementForGeometry(int GeomType) const;
434 
435  virtual int DofForGeometry(int GeomType) const;
436 
437  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
438 
439  virtual const char * Name() const { return "QuadraticPos"; }
440 };
441 
442 /// Piecewise-(bi)cubic continuous finite elements.
444 {
445 private:
446  const PointFiniteElement PointFE;
447  const Cubic1DFiniteElement SegmentFE;
448  const Cubic2DFiniteElement TriangleFE;
449  const BiCubic2DFiniteElement QuadrilateralFE;
450  const Cubic3DFiniteElement TetrahedronFE;
451  const LagrangeHexFiniteElement ParallelepipedFE;
452 
453 public:
454  CubicFECollection() : ParallelepipedFE(3) { }
455 
456  virtual const FiniteElement *
457  FiniteElementForGeometry(int GeomType) const;
458 
459  virtual int DofForGeometry(int GeomType) const;
460 
461  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
462 
463  virtual const char * Name() const { return "Cubic"; }
464 };
465 
466 /// Crouzeix-Raviart nonconforming elements in 2D.
468 {
469 private:
470  const P0SegmentFiniteElement SegmentFE;
471  const CrouzeixRaviartFiniteElement TriangleFE;
472  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
473 public:
474  CrouzeixRaviartFECollection() : SegmentFE(1) { }
475 
476  virtual const FiniteElement *
477  FiniteElementForGeometry(int GeomType) const;
478 
479  virtual int DofForGeometry(int GeomType) const;
480 
481  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
482 
483  virtual const char * Name() const { return "CrouzeixRaviart"; }
484 };
485 
486 /// Piecewise-linear nonconforming finite elements in 3D.
488 {
489 private:
490  const P0TriangleFiniteElement TriangleFE;
491  const P1TetNonConfFiniteElement TetrahedronFE;
492  const P0QuadFiniteElement QuadrilateralFE;
493  const RotTriLinearHexFiniteElement ParallelepipedFE;
494 
495 public:
497 
498  virtual const FiniteElement *
499  FiniteElementForGeometry(int GeomType) const;
500 
501  virtual int DofForGeometry(int GeomType) const;
502 
503  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
504 
505  virtual const char * Name() const { return "LinearNonConf3D"; }
506 };
507 
508 
509 /** First order Raviart-Thomas finite elements in 2D. This class is kept only
510  for backward compatibility, consider using RT_FECollection instead. */
512 {
513 private:
514  const P0SegmentFiniteElement SegmentFE; // normal component on edge
515  const RT0TriangleFiniteElement TriangleFE;
516  const RT0QuadFiniteElement QuadrilateralFE;
517 public:
518  RT0_2DFECollection() : SegmentFE(0) { }
519 
520  virtual const FiniteElement *
521  FiniteElementForGeometry(int GeomType) const;
522 
523  virtual int DofForGeometry(int GeomType) const;
524 
525  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
526 
527  virtual const char * Name() const { return "RT0_2D"; }
528 };
529 
530 /** Second order Raviart-Thomas finite elements in 2D. This class is kept only
531  for backward compatibility, consider using RT_FECollection instead. */
533 {
534 private:
535  const P1SegmentFiniteElement SegmentFE; // normal component on edge
536  const RT1TriangleFiniteElement TriangleFE;
537  const RT1QuadFiniteElement QuadrilateralFE;
538 public:
540 
541  virtual const FiniteElement *
542  FiniteElementForGeometry(int GeomType) const;
543 
544  virtual int DofForGeometry(int GeomType) const;
545 
546  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
547 
548  virtual const char * Name() const { return "RT1_2D"; }
549 };
550 
551 /** Third order Raviart-Thomas finite elements in 2D. This class is kept only
552  for backward compatibility, consider using RT_FECollection instead. */
554 {
555 private:
556  const P2SegmentFiniteElement SegmentFE; // normal component on edge
557  const RT2TriangleFiniteElement TriangleFE;
558  const RT2QuadFiniteElement QuadrilateralFE;
559 public:
561 
562  virtual const FiniteElement *
563  FiniteElementForGeometry(int GeomType) const;
564 
565  virtual int DofForGeometry(int GeomType) const;
566 
567  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
568 
569  virtual const char * Name() const { return "RT2_2D"; }
570 };
571 
572 /** Piecewise-constant discontinuous finite elements in 2D. This class is kept
573  only for backward compatibility, consider using L2_FECollection instead. */
575 {
576 private:
577  const P0TriangleFiniteElement TriangleFE;
578  const P0QuadFiniteElement QuadrilateralFE;
579 public:
581 
582  virtual const FiniteElement *
583  FiniteElementForGeometry(int GeomType) const;
584 
585  virtual int DofForGeometry(int GeomType) const;
586 
587  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
588 
589  virtual const char * Name() const { return "Const2D"; }
590 };
591 
592 /** Piecewise-linear discontinuous finite elements in 2D. This class is kept
593  only for backward compatibility, consider using L2_FECollection instead. */
595 {
596 private:
597  const Linear2DFiniteElement TriangleFE;
598  const BiLinear2DFiniteElement QuadrilateralFE;
599 
600 public:
602 
603  virtual const FiniteElement *
604  FiniteElementForGeometry(int GeomType) const;
605 
606  virtual int DofForGeometry(int GeomType) const;
607 
608  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
609 
610  virtual const char * Name() const { return "LinearDiscont2D"; }
611 };
612 
613 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
615 {
616 private:
617  // const CrouzeixRaviartFiniteElement TriangleFE;
618  const GaussLinear2DFiniteElement TriangleFE;
619  const GaussBiLinear2DFiniteElement QuadrilateralFE;
620 
621 public:
623 
624  virtual const FiniteElement *
625  FiniteElementForGeometry(int GeomType) const;
626 
627  virtual int DofForGeometry(int GeomType) const;
628 
629  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
630 
631  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
632 };
633 
634 /// Linear (P1) finite elements on quadrilaterals.
636 {
637 private:
638  const P1OnQuadFiniteElement QuadrilateralFE;
639 public:
641  virtual const FiniteElement *
642  FiniteElementForGeometry(int GeomType) const;
643  virtual int DofForGeometry(int GeomType) const;
644  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
645  virtual const char * Name() const { return "P1OnQuad"; }
646 };
647 
648 /** Piecewise-quadratic discontinuous finite elements in 2D. This class is kept
649  only for backward compatibility, consider using L2_FECollection instead. */
651 {
652 private:
653  const Quad2DFiniteElement TriangleFE;
654  const BiQuad2DFiniteElement QuadrilateralFE;
655 
656 public:
658 
659  virtual const FiniteElement *
660  FiniteElementForGeometry(int GeomType) const;
661 
662  virtual int DofForGeometry(int GeomType) const;
663 
664  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
665 
666  virtual const char * Name() const { return "QuadraticDiscont2D"; }
667 };
668 
669 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
671 {
672 private:
673  const BiQuadPos2DFiniteElement QuadrilateralFE;
674 
675 public:
677  virtual const FiniteElement *
678  FiniteElementForGeometry(int GeomType) const;
679  virtual int DofForGeometry(int GeomType) const;
680  virtual int * DofOrderForOrientation(int GeomType, int Or) const
681  { return NULL; }
682  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
683 };
684 
685 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
687 {
688 private:
689  // const Quad2DFiniteElement TriangleFE;
690  const GaussQuad2DFiniteElement TriangleFE;
691  const GaussBiQuad2DFiniteElement QuadrilateralFE;
692 
693 public:
695 
696  virtual const FiniteElement *
697  FiniteElementForGeometry(int GeomType) const;
698 
699  virtual int DofForGeometry(int GeomType) const;
700 
701  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
702 
703  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
704 };
705 
706 /** Piecewise-cubic discontinuous finite elements in 2D. This class is kept
707  only for backward compatibility, consider using L2_FECollection instead. */
709 {
710 private:
711  const Cubic2DFiniteElement TriangleFE;
712  const BiCubic2DFiniteElement QuadrilateralFE;
713 
714 public:
716 
717  virtual const FiniteElement *
718  FiniteElementForGeometry(int GeomType) const;
719 
720  virtual int DofForGeometry(int GeomType) const;
721 
722  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
723 
724  virtual const char * Name() const { return "CubicDiscont2D"; }
725 };
726 
727 /** Piecewise-constant discontinuous finite elements in 3D. This class is kept
728  only for backward compatibility, consider using L2_FECollection instead. */
730 {
731 private:
732  const P0TetFiniteElement TetrahedronFE;
733  const P0HexFiniteElement ParallelepipedFE;
734 
735 public:
737 
738  virtual const FiniteElement *
739  FiniteElementForGeometry(int GeomType) const;
740 
741  virtual int DofForGeometry(int GeomType) const;
742 
743  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
744 
745  virtual const char * Name() const { return "Const3D"; }
746 };
747 
748 /** Piecewise-linear discontinuous finite elements in 3D. This class is kept
749  only for backward compatibility, consider using L2_FECollection instead. */
751 {
752 private:
753  const Linear3DFiniteElement TetrahedronFE;
754  const TriLinear3DFiniteElement ParallelepipedFE;
755 
756 public:
758 
759  virtual const FiniteElement *
760  FiniteElementForGeometry(int GeomType) const;
761 
762  virtual int DofForGeometry(int GeomType) const;
763 
764  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
765 
766  virtual const char * Name() const { return "LinearDiscont3D"; }
767 };
768 
769 /** Piecewise-quadratic discontinuous finite elements in 3D. This class is kept
770  only for backward compatibility, consider using L2_FECollection instead. */
772 {
773 private:
774  const Quadratic3DFiniteElement TetrahedronFE;
775  const LagrangeHexFiniteElement ParallelepipedFE;
776 
777 public:
778  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
779 
780  virtual const FiniteElement *
781  FiniteElementForGeometry(int GeomType) const;
782 
783  virtual int DofForGeometry(int GeomType) const;
784 
785  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
786 
787  virtual const char * Name() const { return "QuadraticDiscont3D"; }
788 };
789 
790 /// Finite element collection on a macro-element.
792 {
793 private:
794  const PointFiniteElement PointFE;
795  const RefinedLinear1DFiniteElement SegmentFE;
796  const RefinedLinear2DFiniteElement TriangleFE;
797  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
798  const RefinedLinear3DFiniteElement TetrahedronFE;
799  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
800 
801 public:
803 
804  virtual const FiniteElement *
805  FiniteElementForGeometry(int GeomType) const;
806 
807  virtual int DofForGeometry(int GeomType) const;
808 
809  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
810 
811  virtual const char * Name() const { return "RefinedLinear"; }
812 };
813 
814 /** Lowest order Nedelec finite elements in 3D. This class is kept only for
815  backward compatibility, consider using the new ND_FECollection instead. */
817 {
818 private:
819  const Nedelec1HexFiniteElement HexahedronFE;
820  const Nedelec1TetFiniteElement TetrahedronFE;
821 
822 public:
824 
825  virtual const FiniteElement *
826  FiniteElementForGeometry(int GeomType) const;
827 
828  virtual int DofForGeometry(int GeomType) const;
829 
830  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
831 
832  virtual const char * Name() const { return "ND1_3D"; }
833 };
834 
835 /** First order Raviart-Thomas finite elements in 3D. This class is kept only
836  for backward compatibility, consider using RT_FECollection instead. */
838 {
839 private:
840  const P0TriangleFiniteElement TriangleFE;
841  const P0QuadFiniteElement QuadrilateralFE;
842  const RT0HexFiniteElement HexahedronFE;
843  const RT0TetFiniteElement TetrahedronFE;
844 public:
846 
847  virtual const FiniteElement *
848  FiniteElementForGeometry(int GeomType) const;
849 
850  virtual int DofForGeometry(int GeomType) const;
851 
852  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
853 
854  virtual const char * Name() const { return "RT0_3D"; }
855 };
856 
857 /** Second order Raviart-Thomas finite elements in 3D. This class is kept only
858  for backward compatibility, consider using RT_FECollection instead. */
860 {
861 private:
862  const Linear2DFiniteElement TriangleFE;
863  const BiLinear2DFiniteElement QuadrilateralFE;
864  const RT1HexFiniteElement HexahedronFE;
865 public:
867 
868  virtual const FiniteElement *
869  FiniteElementForGeometry(int GeomType) const;
870 
871  virtual int DofForGeometry(int GeomType) const;
872 
873  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
874 
875  virtual const char * Name() const { return "RT1_3D"; }
876 };
877 
878 /// Discontinuous collection defined locally by a given finite element.
880 {
881 private:
882  char d_name[32];
883  int GeomType;
884  FiniteElement *Local_Element;
885 
886 public:
887  Local_FECollection(const char *fe_name);
888 
889  virtual const FiniteElement *FiniteElementForGeometry(int _GeomType) const
890  { return (GeomType == _GeomType) ? Local_Element : NULL; }
891  virtual int DofForGeometry(int _GeomType) const
892  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
893  virtual int *DofOrderForOrientation(int GeomType, int Or) const
894  { return NULL; }
895  virtual const char *Name() const { return d_name; }
896 
897  virtual ~Local_FECollection() { delete Local_Element; }
898 };
899 
900 }
901 
902 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:47
int GetOrder() const
Definition: fe_coll.hpp:348
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:333
virtual const char * Name() const
Definition: fe_coll.hpp:463
virtual const char * Name() const
Definition: fe_coll.hpp:787
virtual const char * Name() const
Definition: fe_coll.hpp:703
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1243
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:686
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:245
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:301
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1348
void UpdateOrder(int Order)
Change the order of the collection.
Definition: fe_coll.hpp:351
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1047
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:1057
virtual const char * Name() const
Definition: fe_coll.hpp:811
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_coll.hpp:37
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:38
virtual const char * Name() const
Definition: fe_coll.hpp:745
void SubDofOrder(int Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:353
RT_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:1938
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:615
virtual const char * Name() const
Definition: fe_coll.hpp:610
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2468
static const int NumGeom
Definition: geom.hpp:34
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2462
virtual const char * Name() const
Definition: fe_coll.hpp:419
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1411
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1171
H1_FECollection(const int p, const int dim=3, const int type=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1461
const Geometry::Type geom
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1089
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1865
virtual const char * Name() const
Definition: fe_coll.hpp:645
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:704
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1332
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:934
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:992
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:970
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:900
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:125
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:376
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1439
void Reset() const
Definition: fe.hpp:2298
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1271
Class for quadratic FE on triangle.
Definition: fe.hpp:631
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:810
Class for quadratic FE on interval.
Definition: fe.hpp:602
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1278
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
Definition: fe_coll.hpp:44
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1377
virtual const char * Name() const
Definition: fe_coll.hpp:268
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2333
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:265
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1123
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1109
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:892
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:1951
H1_Trace_FECollection(const int p, const int dim, const int type=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1685
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2086
virtual const char * Name() const
Definition: fe_coll.hpp:367
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1310
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2134
virtual const char * Name() const
Definition: fe_coll.hpp:527
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe_coll.hpp:57
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1207
virtual const char * Name() const
Definition: fe_coll.hpp:724
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:837
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1060
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Definition: fe_coll.hpp:73
Finite element collection on a macro-element.
Definition: fe_coll.hpp:791
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:783
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:539
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:562
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1186
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe_coll.hpp:23
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:310
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:824
Bernstein polynomials.
Definition: fe_coll.hpp:30
virtual const char * Name() const
Definition: fe_coll.hpp:315
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:962
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2102
virtual int DofForGeometry(int _GeomType) const
Definition: fe_coll.hpp:891
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:878
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2069
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1393
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:443
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:926
virtual const char * Name() const
Definition: fe_coll.hpp:548
int dim
Definition: ex3.cpp:47
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:577
Class for refined linear FE on interval.
Definition: fe.hpp:1133
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1424
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:153
Class for refined linear FE on triangle.
Definition: fe.hpp:1153
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1173
aka open Newton-Cotes
Definition: intrules.hpp:266
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:244
Class for constant FE on triangle.
Definition: fe.hpp:770
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:278
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:628
int GetBasisType() const
Definition: fe_coll.hpp:230
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:912
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:555
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:979
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1183
int GetBasisType() const
Definition: fe_coll.hpp:168
ND_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2365
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2154
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1206
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:879
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:603
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:670
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2353
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1199
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:741
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:825
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:893
For scalar fields; preserves point values.
Definition: fe.hpp:85
const int * GetDofMap(int GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1647
Nodes: x_i = (i+1/2)/n, i=0,...,n-1.
Definition: fe_coll.hpp:33
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1012
For scalar fields; preserves volume integrals.
Definition: fe.hpp:86
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:866
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:423
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1131
virtual ~Local_FECollection()
Definition: fe_coll.hpp:897
virtual const char * Name() const
Definition: fe_coll.hpp:631
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1294
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1673
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:523
Class for linear FE on tetrahedron.
Definition: fe.hpp:800
Class for linear FE on interval.
Definition: fe.hpp:500
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:873
Class for linear FE on triangle.
Definition: fe.hpp:520
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1363
Class for quadratic FE on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:655
virtual const char * Name() const
Definition: fe_coll.hpp:165
virtual const char * Name() const
Definition: fe_coll.hpp:766
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1628
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:152
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:487
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:670
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:467
NURBSFECollection(int Order)
Definition: fe_coll.hpp:346
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:239
L2_FECollection DG_FECollection
Definition: fe_coll.hpp:236
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:731
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:180
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:717
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:288
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:796
virtual const char * Name() const
Definition: fe_coll.hpp:483
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:213
virtual const char * Name() const
Definition: fe_coll.hpp:395
virtual const char * Name() const
Definition: fe_coll.hpp:569
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:224
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1219
virtual const char * Name() const
Definition: fe_coll.hpp:505
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1026
virtual const char * Name() const
Definition: fe_coll.hpp:832
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:122
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:671
Class for tri-linear FE on cube.
Definition: fe.hpp:838
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_coll.hpp:32
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1034
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:594
Nodes: x_i = (i+1)/(n+1), i=0,...,n-1.
Definition: fe_coll.hpp:31
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:947
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:263
virtual const char * Name() const
Definition: fe_coll.hpp:666
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:542
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1255
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1097
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1074
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:614
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:2443
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:222
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:44
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:121
Class for linear FE on triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:566
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
Definition: fe_coll.hpp:67
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1235
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1320
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:312
virtual const char * Name() const
Definition: fe_coll.hpp:117
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:578
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:160
aka closed Newton-Cotes
Definition: intrules.hpp:267
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:710
virtual const char * Name() const
Definition: fe_coll.hpp:589
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:851
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2310
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1611
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:768
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:162
virtual const char * Name() const
Definition: fe_coll.hpp:854
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2114
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:1112
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:399
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:211
aka &quot;open half&quot; Newton-Cotes
Definition: intrules.hpp:268
int HasFaceDofs(int GeomType) const
Definition: fe_coll.cpp:25
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:297
void Reset() const
Definition: fe_coll.hpp:353
Class for cubic FE on tetrahedron.
Definition: fe.hpp:757
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, int &fg, int &fo, const int face_info)
Definition: fe_coll.cpp:308
virtual const char * Name() const
Definition: fe_coll.hpp:895
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1161
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1850
L2_FECollection(const int p, const int dim, const int type=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1705
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
virtual const char * Name() const
Definition: fe_coll.hpp:875
virtual const char * Name() const
Definition: fe_coll.hpp:682
virtual const FiniteElement * FiniteElementForGeometry(int _GeomType) const
Definition: fe_coll.hpp:889
virtual int DofForGeometry(int GeomType) const =0
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:680
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:2456
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2384
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1045
virtual const char * Name() const
Definition: fe_coll.hpp:439
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1145
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:635
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:654
virtual ~NURBSFECollection()
Definition: fe_coll.hpp:371
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:302
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:754
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:861
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:638
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1000