MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_coll.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 #include "fem.hpp"
13 #include <cstdlib>
14 #include <cstring>
15 #include <cstdio>
16 #ifdef _WIN32
17 #define snprintf _snprintf_s
18 #endif
19 
20 namespace mfem
21 {
22 
23 using namespace std;
24 
26 {
27  switch (geom)
28  {
30  return GetNumDof(Geometry::TRIANGLE, p);
31  case Geometry::CUBE:
32  return GetNumDof(Geometry::SQUARE, p);
33  case Geometry::PRISM:
34  return max(GetNumDof(Geometry::TRIANGLE, p),
35  GetNumDof(Geometry::SQUARE, p));
36  case Geometry::PYRAMID:
37  return max(GetNumDof(Geometry::TRIANGLE, p),
38  GetNumDof(Geometry::SQUARE, p));
39  default:
40  MFEM_ABORT("unknown geometry type");
41  }
42  return 0;
43 }
44 
46 {
47  MFEM_ABORT("this method is not implemented in this derived class!");
48  return NULL;
49 }
50 
52 {
53  FiniteElementCollection *fec = NULL;
54 
55  if (!strcmp(name, "Linear"))
56  {
57  fec = new LinearFECollection;
58  }
59  else if (!strcmp(name, "Quadratic"))
60  {
61  fec = new QuadraticFECollection;
62  }
63  else if (!strcmp(name, "QuadraticPos"))
64  {
65  fec = new QuadraticPosFECollection;
66  }
67  else if (!strcmp(name, "Cubic"))
68  {
69  fec = new CubicFECollection;
70  }
71  else if (!strcmp(name, "Const3D"))
72  {
73  fec = new Const3DFECollection;
74  }
75  else if (!strcmp(name, "Const2D"))
76  {
77  fec = new Const2DFECollection;
78  }
79  else if (!strcmp(name, "LinearDiscont2D"))
80  {
82  }
83  else if (!strcmp(name, "GaussLinearDiscont2D"))
84  {
86  }
87  else if (!strcmp(name, "P1OnQuad"))
88  {
89  fec = new P1OnQuadFECollection;
90  }
91  else if (!strcmp(name, "QuadraticDiscont2D"))
92  {
94  }
95  else if (!strcmp(name, "QuadraticPosDiscont2D"))
96  {
98  }
99  else if (!strcmp(name, "GaussQuadraticDiscont2D"))
100  {
102  }
103  else if (!strcmp(name, "CubicDiscont2D"))
104  {
105  fec = new CubicDiscont2DFECollection;
106  }
107  else if (!strcmp(name, "LinearDiscont3D"))
108  {
109  fec = new LinearDiscont3DFECollection;
110  }
111  else if (!strcmp(name, "QuadraticDiscont3D"))
112  {
114  }
115  else if (!strcmp(name, "LinearNonConf3D"))
116  {
117  fec = new LinearNonConf3DFECollection;
118  }
119  else if (!strcmp(name, "CrouzeixRaviart"))
120  {
121  fec = new CrouzeixRaviartFECollection;
122  }
123  else if (!strcmp(name, "ND1_3D"))
124  {
125  fec = new ND1_3DFECollection;
126  }
127  else if (!strcmp(name, "RT0_2D"))
128  {
129  fec = new RT0_2DFECollection;
130  }
131  else if (!strcmp(name, "RT1_2D"))
132  {
133  fec = new RT1_2DFECollection;
134  }
135  else if (!strcmp(name, "RT2_2D"))
136  {
137  fec = new RT2_2DFECollection;
138  }
139  else if (!strcmp(name, "RT0_3D"))
140  {
141  fec = new RT0_3DFECollection;
142  }
143  else if (!strcmp(name, "RT1_3D"))
144  {
145  fec = new RT1_3DFECollection;
146  }
147  else if (!strncmp(name, "H1_Trace_", 9))
148  {
149  fec = new H1_Trace_FECollection(atoi(name + 13), atoi(name + 9));
150  }
151  else if (!strncmp(name, "H1_Trace@", 9))
152  {
153  fec = new H1_Trace_FECollection(atoi(name + 15), atoi(name + 11),
154  BasisType::GetType(name[9]));
155  }
156  else if (!strncmp(name, "H1_", 3))
157  {
158  fec = new H1_FECollection(atoi(name + 7), atoi(name + 3));
159  }
160  else if (!strncmp(name, "H1Pos_Trace_", 12))
161  {
162  fec = new H1_Trace_FECollection(atoi(name + 16), atoi(name + 12),
164  }
165  else if (!strncmp(name, "H1Pos_", 6))
166  {
167  fec = new H1Pos_FECollection(atoi(name + 10), atoi(name + 6));
168  }
169  else if (!strncmp(name, "H1Ser_", 6))
170  {
171  fec = new H1Ser_FECollection(atoi(name + 10), atoi(name + 6));
172  }
173  else if (!strncmp(name, "H1@", 3))
174  {
175  fec = new H1_FECollection(atoi(name + 9), atoi(name + 5),
176  BasisType::GetType(name[3]));
177  }
178  else if (!strncmp(name, "L2_T", 4))
179  fec = new L2_FECollection(atoi(name + 10), atoi(name + 6),
180  atoi(name + 4));
181  else if (!strncmp(name, "L2_", 3))
182  {
183  fec = new L2_FECollection(atoi(name + 7), atoi(name + 3));
184  }
185  else if (!strncmp(name, "L2Int_T", 7))
186  {
187  fec = new L2_FECollection(atoi(name + 13), atoi(name + 9),
188  atoi(name + 7), FiniteElement::INTEGRAL);
189  }
190  else if (!strncmp(name, "L2Int_", 6))
191  {
192  fec = new L2_FECollection(atoi(name + 10), atoi(name + 6),
195  }
196  else if (!strncmp(name, "RT_Trace_", 9))
197  {
198  fec = new RT_Trace_FECollection(atoi(name + 13), atoi(name + 9));
199  }
200  else if (!strncmp(name, "RT_ValTrace_", 12))
201  {
202  fec = new RT_Trace_FECollection(atoi(name + 16), atoi(name + 12),
204  }
205  else if (!strncmp(name, "RT_Trace@", 9))
206  {
207  fec = new RT_Trace_FECollection(atoi(name + 15), atoi(name + 11),
209  BasisType::GetType(name[9]));
210  }
211  else if (!strncmp(name, "RT_ValTrace@", 12))
212  {
213  fec = new RT_Trace_FECollection(atoi(name + 18), atoi(name + 14),
215  BasisType::GetType(name[12]));
216  }
217  else if (!strncmp(name, "DG_Iface_", 9))
218  {
219  fec = new DG_Interface_FECollection(atoi(name + 13), atoi(name + 9));
220  }
221  else if (!strncmp(name, "DG_Iface@", 9))
222  {
223  fec = new DG_Interface_FECollection(atoi(name + 15), atoi(name + 11),
225  BasisType::GetType(name[9]));
226  }
227  else if (!strncmp(name, "DG_IntIface_", 12))
228  {
229  fec = new DG_Interface_FECollection(atoi(name + 16), atoi(name + 12),
231  }
232  else if (!strncmp(name, "DG_IntIface@", 12))
233  {
234  fec = new DG_Interface_FECollection(atoi(name + 18), atoi(name + 14),
236  BasisType::GetType(name[12]));
237  }
238  else if (!strncmp(name, "RT_", 3))
239  {
240  fec = new RT_FECollection(atoi(name + 7), atoi(name + 3));
241  }
242  else if (!strncmp(name, "RT@", 3))
243  {
244  fec = new RT_FECollection(atoi(name + 10), atoi(name + 6),
245  BasisType::GetType(name[3]),
246  BasisType::GetType(name[4]));
247  }
248  else if (!strncmp(name, "ND_Trace_", 9))
249  {
250  fec = new ND_Trace_FECollection(atoi(name + 13), atoi(name + 9));
251  }
252  else if (!strncmp(name, "ND_Trace@", 9))
253  {
254  fec = new ND_Trace_FECollection(atoi(name + 16), atoi(name + 12),
255  BasisType::GetType(name[9]),
256  BasisType::GetType(name[10]));
257  }
258  else if (!strncmp(name, "ND_", 3))
259  {
260  fec = new ND_FECollection(atoi(name + 7), atoi(name + 3));
261  }
262  else if (!strncmp(name, "ND@", 3))
263  {
264  fec = new ND_FECollection(atoi(name + 10), atoi(name + 6),
265  BasisType::GetType(name[3]),
266  BasisType::GetType(name[4]));
267  }
268  else if (!strncmp(name, "Local_", 6))
269  {
270  fec = new Local_FECollection(name + 6);
271  }
272  else if (!strncmp(name, "NURBS", 5))
273  {
274  if (name[5] != '\0')
275  {
276  // "NURBS" + "number" --> fixed order nurbs collection
277  fec = new NURBSFECollection(atoi(name + 5));
278  }
279  else
280  {
281  // "NURBS" --> variable order nurbs collection
282  fec = new NURBSFECollection();
283  }
284  }
285  else
286  {
287  MFEM_ABORT("unknown FiniteElementCollection: " << name);
288  }
289  MFEM_VERIFY(!strcmp(fec->Name(), name), "input name: \"" << name
290  << "\" does not match the created collection name: \""
291  << fec->Name() << '"');
292 
293  return fec;
294 }
295 
297 {
298  // default implementation for collections that don't care about variable p
299  MFEM_ABORT("Collection " << Name() << " does not support variable orders.");
300  (void) p;
301  return NULL;
302 }
303 
305 {
306  if (p >= var_orders.Size())
307  {
308  var_orders.SetSize(p+1, NULL);
309  }
310  var_orders[p] = Clone(p);
311 }
312 
314 {
315  for (int i = 0; i < var_orders.Size(); i++)
316  {
317  delete var_orders[i];
318  }
319 }
320 
321 template <Geometry::Type geom>
322 inline void FiniteElementCollection::GetNVE(int &nv, int &ne)
323 {
324  typedef typename Geometry::Constants<geom> g_consts;
325 
326  nv = g_consts::NumVert;
327  ne = g_consts::NumEdges;
328 }
329 
330 template <Geometry::Type geom, typename v_t>
331 inline void FiniteElementCollection::
332 GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
333 {
334  typedef typename Geometry::Constants<Geometry::SEGMENT> e_consts;
335  typedef typename Geometry::Constants<geom> g_consts;
336 
337  nv = e_consts::NumVert;
338  ne = 1;
339  e = edge_info/64;
340  eo = edge_info%64;
341  MFEM_ASSERT(0 <= e && e < g_consts::NumEdges, "");
342  MFEM_ASSERT(0 <= eo && eo < e_consts::NumOrient, "");
343  v[0] = e_consts::Orient[eo][0];
344  v[1] = e_consts::Orient[eo][1];
345  v[0] = g_consts::Edges[e][v[0]];
346  v[1] = g_consts::Edges[e][v[1]];
347 }
348 
349 template <Geometry::Type geom, Geometry::Type f_geom,
350  typename v_t, typename e_t, typename eo_t>
351 inline void FiniteElementCollection::
352 GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo,
353  int &nf, int &f, Geometry::Type &fg, int &fo, const int face_info)
354 {
355  typedef typename Geometry::Constants< geom> g_consts;
356  typedef typename Geometry::Constants<f_geom> f_consts;
357 
358  nv = f_consts::NumVert;
359  nf = 1;
360  f = face_info/64;
361  fg = f_geom;
362  fo = face_info%64;
363  MFEM_ASSERT(0 <= f && f < g_consts::NumFaces, "");
364  MFEM_ASSERT(0 <= fo && fo < f_consts::NumOrient, "");
365  for (int i = 0; i < f_consts::NumVert; i++)
366  {
367  v[i] = f_consts::Orient[fo][i];
368  v[i] = g_consts::FaceVert[f][v[i]];
369  }
370  ne = f_consts::NumEdges;
371  for (int i = 0; i < f_consts::NumEdges; i++)
372  {
373  int v0 = v[f_consts::Edges[i][0]];
374  int v1 = v[f_consts::Edges[i][1]];
375  int eor = 0;
376  if (v0 > v1) { swap(v0, v1); eor = 1; }
377  for (int j = g_consts::VertToVert::I[v0]; true; j++)
378  {
379  MFEM_ASSERT(j < g_consts::VertToVert::I[v0+1],
380  "internal error, edge not found");
381  if (v1 == g_consts::VertToVert::J[j][0])
382  {
383  int en = g_consts::VertToVert::J[j][1];
384  if (en < 0)
385  {
386  en = -1-en;
387  eor = 1-eor;
388  }
389  e[i] = en;
390  eo[i] = eor;
391  break;
392  }
393  }
394  }
395 }
396 
398  int Info,
399  Array<int> &dofs) const
400 {
401  // Info = 64 * SubIndex + SubOrientation
402  MFEM_ASSERT(0 <= Geom && Geom < Geometry::NumGeom,
403  "invalid Geom = " << Geom);
404  MFEM_ASSERT(0 <= SDim && SDim <= Geometry::Dimension[Geom],
405  "invalid SDim = " << SDim <<
406  " for Geom = " << Geometry::Name[Geom]);
407 
408  const int nvd = DofForGeometry(Geometry::POINT);
409  if (SDim == 0) // vertex
410  {
411  const int off = nvd*(Info/64);
412  dofs.SetSize(nvd);
413  for (int i = 0; i < nvd; i++)
414  {
415  dofs[i] = off + i;
416  }
417  }
418  else
419  {
420  int v[4], e[4], eo[4], f[1], fo[1];
421  int av = 0, nv = 0, ae = 0, ne = 0, nf = 0;
422  Geometry::Type fg[1];
423 
424  switch (Geom)
425  {
426  case Geometry::SEGMENT:
427  {
428  GetNVE<Geometry::SEGMENT>(av, ae);
429  GetEdge<Geometry::SEGMENT>(nv, v, ne, e[0], eo[0], Info);
430  break;
431  }
432 
433  case Geometry::TRIANGLE:
434  {
435  GetNVE<Geometry::TRIANGLE>(av, ae);
436  switch (SDim)
437  {
438  case 1:
439  GetEdge<Geometry::TRIANGLE>(nv, v, ne, e[0], eo[0], Info);
440  break;
441  case 2:
442  GetFace<Geometry::TRIANGLE,Geometry::TRIANGLE>(
443  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
444  break;
445  default:
446  goto not_supp;
447  }
448  break;
449  }
450 
451  case Geometry::SQUARE:
452  {
453  GetNVE<Geometry::SQUARE>(av, ae);
454  switch (SDim)
455  {
456  case 1:
457  GetEdge<Geometry::SQUARE>(nv, v, ne, e[0], eo[0], Info);
458  break;
459  case 2:
460  GetFace<Geometry::SQUARE,Geometry::SQUARE>(
461  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
462  break;
463  default:
464  goto not_supp;
465  }
466  break;
467  }
468 
470  {
471  GetNVE<Geometry::TETRAHEDRON>(av, ae);
472  switch (SDim)
473  {
474  case 1:
475  GetEdge<Geometry::TETRAHEDRON>(nv, v, ne, e[0], eo[0], Info);
476  break;
477  case 2:
478  GetFace<Geometry::TETRAHEDRON,Geometry::TRIANGLE>(
479  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
480  break;
481  default:
482  goto not_supp;
483  }
484  break;
485  }
486 
487  case Geometry::CUBE:
488  {
489  GetNVE<Geometry::CUBE>(av, ae);
490  switch (SDim)
491  {
492  case 1:
493  GetEdge<Geometry::CUBE>(nv, v, ne, e[0], eo[0], Info);
494  break;
495  case 2:
496  GetFace<Geometry::CUBE,Geometry::SQUARE>(
497  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
498  break;
499  default:
500  goto not_supp;
501  }
502  break;
503  }
504 
505  default:
506  MFEM_ABORT("invalid Geom = " << Geom);
507  }
508 
509  int ned = (ne > 0) ? DofForGeometry(Geometry::SEGMENT) : 0;
510 
511  // add vertex dofs
512  dofs.SetSize(nv*nvd+ne*ned);
513  for (int i = 0; i < nv; i++)
514  {
515  for (int j = 0; j < nvd; j++)
516  {
517  dofs[i*nvd+j] = v[i]*nvd+j;
518  }
519  }
520  int l_off = nv*nvd, g_off = av*nvd;
521 
522  // add edge dofs
523  if (ned > 0)
524  {
525  for (int i = 0; i < ne; i++)
526  {
527  const int *ed = DofOrderForOrientation(Geometry::SEGMENT,
528  eo[i] ? -1 : 1);
529  for (int j = 0; j < ned; j++)
530  {
531  dofs[l_off+i*ned+j] =
532  ed[j] >= 0 ?
533  g_off+e[i]*ned+ed[j] :
534  -1-(g_off+e[i]*ned+(-1-ed[j]));
535  }
536  }
537  l_off += ne*ned;
538  g_off += ae*ned;
539  }
540 
541  // add face dofs
542  if (nf > 0)
543  {
544  const int nfd = DofForGeometry(fg[0]); // assume same face geometry
545  dofs.SetSize(dofs.Size()+nf*nfd);
546  for (int i = 0; i < nf; i++)
547  {
548  const int *fd = DofOrderForOrientation(fg[i], fo[i]);
549  for (int j = 0; j < nfd; j++)
550  {
551  dofs[l_off+i*nfd+j] =
552  fd[j] >= 0 ?
553  g_off+f[i]*nfd+fd[j] :
554  -1-(g_off+f[i]*nfd+(-1-fd[j]));
555  }
556  }
557  }
558 
559  // add volume dofs ...
560  }
561  return;
562 
563 not_supp:
564  MFEM_ABORT("Geom = " << Geometry::Name[Geom] <<
565  ", SDim = " << SDim << " is not supported");
566 }
567 
568 const FiniteElement *
570 {
571  switch (GeomType)
572  {
573  case Geometry::POINT: return &PointFE;
574  case Geometry::SEGMENT: return &SegmentFE;
575  case Geometry::TRIANGLE: return &TriangleFE;
576  case Geometry::SQUARE: return &QuadrilateralFE;
577  case Geometry::TETRAHEDRON: return &TetrahedronFE;
578  case Geometry::CUBE: return &ParallelepipedFE;
579  case Geometry::PRISM: return &WedgeFE;
580  case Geometry::PYRAMID: return &PyramidFE;
581  default:
582  mfem_error ("LinearFECollection: unknown geometry type.");
583  }
584  return &SegmentFE; // Make some compilers happy
585 }
586 
588 {
589  switch (GeomType)
590  {
591  case Geometry::POINT: return 1;
592  case Geometry::SEGMENT: return 0;
593  case Geometry::TRIANGLE: return 0;
594  case Geometry::SQUARE: return 0;
595  case Geometry::TETRAHEDRON: return 0;
596  case Geometry::CUBE: return 0;
597  case Geometry::PRISM: return 0;
598  case Geometry::PYRAMID: return 0;
599  default:
600  mfem_error ("LinearFECollection: unknown geometry type.");
601  }
602  return 0; // Make some compilers happy
603 }
604 
606  int Or) const
607 {
608  return NULL;
609 }
610 
611 
612 const FiniteElement *
614 {
615  switch (GeomType)
616  {
617  case Geometry::POINT: return &PointFE;
618  case Geometry::SEGMENT: return &SegmentFE;
619  case Geometry::TRIANGLE: return &TriangleFE;
620  case Geometry::SQUARE: return &QuadrilateralFE;
621  case Geometry::TETRAHEDRON: return &TetrahedronFE;
622  case Geometry::CUBE: return &ParallelepipedFE;
623  case Geometry::PRISM: return &WedgeFE;
624  default:
625  mfem_error ("QuadraticFECollection: unknown geometry type.");
626  }
627  return &SegmentFE; // Make some compilers happy
628 }
629 
631 {
632  switch (GeomType)
633  {
634  case Geometry::POINT: return 1;
635  case Geometry::SEGMENT: return 1;
636  case Geometry::TRIANGLE: return 0;
637  case Geometry::SQUARE: return 1;
638  case Geometry::TETRAHEDRON: return 0;
639  case Geometry::CUBE: return 1;
640  case Geometry::PRISM: return 0;
641  default:
642  mfem_error ("QuadraticFECollection: unknown geometry type.");
643  }
644  return 0; // Make some compilers happy
645 }
646 
648  Geometry::Type GeomType, int Or) const
649 {
650  static int indexes[] = { 0 };
651 
652  return indexes;
653 }
654 
655 
656 const FiniteElement *
658  Geometry::Type GeomType) const
659 {
660  switch (GeomType)
661  {
662  case Geometry::SEGMENT: return &SegmentFE;
663  case Geometry::SQUARE: return &QuadrilateralFE;
664  default:
665  mfem_error ("QuadraticPosFECollection: unknown geometry type.");
666  }
667  return NULL; // Make some compilers happy
668 }
669 
671 {
672  switch (GeomType)
673  {
674  case Geometry::POINT: return 1;
675  case Geometry::SEGMENT: return 1;
676  case Geometry::SQUARE: return 1;
677  default:
678  mfem_error ("QuadraticPosFECollection: unknown geometry type.");
679  }
680  return 0; // Make some compilers happy
681 }
682 
684  Geometry::Type GeomType, int Or) const
685 {
686  static int indexes[] = { 0 };
687 
688  return indexes;
689 }
690 
691 
692 const FiniteElement *
694 {
695  switch (GeomType)
696  {
697  case Geometry::POINT: return &PointFE;
698  case Geometry::SEGMENT: return &SegmentFE;
699  case Geometry::TRIANGLE: return &TriangleFE;
700  case Geometry::SQUARE: return &QuadrilateralFE;
701  case Geometry::TETRAHEDRON: return &TetrahedronFE;
702  case Geometry::CUBE: return &ParallelepipedFE;
703  case Geometry::PRISM: return &WedgeFE;
704  default:
705  mfem_error ("CubicFECollection: unknown geometry type.");
706  }
707  return &SegmentFE; // Make some compilers happy
708 }
709 
711 {
712  switch (GeomType)
713  {
714  case Geometry::POINT: return 1;
715  case Geometry::SEGMENT: return 2;
716  case Geometry::TRIANGLE: return 1;
717  case Geometry::SQUARE: return 4;
718  case Geometry::TETRAHEDRON: return 0;
719  case Geometry::CUBE: return 8;
720  case Geometry::PRISM: return 2;
721  default:
722  mfem_error ("CubicFECollection: unknown geometry type.");
723  }
724  return 0; // Make some compilers happy
725 }
726 
728  int Or) const
729 {
730  if (GeomType == Geometry::SEGMENT)
731  {
732  static int ind_pos[] = { 0, 1 };
733  static int ind_neg[] = { 1, 0 };
734 
735  if (Or < 0)
736  {
737  return ind_neg;
738  }
739  return ind_pos;
740  }
741  else if (GeomType == Geometry::TRIANGLE)
742  {
743  static int indexes[] = { 0 };
744 
745  return indexes;
746  }
747  else if (GeomType == Geometry::SQUARE)
748  {
749  static int sq_ind[8][4] = {{0, 1, 2, 3}, {0, 2, 1, 3},
750  {2, 0, 3, 1}, {1, 0, 3, 2},
751  {3, 2, 1, 0}, {3, 1, 2, 0},
752  {1, 3, 0, 2}, {2, 3, 0, 1}
753  };
754  return sq_ind[Or];
755  }
756 
757  return NULL;
758 }
759 
760 
761 const FiniteElement *
763  Geometry::Type GeomType) const
764 {
765  switch (GeomType)
766  {
767  case Geometry::SEGMENT: return &SegmentFE;
768  case Geometry::TRIANGLE: return &TriangleFE;
769  case Geometry::SQUARE: return &QuadrilateralFE;
770  default:
771  mfem_error ("CrouzeixRaviartFECollection: unknown geometry type.");
772  }
773  return &SegmentFE; // Make some compilers happy
774 }
775 
777 {
778  switch (GeomType)
779  {
780  case Geometry::POINT: return 0;
781  case Geometry::SEGMENT: return 1;
782  case Geometry::TRIANGLE: return 0;
783  case Geometry::SQUARE: return 0;
784  default:
785  mfem_error ("CrouzeixRaviartFECollection: unknown geometry type.");
786  }
787  return 0; // Make some compilers happy
788 }
789 
791  Geometry::Type GeomType, int Or) const
792 {
793  static int indexes[] = { 0 };
794 
795  return indexes;
796 }
797 
798 
799 const FiniteElement *
801 {
802  switch (GeomType)
803  {
804  case Geometry::SEGMENT: return &SegmentFE;
805  case Geometry::TRIANGLE: return &TriangleFE;
806  case Geometry::SQUARE: return &QuadrilateralFE;
807  default:
808  mfem_error ("RT0_2DFECollection: unknown geometry type.");
809  }
810  return &SegmentFE; // Make some compilers happy
811 }
812 
814 {
815  switch (GeomType)
816  {
817  case Geometry::POINT: return 0;
818  case Geometry::SEGMENT: return 1;
819  case Geometry::TRIANGLE: return 0;
820  case Geometry::SQUARE: return 0;
821  default:
822  mfem_error ("RT0_2DFECollection: unknown geometry type.");
823  }
824  return 0; // Make some compilers happy
825 }
826 
828  int Or) const
829 {
830  static int ind_pos[] = { 0 };
831  static int ind_neg[] = { -1 };
832 
833  if (Or > 0)
834  {
835  return ind_pos;
836  }
837  return ind_neg;
838 }
839 
840 
841 const FiniteElement *
843 {
844  switch (GeomType)
845  {
846  case Geometry::SEGMENT: return &SegmentFE;
847  case Geometry::TRIANGLE: return &TriangleFE;
848  case Geometry::SQUARE: return &QuadrilateralFE;
849  default:
850  mfem_error ("RT1_2DFECollection: unknown geometry type.");
851  }
852  return &SegmentFE; // Make some compilers happy
853 }
854 
856 {
857  switch (GeomType)
858  {
859  case Geometry::POINT: return 0;
860  case Geometry::SEGMENT: return 2;
861  case Geometry::TRIANGLE: return 2;
862  case Geometry::SQUARE: return 4;
863  default:
864  mfem_error ("RT1_2DFECollection: unknown geometry type.");
865  }
866  return 0; // Make some compilers happy
867 }
868 
870  int Or) const
871 {
872  static int ind_pos[] = { 0, 1 };
873  static int ind_neg[] = { -2, -1 };
874 
875  if (Or > 0)
876  {
877  return ind_pos;
878  }
879  return ind_neg;
880 }
881 
882 const FiniteElement *
884 {
885  switch (GeomType)
886  {
887  case Geometry::SEGMENT: return &SegmentFE;
888  case Geometry::TRIANGLE: return &TriangleFE;
889  case Geometry::SQUARE: return &QuadrilateralFE;
890  default:
891  mfem_error ("RT2_2DFECollection: unknown geometry type.");
892  }
893  return &SegmentFE; // Make some compilers happy
894 }
895 
897 {
898  switch (GeomType)
899  {
900  case Geometry::POINT: return 0;
901  case Geometry::SEGMENT: return 3;
902  case Geometry::TRIANGLE: return 6;
903  case Geometry::SQUARE: return 12;
904  default:
905  mfem_error ("RT2_2DFECollection: unknown geometry type.");
906  }
907  return 0; // Make some compilers happy
908 }
909 
911  int Or) const
912 {
913  static int ind_pos[] = { 0, 1, 2 };
914  static int ind_neg[] = { -3, -2, -1 };
915 
916  if (Or > 0)
917  {
918  return ind_pos;
919  }
920  return ind_neg;
921 }
922 
923 
924 const FiniteElement *
926 {
927  switch (GeomType)
928  {
929  case Geometry::TRIANGLE: return &TriangleFE;
930  case Geometry::SQUARE: return &QuadrilateralFE;
931  default:
932  mfem_error ("Const2DFECollection: unknown geometry type.");
933  }
934  return &TriangleFE; // Make some compilers happy
935 }
936 
938 {
939  switch (GeomType)
940  {
941  case Geometry::POINT: return 0;
942  case Geometry::SEGMENT: return 0;
943  case Geometry::TRIANGLE: return 1;
944  case Geometry::SQUARE: return 1;
945  default:
946  mfem_error ("Const2DFECollection: unknown geometry type.");
947  }
948  return 0; // Make some compilers happy
949 }
950 
952  int Or) const
953 {
954  return NULL;
955 }
956 
957 
958 const FiniteElement *
960  Geometry::Type GeomType) const
961 {
962  switch (GeomType)
963  {
964  case Geometry::TRIANGLE: return &TriangleFE;
965  case Geometry::SQUARE: return &QuadrilateralFE;
966  default:
967  mfem_error ("LinearDiscont2DFECollection: unknown geometry type.");
968  }
969  return &TriangleFE; // Make some compilers happy
970 }
971 
973 {
974  switch (GeomType)
975  {
976  case Geometry::POINT: return 0;
977  case Geometry::SEGMENT: return 0;
978  case Geometry::TRIANGLE: return 3;
979  case Geometry::SQUARE: return 4;
980  default:
981  mfem_error ("LinearDiscont2DFECollection: unknown geometry type.");
982  }
983  return 0; // Make some compilers happy
984 }
985 
987  Geometry::Type GeomType, int Or) const
988 {
989  return NULL;
990 }
991 
992 
993 const FiniteElement *
995  Geometry::Type GeomType) const
996 {
997  switch (GeomType)
998  {
999  case Geometry::TRIANGLE: return &TriangleFE;
1000  case Geometry::SQUARE: return &QuadrilateralFE;
1001  default:
1002  mfem_error ("GaussLinearDiscont2DFECollection:"
1003  " unknown geometry type.");
1004  }
1005  return &TriangleFE; // Make some compilers happy
1006 }
1007 
1009 const
1010 {
1011  switch (GeomType)
1012  {
1013  case Geometry::POINT: return 0;
1014  case Geometry::SEGMENT: return 0;
1015  case Geometry::TRIANGLE: return 3;
1016  case Geometry::SQUARE: return 4;
1017  default:
1018  mfem_error ("GaussLinearDiscont2DFECollection:"
1019  " unknown geometry type.");
1020  }
1021  return 0; // Make some compilers happy
1022 }
1023 
1025  Geometry::Type GeomType, int Or) const
1026 {
1027  return NULL;
1028 }
1029 
1030 
1031 const FiniteElement *
1033 {
1034  if (GeomType != Geometry::SQUARE)
1035  {
1036  mfem_error ("P1OnQuadFECollection: unknown geometry type.");
1037  }
1038  return &QuadrilateralFE;
1039 }
1040 
1042 {
1043  switch (GeomType)
1044  {
1045  case Geometry::POINT: return 0;
1046  case Geometry::SEGMENT: return 0;
1047  case Geometry::SQUARE: return 3;
1048  default:
1049  mfem_error ("P1OnQuadFECollection: unknown geometry type.");
1050  }
1051  return 0; // Make some compilers happy
1052 }
1053 
1055  Geometry::Type GeomType, int Or) const
1056 {
1057  return NULL;
1058 }
1059 
1060 
1061 const FiniteElement *
1063  Geometry::Type GeomType) const
1064 {
1065  switch (GeomType)
1066  {
1067  case Geometry::TRIANGLE: return &TriangleFE;
1068  case Geometry::SQUARE: return &QuadrilateralFE;
1069  default:
1070  mfem_error ("QuadraticDiscont2DFECollection: unknown geometry type.");
1071  }
1072  return &TriangleFE; // Make some compilers happy
1073 }
1074 
1076 const
1077 {
1078  switch (GeomType)
1079  {
1080  case Geometry::POINT: return 0;
1081  case Geometry::SEGMENT: return 0;
1082  case Geometry::TRIANGLE: return 6;
1083  case Geometry::SQUARE: return 9;
1084  default:
1085  mfem_error ("QuadraticDiscont2DFECollection: unknown geometry type.");
1086  }
1087  return 0; // Make some compilers happy
1088 }
1089 
1091  Geometry::Type GeomType, int Or) const
1092 {
1093  return NULL;
1094 }
1095 
1096 
1097 const FiniteElement *
1099  Geometry::Type GeomType) const
1100 {
1101  switch (GeomType)
1102  {
1103  case Geometry::SQUARE: return &QuadrilateralFE;
1104  default:
1105  mfem_error ("QuadraticPosDiscont2DFECollection: unknown geometry type.");
1106  }
1107  return NULL; // Make some compilers happy
1108 }
1109 
1111 const
1112 {
1113  switch (GeomType)
1114  {
1115  case Geometry::POINT: return 0;
1116  case Geometry::SEGMENT: return 0;
1117  case Geometry::SQUARE: return 9;
1118  default:
1119  mfem_error ("QuadraticPosDiscont2DFECollection: unknown geometry type.");
1120  }
1121  return 0; // Make some compilers happy
1122 }
1123 
1124 
1125 const FiniteElement *
1127  Geometry::Type GeomType)
1128 const
1129 {
1130  switch (GeomType)
1131  {
1132  case Geometry::TRIANGLE: return &TriangleFE;
1133  case Geometry::SQUARE: return &QuadrilateralFE;
1134  default:
1135  mfem_error ("GaussQuadraticDiscont2DFECollection:"
1136  " unknown geometry type.");
1137  }
1138  return &QuadrilateralFE; // Make some compilers happy
1139 }
1140 
1142  Geometry::Type GeomType) const
1143 {
1144  switch (GeomType)
1145  {
1146  case Geometry::POINT: return 0;
1147  case Geometry::SEGMENT: return 0;
1148  case Geometry::TRIANGLE: return 6;
1149  case Geometry::SQUARE: return 9;
1150  default:
1151  mfem_error ("GaussQuadraticDiscont2DFECollection:"
1152  " unknown geometry type.");
1153  }
1154  return 0; // Make some compilers happy
1155 }
1156 
1158  Geometry::Type GeomType, int Or) const
1159 {
1160  return NULL;
1161 }
1162 
1163 
1164 const FiniteElement *
1166  Geometry::Type GeomType) const
1167 {
1168  switch (GeomType)
1169  {
1170  case Geometry::TRIANGLE: return &TriangleFE;
1171  case Geometry::SQUARE: return &QuadrilateralFE;
1172  default:
1173  mfem_error ("CubicDiscont2DFECollection: unknown geometry type.");
1174  }
1175  return &TriangleFE; // Make some compilers happy
1176 }
1177 
1179 {
1180  switch (GeomType)
1181  {
1182  case Geometry::POINT: return 0;
1183  case Geometry::SEGMENT: return 0;
1184  case Geometry::TRIANGLE: return 10;
1185  case Geometry::SQUARE: return 16;
1186  default:
1187  mfem_error ("CubicDiscont2DFECollection: unknown geometry type.");
1188  }
1189  return 0; // Make some compilers happy
1190 }
1191 
1193  Geometry::Type GeomType, int Or) const
1194 {
1195  return NULL;
1196 }
1197 
1198 
1199 const FiniteElement *
1201  Geometry::Type GeomType) const
1202 {
1203  switch (GeomType)
1204  {
1205  case Geometry::TRIANGLE: return &TriangleFE;
1206  case Geometry::SQUARE: return &QuadrilateralFE;
1207  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1208  case Geometry::CUBE: return &ParallelepipedFE;
1209  default:
1210  mfem_error ("LinearNonConf3DFECollection: unknown geometry type.");
1211  }
1212  return &TriangleFE; // Make some compilers happy
1213 }
1214 
1216 {
1217  switch (GeomType)
1218  {
1219  case Geometry::POINT: return 0;
1220  case Geometry::SEGMENT: return 0;
1221  case Geometry::TRIANGLE: return 1;
1222  case Geometry::SQUARE: return 1;
1223  case Geometry::TETRAHEDRON: return 0;
1224  case Geometry::CUBE: return 0;
1225  default:
1226  mfem_error ("LinearNonConf3DFECollection: unknown geometry type.");
1227  }
1228  return 0; // Make some compilers happy
1229 }
1230 
1232  Geometry::Type GeomType, int Or) const
1233 {
1234  static int indexes[] = { 0 };
1235 
1236  return indexes;
1237 }
1238 
1239 
1240 const FiniteElement *
1242 {
1243  switch (GeomType)
1244  {
1245  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1246  case Geometry::CUBE: return &ParallelepipedFE;
1247  case Geometry::PRISM: return &WedgeFE;
1248  case Geometry::PYRAMID: return &PyramidFE;
1249  default:
1250  mfem_error ("Const3DFECollection: unknown geometry type.");
1251  }
1252  return &TetrahedronFE; // Make some compilers happy
1253 }
1254 
1256 {
1257  switch (GeomType)
1258  {
1259  case Geometry::POINT: return 0;
1260  case Geometry::SEGMENT: return 0;
1261  case Geometry::TRIANGLE: return 0;
1262  case Geometry::SQUARE: return 0;
1263  case Geometry::TETRAHEDRON: return 1;
1264  case Geometry::CUBE: return 1;
1265  case Geometry::PRISM: return 1;
1266  case Geometry::PYRAMID: return 1;
1267  default:
1268  mfem_error ("Const3DFECollection: unknown geometry type.");
1269  }
1270  return 0; // Make some compilers happy
1271 }
1272 
1274  int Or) const
1275 {
1276  return NULL;
1277 }
1278 
1279 
1280 const FiniteElement *
1282  Geometry::Type GeomType) const
1283 {
1284  switch (GeomType)
1285  {
1286  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1287  case Geometry::PYRAMID: return &PyramidFE;
1288  case Geometry::PRISM: return &WedgeFE;
1289  case Geometry::CUBE: return &ParallelepipedFE;
1290  default:
1291  mfem_error ("LinearDiscont3DFECollection: unknown geometry type.");
1292  }
1293  return &TetrahedronFE; // Make some compilers happy
1294 }
1295 
1297 {
1298  switch (GeomType)
1299  {
1300  case Geometry::POINT: return 0;
1301  case Geometry::SEGMENT: return 0;
1302  case Geometry::TRIANGLE: return 0;
1303  case Geometry::SQUARE: return 0;
1304  case Geometry::TETRAHEDRON: return 4;
1305  case Geometry::PYRAMID: return 5;
1306  case Geometry::PRISM: return 6;
1307  case Geometry::CUBE: return 8;
1308  default:
1309  mfem_error ("LinearDiscont3DFECollection: unknown geometry type.");
1310  }
1311  return 0; // Make some compilers happy
1312 }
1313 
1315  Geometry::Type GeomType, int Or) const
1316 {
1317  return NULL;
1318 }
1319 
1320 
1321 const FiniteElement *
1323  Geometry::Type GeomType) const
1324 {
1325  switch (GeomType)
1326  {
1327  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1328  case Geometry::CUBE: return &ParallelepipedFE;
1329  default:
1330  mfem_error ("QuadraticDiscont3DFECollection: unknown geometry type.");
1331  }
1332  return &TetrahedronFE; // Make some compilers happy
1333 }
1334 
1336 const
1337 {
1338  switch (GeomType)
1339  {
1340  case Geometry::POINT: return 0;
1341  case Geometry::SEGMENT: return 0;
1342  case Geometry::TRIANGLE: return 0;
1343  case Geometry::SQUARE: return 0;
1344  case Geometry::TETRAHEDRON: return 10;
1345  case Geometry::CUBE: return 27;
1346  default:
1347  mfem_error ("QuadraticDiscont3DFECollection: unknown geometry type.");
1348  }
1349  return 0; // Make some compilers happy
1350 }
1351 
1353  Geometry::Type GeomType, int Or) const
1354 {
1355  return NULL;
1356 }
1357 
1358 const FiniteElement *
1360  Geometry::Type GeomType) const
1361 {
1362  switch (GeomType)
1363  {
1364  case Geometry::POINT: return &PointFE;
1365  case Geometry::SEGMENT: return &SegmentFE;
1366  case Geometry::TRIANGLE: return &TriangleFE;
1367  case Geometry::SQUARE: return &QuadrilateralFE;
1368  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1369  case Geometry::CUBE: return &ParallelepipedFE;
1370  default:
1371  mfem_error ("RefinedLinearFECollection: unknown geometry type.");
1372  }
1373  return &SegmentFE; // Make some compilers happy
1374 }
1375 
1377 {
1378  switch (GeomType)
1379  {
1380  case Geometry::POINT: return 1;
1381  case Geometry::SEGMENT: return 1;
1382  case Geometry::TRIANGLE: return 0;
1383  case Geometry::SQUARE: return 1;
1384  case Geometry::TETRAHEDRON: return 0;
1385  case Geometry::CUBE: return 1;
1386  default:
1387  mfem_error ("RefinedLinearFECollection: unknown geometry type.");
1388  }
1389  return 0; // Make some compilers happy
1390 }
1391 
1393  Geometry::Type GeomType, int Or) const
1394 {
1395  static int indexes[] = { 0 };
1396 
1397  return indexes;
1398 }
1399 
1400 
1401 const FiniteElement *
1403 {
1404  switch (GeomType)
1405  {
1406  case Geometry::CUBE: return &HexahedronFE;
1407  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1408  case Geometry::PRISM: return &WedgeFE;
1409  case Geometry::PYRAMID: return &PyramidFE;
1410  default:
1411  mfem_error ("ND1_3DFECollection: unknown geometry type.");
1412  }
1413  return &HexahedronFE; // Make some compilers happy
1414 }
1415 
1417 {
1418  switch (GeomType)
1419  {
1420  case Geometry::POINT: return 0;
1421  case Geometry::SEGMENT: return 1;
1422  case Geometry::TRIANGLE: return 0;
1423  case Geometry::SQUARE: return 0;
1424  case Geometry::TETRAHEDRON: return 0;
1425  case Geometry::CUBE: return 0;
1426  case Geometry::PRISM: return 0;
1427  case Geometry::PYRAMID: return 0;
1428  default:
1429  mfem_error ("ND1_3DFECollection: unknown geometry type.");
1430  }
1431  return 0; // Make some compilers happy
1432 }
1433 
1435  int Or) const
1436 {
1437  static int ind_pos[] = { 0 };
1438  static int ind_neg[] = { -1 };
1439 
1440  if (Or > 0)
1441  {
1442  return ind_pos;
1443  }
1444  return ind_neg;
1445 }
1446 
1447 
1448 const FiniteElement *
1450 {
1451  switch (GeomType)
1452  {
1453  case Geometry::TRIANGLE: return &TriangleFE;
1454  case Geometry::SQUARE: return &QuadrilateralFE;
1455  case Geometry::CUBE: return &HexahedronFE;
1456  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1457  case Geometry::PRISM: return &WedgeFE;
1458  case Geometry::PYRAMID: return &PyramidFE;
1459  default:
1460  mfem_error ("RT0_3DFECollection: unknown geometry type.");
1461  }
1462  return &HexahedronFE; // Make some compilers happy
1463 }
1464 
1466 {
1467  switch (GeomType)
1468  {
1469  case Geometry::POINT: return 0;
1470  case Geometry::SEGMENT: return 0;
1471  case Geometry::TRIANGLE: return 1;
1472  case Geometry::SQUARE: return 1;
1473  case Geometry::TETRAHEDRON: return 0;
1474  case Geometry::CUBE: return 0;
1475  case Geometry::PRISM: return 0;
1476  case Geometry::PYRAMID: return 0;
1477  default:
1478  mfem_error ("RT0_3DFECollection: unknown geometry type.");
1479  }
1480  return 0; // Make some compilers happy
1481 }
1482 
1484  int Or) const
1485 {
1486  static int ind_pos[] = { 0 };
1487  static int ind_neg[] = { -1 };
1488 
1489  if ((GeomType == Geometry::TRIANGLE) || (GeomType == Geometry::SQUARE))
1490  {
1491  if (Or % 2 == 0)
1492  {
1493  return ind_pos;
1494  }
1495  return ind_neg;
1496  }
1497  return NULL;
1498 }
1499 
1500 const FiniteElement *
1502 {
1503  switch (GeomType)
1504  {
1505  case Geometry::TRIANGLE: return &TriangleFE;
1506  case Geometry::SQUARE: return &QuadrilateralFE;
1507  case Geometry::CUBE: return &HexahedronFE;
1508  default:
1509  mfem_error ("RT1_3DFECollection: unknown geometry type.");
1510  }
1511  return &HexahedronFE; // Make some compilers happy
1512 }
1513 
1515 {
1516  switch (GeomType)
1517  {
1518  case Geometry::POINT: return 0;
1519  case Geometry::SEGMENT: return 0;
1520  case Geometry::TRIANGLE: return 2;
1521  case Geometry::SQUARE: return 4;
1522  case Geometry::CUBE: return 12;
1523  default:
1524  mfem_error ("RT1_3DFECollection: unknown geometry type.");
1525  }
1526  return 0; // Make some compilers happy
1527 }
1528 
1530  int Or) const
1531 {
1532  if (GeomType == Geometry::SQUARE)
1533  {
1534  static int sq_ind[8][4] =
1535  {
1536  {0, 1, 2, 3}, {-1, -3, -2, -4},
1537  {2, 0, 3, 1}, {-2, -1, -4, -3},
1538  {3, 2, 1, 0}, {-4, -2, -3, -1},
1539  {1, 3, 0, 2}, {-3, -4, -1, -2}
1540  };
1541 
1542  return sq_ind[Or];
1543  }
1544  else
1545  {
1546  return NULL;
1547  }
1548 }
1549 
1550 
1551 H1_FECollection::H1_FECollection(const int p, const int dim, const int btype)
1553  , dim(dim)
1554 {
1555  MFEM_VERIFY(p >= 1, "H1_FECollection requires order >= 1.");
1556  MFEM_VERIFY(dim >= 0 && dim <= 3, "H1_FECollection requires 0 <= dim <= 3.");
1557 
1558  const int pm1 = p - 1, pm2 = pm1 - 1, pm3 = pm2 - 1, pm4 = pm3 - 1;
1559 
1560  int pt_type = BasisType::GetQuadrature1D(btype);
1561  b_type = BasisType::Check(btype);
1562  switch (btype)
1563  {
1565  {
1566  snprintf(h1_name, 32, "H1_%dD_P%d", dim, p);
1567  break;
1568  }
1569  case BasisType::Positive:
1570  {
1571  snprintf(h1_name, 32, "H1Pos_%dD_P%d", dim, p);
1572  break;
1573  }
1575  {
1576  snprintf(h1_name, 32, "H1Ser_%dD_P%d", dim, p);
1577  break;
1578  }
1579  default:
1580  {
1581  MFEM_VERIFY(Quadrature1D::CheckClosed(pt_type) !=
1583  "unsupported BasisType: " << BasisType::Name(btype));
1584 
1585  snprintf(h1_name, 32, "H1@%c_%dD_P%d",
1586  (int)BasisType::GetChar(btype), dim, p);
1587  }
1588  }
1589 
1590  for (int g = 0; g < Geometry::NumGeom; g++)
1591  {
1592  H1_dof[g] = 0;
1593  H1_Elements[g] = NULL;
1594  }
1595  for (int i = 0; i < 2; i++)
1596  {
1597  SegDofOrd[i] = NULL;
1598  }
1599  for (int i = 0; i < 6; i++)
1600  {
1601  TriDofOrd[i] = NULL;
1602  }
1603  for (int i = 0; i < 8; i++)
1604  {
1605  QuadDofOrd[i] = NULL;
1606  }
1607  for (int i = 0; i < 24; i++)
1608  {
1609  TetDofOrd[i] = NULL;
1610  }
1611 
1612  H1_dof[Geometry::POINT] = 1;
1614 
1615  if (dim >= 1)
1616  {
1617  H1_dof[Geometry::SEGMENT] = pm1;
1618  if (b_type == BasisType::Positive)
1619  {
1621  }
1622  else
1623  {
1625  }
1626 
1627  SegDofOrd[0] = new int[2*pm1];
1628  SegDofOrd[1] = SegDofOrd[0] + pm1;
1629  for (int i = 0; i < pm1; i++)
1630  {
1631  SegDofOrd[0][i] = i;
1632  SegDofOrd[1][i] = pm2 - i;
1633  }
1634  }
1635 
1636  if (dim >= 2)
1637  {
1638  H1_dof[Geometry::TRIANGLE] = (pm1*pm2)/2;
1639  H1_dof[Geometry::SQUARE] = pm1*pm1;
1640  if (b_type == BasisType::Positive)
1641  {
1644  }
1645  else if (b_type == BasisType::Serendipity)
1646  {
1647  // Note: in fe_coll.hpp the DofForGeometry(Geometry::Type) method
1648  // returns H1_dof[GeomType], so we need to fix the value of H1_dof here
1649  // for the serendipity case.
1650 
1651  // formula for number of interior serendipity DoFs (when p>1)
1652  H1_dof[Geometry::SQUARE] = (pm3*pm2)/2;
1654  // allows for mixed tri/quad meshes
1656  }
1657  else
1658  {
1661  }
1662 
1663  const int &TriDof = H1_dof[Geometry::TRIANGLE];
1664  const int &QuadDof = H1_dof[Geometry::SQUARE];
1665  TriDofOrd[0] = new int[6*TriDof];
1666  for (int i = 1; i < 6; i++)
1667  {
1668  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1669  }
1670  // see Mesh::GetTriOrientation in mesh/mesh.cpp
1671  for (int j = 0; j < pm2; j++)
1672  {
1673  for (int i = 0; i + j < pm2; i++)
1674  {
1675  int o = TriDof - ((pm1 - j)*(pm2 - j))/2 + i;
1676  int k = pm3 - j - i;
1677  TriDofOrd[0][o] = o; // (0,1,2)
1678  TriDofOrd[1][o] = TriDof - ((pm1-j)*(pm2-j))/2 + k; // (1,0,2)
1679  TriDofOrd[2][o] = TriDof - ((pm1-i)*(pm2-i))/2 + k; // (2,0,1)
1680  TriDofOrd[3][o] = TriDof - ((pm1-k)*(pm2-k))/2 + i; // (2,1,0)
1681  TriDofOrd[4][o] = TriDof - ((pm1-k)*(pm2-k))/2 + j; // (1,2,0)
1682  TriDofOrd[5][o] = TriDof - ((pm1-i)*(pm2-i))/2 + j; // (0,2,1)
1683  }
1684  }
1685 
1686  QuadDofOrd[0] = new int[8*QuadDof];
1687  for (int i = 1; i < 8; i++)
1688  {
1689  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
1690  }
1691 
1692  // For serendipity order >=4, the QuadDofOrd array must be re-defined. We
1693  // do this by computing the corresponding tensor product QuadDofOrd array
1694  // or two orders less, which contains enough DoFs for their serendipity
1695  // basis. This could be optimized.
1697  {
1698  if (p < 4)
1699  {
1700  // no face dofs --> don't need to adjust QuadDofOrd
1701  }
1702  else // p >= 4 --> have face dofs
1703  {
1704  // Exactly the same as tensor product case, but with all orders
1705  // reduced by 2 e.g. in case p=5 it builds a 2x2 array, even though
1706  // there are only 3 serendipity dofs.
1707  // In the tensor product case, the i and j index tensor directions,
1708  // and o index from 0 to (pm1)^2,
1709 
1710  for (int j = 0; j < pm3; j++) // pm3 instead of pm1, etc
1711  {
1712  for (int i = 0; i < pm3; i++)
1713  {
1714  int o = i + j*pm3;
1715  QuadDofOrd[0][o] = i + j*pm3; // (0,1,2,3)
1716  QuadDofOrd[1][o] = j + i*pm3; // (0,3,2,1)
1717  QuadDofOrd[2][o] = j + (pm4 - i)*pm3; // (1,2,3,0)
1718  QuadDofOrd[3][o] = (pm4 - i) + j*pm3; // (1,0,3,2)
1719  QuadDofOrd[4][o] = (pm4 - i) + (pm4 - j)*pm3; // (2,3,0,1)
1720  QuadDofOrd[5][o] = (pm4 - j) + (pm4 - i)*pm3; // (2,1,0,3)
1721  QuadDofOrd[6][o] = (pm4 - j) + i*pm3; // (3,0,1,2)
1722  QuadDofOrd[7][o] = i + (pm4 - j)*pm3; // (3,2,1,0)
1723  }
1724  }
1725 
1726  }
1727  }
1728  else // not serendipity
1729  {
1730  for (int j = 0; j < pm1; j++)
1731  {
1732  for (int i = 0; i < pm1; i++)
1733  {
1734  int o = i + j*pm1;
1735  QuadDofOrd[0][o] = i + j*pm1; // (0,1,2,3)
1736  QuadDofOrd[1][o] = j + i*pm1; // (0,3,2,1)
1737  QuadDofOrd[2][o] = j + (pm2 - i)*pm1; // (1,2,3,0)
1738  QuadDofOrd[3][o] = (pm2 - i) + j*pm1; // (1,0,3,2)
1739  QuadDofOrd[4][o] = (pm2 - i) + (pm2 - j)*pm1; // (2,3,0,1)
1740  QuadDofOrd[5][o] = (pm2 - j) + (pm2 - i)*pm1; // (2,1,0,3)
1741  QuadDofOrd[6][o] = (pm2 - j) + i*pm1; // (3,0,1,2)
1742  QuadDofOrd[7][o] = i + (pm2 - j)*pm1; // (3,2,1,0)
1743  }
1744  }
1745  }
1746 
1747  if (dim >= 3)
1748  {
1749  H1_dof[Geometry::TETRAHEDRON] = (TriDof*pm3)/3;
1750  H1_dof[Geometry::CUBE] = QuadDof*pm1;
1751  H1_dof[Geometry::PRISM] = TriDof*pm1;
1753  if (b_type == BasisType::Positive)
1754  {
1758  }
1759  else
1760  {
1762  new H1_TetrahedronElement(p, btype);
1764  H1_Elements[Geometry::PRISM] = new H1_WedgeElement(p, btype);
1765  }
1767 
1768  const int &TetDof = H1_dof[Geometry::TETRAHEDRON];
1769  TetDofOrd[0] = new int[24*TetDof];
1770  for (int i = 1; i < 24; i++)
1771  {
1772  TetDofOrd[i] = TetDofOrd[i-1] + TetDof;
1773  }
1774  // see Mesh::GetTetOrientation in mesh/mesh.cpp
1775  for (int k = 0; k < pm3; k++)
1776  {
1777  for (int j = 0; j + k < pm3; j++)
1778  {
1779  for (int i = 0; i + j + k < pm3; i++)
1780  {
1781  int l = pm4 - k - j - i;
1782  int o = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1783  + (j * (2 * p - 5 - j - 2 * k)) / 2 + i;
1784  int o1 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1785  + (k * (2 * p - 5 - k - 2 * j)) / 2 + i;
1786  int o2 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1787  + (k * (2 * p - 5 - k - 2 * i)) / 2 + j;
1788  int o3 = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1789  + (i * (2 * p - 5 - i - 2 * k)) / 2 + j;
1790  int o4 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1791  + (i * (2 * p - 5 - i - 2 * j)) / 2 + k;
1792  int o5 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1793  + (j * (2 * p - 5 - j - 2 * i)) / 2 + k;
1794  int o6 = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1795  + (l * (2 * p - 5 - l - 2 * k)) / 2 + j;
1796  int o7 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1797  + (k * (2 * p - 5 - k - 2 * l)) / 2 + j;
1798  int o8 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1799  + (j * (2 * p - 5 - j - 2 * l)) / 2 + k;
1800  int o9 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1801  + (l * (2 * p - 5 - l - 2 * j)) / 2 + k;
1802  int o10 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1803  + (k * (2 * p - 5 - k - 2 * j)) / 2 + l;
1804  int o11 = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1805  + (j * (2 * p - 5 - j - 2 * k)) / 2 + l;
1806  int o12 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1807  + (l * (2 * p - 5 - l - 2 * i)) / 2 + k;
1808  int o13 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1809  + (i * (2 * p - 5 - i - 2 * l)) / 2 + k;
1810  int o14 = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1811  + (i * (2 * p - 5 - i - 2 * k)) / 2 + l;
1812  int o15 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1813  + (k * (2 * p - 5 - k - 2 * i)) / 2 + l;
1814  int o16 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1815  + (k * (2 * p - 5 - k - 2 * l)) / 2 + i;
1816  int o17 = TetDof - ((pm1 - k) * (pm2 - k) * (pm3 - k)) / 6
1817  + (l * (2 * p - 5 - l - 2 * k)) / 2 + i;
1818  int o18 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1819  + (j * (2 * p - 5 - j - 2 * i)) / 2 + l;
1820  int o19 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1821  + (i * (2 * p - 5 - i - 2 * j)) / 2 + l;
1822  int o20 = TetDof - ((pm1 - j) * (pm2 - j) * (pm3 - j)) / 6
1823  + (l * (2 * p - 5 - l - 2 * j)) / 2 + i;
1824  int o21 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1825  + (j * (2 * p - 5 - j - 2 * l)) / 2 + i;
1826  int o22 = TetDof - ((pm1 - l) * (pm2 - l) * (pm3 - l)) / 6
1827  + (i * (2 * p - 5 - i - 2 * l)) / 2 + j;
1828  int o23 = TetDof - ((pm1 - i) * (pm2 - i) * (pm3 - i)) / 6
1829  + (l * (2 * p - 5 - l - 2 * i)) / 2 + j;
1830  TetDofOrd[ 0][o] = o; // (0,1,2,3)
1831  TetDofOrd[ 1][o] = o1; // (0,1,3,2)
1832  TetDofOrd[ 2][o] = o2; // (0,2,3,1)
1833  TetDofOrd[ 3][o] = o3; // (0,2,1,3)
1834  TetDofOrd[ 4][o] = o4; // (0,3,1,2)
1835  TetDofOrd[ 5][o] = o5; // (0,3,2,1)
1836  TetDofOrd[ 6][o] = o6; // (1,2,0,3)
1837  TetDofOrd[ 7][o] = o7; // (1,2,3,0)
1838  TetDofOrd[ 8][o] = o8; // (1,3,2,0)
1839  TetDofOrd[ 9][o] = o9; // (1,3,0,2)
1840  TetDofOrd[10][o] = o10; // (1,0,3,2)
1841  TetDofOrd[11][o] = o11; // (1,0,2,3)
1842  TetDofOrd[12][o] = o12; // (2,3,0,1)
1843  TetDofOrd[13][o] = o13; // (2,3,1,0)
1844  TetDofOrd[14][o] = o14; // (2,0,1,3)
1845  TetDofOrd[15][o] = o15; // (2,0,3,1)
1846  TetDofOrd[16][o] = o16; // (2,1,3,0)
1847  TetDofOrd[17][o] = o17; // (2,1,0,3)
1848  TetDofOrd[18][o] = o18; // (3,0,2,1)
1849  TetDofOrd[19][o] = o19; // (3,0,1,2)
1850  TetDofOrd[20][o] = o20; // (3,1,0,2)
1851  TetDofOrd[21][o] = o21; // (3,1,2,0)
1852  TetDofOrd[22][o] = o22; // (3,2,1,0)
1853  TetDofOrd[23][o] = o23; // (3,2,0,1)
1854  }
1855  }
1856  }
1857  }
1858  }
1859 }
1860 
1861 const FiniteElement *
1863 {
1864  if (GeomType != Geometry::PYRAMID || this->GetOrder() == 1)
1865  {
1866  return H1_Elements[GeomType];
1867  }
1868  else
1869  {
1870  MFEM_ABORT("H1 Pyramid basis functions are not yet supported "
1871  "for order > 1.");
1872  return NULL;
1873  }
1874 }
1875 
1877  int Or) const
1878 {
1879  if (GeomType == Geometry::SEGMENT)
1880  {
1881  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
1882  }
1883  else if (GeomType == Geometry::TRIANGLE)
1884  {
1885  return TriDofOrd[Or%6];
1886  }
1887  else if (GeomType == Geometry::SQUARE)
1888  {
1889  return QuadDofOrd[Or%8];
1890  }
1891  else if (GeomType == Geometry::TETRAHEDRON)
1892  {
1893  return TetDofOrd[Or%24];
1894  }
1895  return NULL;
1896 }
1897 
1899 {
1900  int tr_p = H1_dof[Geometry::SEGMENT] + 1;
1901  int tr_dim = -1;
1902  if (!strncmp(h1_name, "H1_", 3))
1903  {
1904  tr_dim = atoi(h1_name + 3);
1905  }
1906  else if (!strncmp(h1_name, "H1Pos_", 6))
1907  {
1908  tr_dim = atoi(h1_name + 6);
1909  }
1910  else if (!strncmp(h1_name, "H1@", 3))
1911  {
1912  tr_dim = atoi(h1_name + 5);
1913  }
1914  return (dim < 0) ? NULL : new H1_Trace_FECollection(tr_p, tr_dim, b_type);
1915 }
1916 
1917 const int *H1_FECollection::GetDofMap(Geometry::Type GeomType) const
1918 {
1919  const int *dof_map = NULL;
1920  const FiniteElement *fe = H1_Elements[GeomType];
1921  const NodalFiniteElement *nodal_fe =
1922  dynamic_cast<const NodalFiniteElement*>(fe);
1923  if (nodal_fe)
1924  {
1925  dof_map = nodal_fe->GetLexicographicOrdering().GetData();
1926  }
1927  else
1928  {
1929  MFEM_ABORT("Geometry type " << Geometry::Name[GeomType] << " is not "
1930  "implemented");
1931  }
1932  return dof_map;
1933 }
1934 
1935 const int *H1_FECollection::GetDofMap(Geometry::Type GeomType, int p) const
1936 {
1937  if (p == base_p) { return GetDofMap(GeomType); }
1938  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
1939  return ((H1_FECollection*) var_orders[p])->GetDofMap(GeomType);
1940 }
1941 
1943 {
1944  delete [] SegDofOrd[0];
1945  delete [] TriDofOrd[0];
1946  delete [] QuadDofOrd[0];
1947  delete [] TetDofOrd[0];
1948  for (int g = 0; g < Geometry::NumGeom; g++)
1949  {
1950  delete H1_Elements[g];
1951  }
1952 }
1953 
1954 
1956  const int btype)
1957  : H1_FECollection(p, dim-1, btype)
1958 {
1959  if (btype == BasisType::GaussLobatto)
1960  {
1961  snprintf(h1_name, 32, "H1_Trace_%dD_P%d", dim, p);
1962  }
1963  else if (btype == BasisType::Positive)
1964  {
1965  snprintf(h1_name, 32, "H1Pos_Trace_%dD_P%d", dim, p);
1966  }
1967  else // base class checks that type is closed
1968  {
1969  snprintf(h1_name, 32, "H1_Trace@%c_%dD_P%d",
1970  (int)BasisType::GetChar(btype), dim, p);
1971  }
1972 }
1973 
1974 
1975 L2_FECollection::L2_FECollection(const int p, const int dim, const int btype,
1976  const int map_type)
1978  , dim(dim)
1979  , m_type(map_type)
1980 {
1981  MFEM_VERIFY(p >= 0, "L2_FECollection requires order >= 0.");
1982 
1983  b_type = BasisType::Check(btype);
1984  const char *prefix = NULL;
1985  switch (map_type)
1986  {
1987  case FiniteElement::VALUE: prefix = "L2"; break;
1988  case FiniteElement::INTEGRAL: prefix = "L2Int"; break;
1989  default:
1990  MFEM_ABORT("invalid map_type: " << map_type);
1991  }
1992  switch (btype)
1993  {
1995  snprintf(d_name, 32, "%s_%dD_P%d", prefix, dim, p);
1996  break;
1997  default:
1998  snprintf(d_name, 32, "%s_T%d_%dD_P%d", prefix, btype, dim, p);
1999  }
2000 
2001  for (int g = 0; g < Geometry::NumGeom; g++)
2002  {
2003  L2_Elements[g] = NULL;
2004  Tr_Elements[g] = NULL;
2005  }
2006  for (int i = 0; i < 2; i++)
2007  {
2008  SegDofOrd[i] = NULL;
2009  }
2010  for (int i = 0; i < 6; i++)
2011  {
2012  TriDofOrd[i] = NULL;
2013  }
2014  for (int i = 0; i < 24; i++)
2015  {
2016  TetDofOrd[i] = NULL;
2017  }
2018  OtherDofOrd = NULL;
2019 
2020  if (dim == 0)
2021  {
2022  L2_Elements[Geometry::POINT] = new PointFiniteElement;
2023  }
2024  else if (dim == 1)
2025  {
2026  if (b_type == BasisType::Positive)
2027  {
2028  L2_Elements[Geometry::SEGMENT] = new L2Pos_SegmentElement(p);
2029  }
2030  else
2031  {
2032  L2_Elements[Geometry::SEGMENT] = new L2_SegmentElement(p, btype);
2033  }
2034  L2_Elements[Geometry::SEGMENT]->SetMapType(map_type);
2035 
2036  Tr_Elements[Geometry::POINT] = new PointFiniteElement;
2037  // No need to set the map_type for Tr_Elements.
2038 
2039  const int pp1 = p + 1;
2040  SegDofOrd[0] = new int[2*pp1];
2041  SegDofOrd[1] = SegDofOrd[0] + pp1;
2042  for (int i = 0; i <= p; i++)
2043  {
2044  SegDofOrd[0][i] = i;
2045  SegDofOrd[1][i] = p - i;
2046  }
2047  }
2048  else if (dim == 2)
2049  {
2050  if (b_type == BasisType::Positive)
2051  {
2052  L2_Elements[Geometry::TRIANGLE] = new L2Pos_TriangleElement(p);
2053  L2_Elements[Geometry::SQUARE] = new L2Pos_QuadrilateralElement(p);
2054  }
2055  else
2056  {
2057  L2_Elements[Geometry::TRIANGLE] = new L2_TriangleElement(p, btype);
2058  L2_Elements[Geometry::SQUARE] = new L2_QuadrilateralElement(p, btype);
2059  }
2060  L2_Elements[Geometry::TRIANGLE]->SetMapType(map_type);
2061  L2_Elements[Geometry::SQUARE]->SetMapType(map_type);
2062  // Trace element use the default Gauss-Legendre nodal points for positive basis
2063  if (b_type == BasisType::Positive)
2064  {
2065  Tr_Elements[Geometry::SEGMENT] = new L2Pos_SegmentElement(p);
2066  }
2067  else
2068  {
2069  Tr_Elements[Geometry::SEGMENT] = new L2_SegmentElement(p, btype);
2070  }
2071 
2072  const int TriDof = L2_Elements[Geometry::TRIANGLE]->GetDof();
2073  TriDofOrd[0] = new int[6*TriDof];
2074  for (int i = 1; i < 6; i++)
2075  {
2076  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
2077  }
2078  const int pp1 = p + 1, pp2 = pp1 + 1;
2079  for (int j = 0; j <= p; j++)
2080  {
2081  for (int i = 0; i + j <= p; i++)
2082  {
2083  int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
2084  int k = p - j - i;
2085  TriDofOrd[0][o] = o; // (0,1,2)
2086  TriDofOrd[1][o] = TriDof - ((pp2-j)*(pp1-j))/2 + k; // (1,0,2)
2087  TriDofOrd[2][o] = TriDof - ((pp2-i)*(pp1-i))/2 + k; // (2,0,1)
2088  TriDofOrd[3][o] = TriDof - ((pp2-k)*(pp1-k))/2 + i; // (2,1,0)
2089  TriDofOrd[4][o] = TriDof - ((pp2-k)*(pp1-k))/2 + j; // (1,2,0)
2090  TriDofOrd[5][o] = TriDof - ((pp2-i)*(pp1-i))/2 + j; // (0,2,1)
2091  }
2092  }
2093  const int QuadDof = L2_Elements[Geometry::SQUARE]->GetDof();
2094  OtherDofOrd = new int[QuadDof];
2095  for (int j = 0; j < QuadDof; j++)
2096  {
2097  OtherDofOrd[j] = j; // for Or == 0
2098  }
2099  }
2100  else if (dim == 3)
2101  {
2102  if (b_type == BasisType::Positive)
2103  {
2104  L2_Elements[Geometry::TETRAHEDRON] = new L2Pos_TetrahedronElement(p);
2105  L2_Elements[Geometry::CUBE] = new L2Pos_HexahedronElement(p);
2106  L2_Elements[Geometry::PRISM] = new L2Pos_WedgeElement(p);
2107  }
2108  else
2109  {
2110  L2_Elements[Geometry::TETRAHEDRON] =
2111  new L2_TetrahedronElement(p, btype);
2112  L2_Elements[Geometry::CUBE] = new L2_HexahedronElement(p, btype);
2113  L2_Elements[Geometry::PRISM] = new L2_WedgeElement(p, btype);
2114  }
2115  L2_Elements[Geometry::PYRAMID] = new P0PyrFiniteElement;
2116 
2117  L2_Elements[Geometry::TETRAHEDRON]->SetMapType(map_type);
2118  L2_Elements[Geometry::CUBE]->SetMapType(map_type);
2119  L2_Elements[Geometry::PRISM]->SetMapType(map_type);
2120  L2_Elements[Geometry::PYRAMID]->SetMapType(map_type);
2121  // Trace element use the default Gauss-Legendre nodal points for positive basis
2122  if (b_type == BasisType::Positive)
2123  {
2124  Tr_Elements[Geometry::TRIANGLE] = new L2Pos_TriangleElement(p);
2125  Tr_Elements[Geometry::SQUARE] = new L2Pos_QuadrilateralElement(p);
2126  }
2127  else
2128  {
2129  Tr_Elements[Geometry::TRIANGLE] = new L2_TriangleElement(p, btype);
2130  Tr_Elements[Geometry::SQUARE] = new L2_QuadrilateralElement(p, btype);
2131  }
2132 
2133  const int TetDof = L2_Elements[Geometry::TETRAHEDRON]->GetDof();
2134  const int HexDof = L2_Elements[Geometry::CUBE]->GetDof();
2135  const int PriDof = L2_Elements[Geometry::PRISM]->GetDof();
2136  const int MaxDof = std::max(TetDof, std::max(PriDof, HexDof));
2137 
2138  TetDofOrd[0] = new int[24*TetDof];
2139  for (int i = 1; i < 24; i++)
2140  {
2141  TetDofOrd[i] = TetDofOrd[i-1] + TetDof;
2142  }
2143  // see Mesh::GetTetOrientation in mesh/mesh.cpp
2144  const int pp1 = p + 1, pp2 = pp1 + 1, pp3 = pp2 + 1;
2145  for (int k = 0; k <= p; k++)
2146  {
2147  for (int j = 0; j + k <= p; j++)
2148  {
2149  for (int i = 0; i + j + k <= p; i++)
2150  {
2151  int l = p - k - j - i;
2152  int o = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2153  + (j * (2 * p + 3 - j - 2 * k)) / 2 + i;
2154  int o1 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2155  + (k * (2 * p + 3 - k - 2 * j)) / 2 + i;
2156  int o2 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2157  + (k * (2 * p + 3 - k - 2 * i)) / 2 + j;
2158  int o3 = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2159  + (i * (2 * p + 3 - i - 2 * k)) / 2 + j;
2160  int o4 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2161  + (i * (2 * p + 3 - i - 2 * j)) / 2 + k;
2162  int o5 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2163  + (j * (2 * p + 3 - j - 2 * i)) / 2 + k;
2164  int o6 = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2165  + (l * (2 * p + 3 - l - 2 * k)) / 2 + j;
2166  int o7 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2167  + (k * (2 * p + 3 - k - 2 * l)) / 2 + j;
2168  int o8 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2169  + (j * (2 * p + 3 - j - 2 * l)) / 2 + k;
2170  int o9 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2171  + (l * (2 * p + 3 - l - 2 * j)) / 2 + k;
2172  int o10 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2173  + (k * (2 * p + 3 - k - 2 * j)) / 2 + l;
2174  int o11 = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2175  + (j * (2 * p + 3 - j - 2 * k)) / 2 + l;
2176  int o12 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2177  + (l * (2 * p + 3 - l - 2 * i)) / 2 + k;
2178  int o13 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2179  + (i * (2 * p + 3 - i - 2 * l)) / 2 + k;
2180  int o14 = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2181  + (i * (2 * p + 3 - i - 2 * k)) / 2 + l;
2182  int o15 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2183  + (k * (2 * p + 3 - k - 2 * i)) / 2 + l;
2184  int o16 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2185  + (k * (2 * p + 3 - k - 2 * l)) / 2 + i;
2186  int o17 = TetDof - ((pp1 - k) * (pp2 - k) * (pp3 - k)) / 6
2187  + (l * (2 * p + 3 - l - 2 * k)) / 2 + i;
2188  int o18 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2189  + (j * (2 * p + 3 - j - 2 * i)) / 2 + l;
2190  int o19 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2191  + (i * (2 * p + 3 - i - 2 * j)) / 2 + l;
2192  int o20 = TetDof - ((pp1 - j) * (pp2 - j) * (pp3 - j)) / 6
2193  + (l * (2 * p + 3 - l - 2 * j)) / 2 + i;
2194  int o21 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2195  + (j * (2 * p + 3 - j - 2 * l)) / 2 + i;
2196  int o22 = TetDof - ((pp1 - l) * (pp2 - l) * (pp3 - l)) / 6
2197  + (i * (2 * p + 3 - i - 2 * l)) / 2 + j;
2198  int o23 = TetDof - ((pp1 - i) * (pp2 - i) * (pp3 - i)) / 6
2199  + (l * (2 * p + 3 - l - 2 * i)) / 2 + j;
2200  TetDofOrd[ 0][o] = o; // (0,1,2,3)
2201  TetDofOrd[ 1][o] = o1; // (0,1,3,2)
2202  TetDofOrd[ 2][o] = o2; // (0,2,3,1)
2203  TetDofOrd[ 3][o] = o3; // (0,2,1,3)
2204  TetDofOrd[ 4][o] = o4; // (0,3,1,2)
2205  TetDofOrd[ 5][o] = o5; // (0,3,2,1)
2206  TetDofOrd[ 6][o] = o6; // (1,2,0,3)
2207  TetDofOrd[ 7][o] = o7; // (1,2,3,0)
2208  TetDofOrd[ 8][o] = o8; // (1,3,2,0)
2209  TetDofOrd[ 9][o] = o9; // (1,3,0,2)
2210  TetDofOrd[10][o] = o10; // (1,0,3,2)
2211  TetDofOrd[11][o] = o11; // (1,0,2,3)
2212  TetDofOrd[12][o] = o12; // (2,3,0,1)
2213  TetDofOrd[13][o] = o13; // (2,3,1,0)
2214  TetDofOrd[14][o] = o14; // (2,0,1,3)
2215  TetDofOrd[15][o] = o15; // (2,0,3,1)
2216  TetDofOrd[16][o] = o16; // (2,1,3,0)
2217  TetDofOrd[17][o] = o17; // (2,1,0,3)
2218  TetDofOrd[18][o] = o18; // (3,0,2,1)
2219  TetDofOrd[19][o] = o19; // (3,0,1,2)
2220  TetDofOrd[20][o] = o20; // (3,1,0,2)
2221  TetDofOrd[21][o] = o21; // (3,1,2,0)
2222  TetDofOrd[22][o] = o22; // (3,2,1,0)
2223  TetDofOrd[23][o] = o23; // (3,2,0,1)
2224  }
2225  }
2226  }
2227  OtherDofOrd = new int[MaxDof];
2228  for (int j = 0; j < MaxDof; j++)
2229  {
2230  OtherDofOrd[j] = j; // for Or == 0
2231  }
2232  }
2233  else
2234  {
2235  mfem::err << "L2_FECollection::L2_FECollection : dim = "
2236  << dim << endl;
2237  mfem_error();
2238  }
2239 }
2240 
2241 const FiniteElement *
2243 {
2244  if (GeomType != Geometry::PYRAMID || this->GetOrder() == 0)
2245  {
2246  return L2_Elements[GeomType];
2247  }
2248  else
2249  {
2250  MFEM_ABORT("L2 Pyramid basis functions are not yet supported "
2251  "for order > 0.");
2252  return NULL;
2253  }
2254 }
2255 
2257  int Or) const
2258 {
2259  switch (GeomType)
2260  {
2261  case Geometry::SEGMENT:
2262  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
2263 
2264  case Geometry::TRIANGLE:
2265  return TriDofOrd[Or%6];
2266 
2267  case Geometry::TETRAHEDRON:
2268  return TetDofOrd[Or%24];
2269 
2270  default:
2271  return (Or == 0) ? OtherDofOrd : NULL;
2272  }
2273 }
2274 
2276 {
2277  delete [] OtherDofOrd;
2278  delete [] SegDofOrd[0];
2279  delete [] TriDofOrd[0];
2280  delete [] TetDofOrd[0];
2281  for (int i = 0; i < Geometry::NumGeom; i++)
2282  {
2283  delete L2_Elements[i];
2284  delete Tr_Elements[i];
2285  }
2286 }
2287 
2288 
2289 RT_FECollection::RT_FECollection(const int order, const int dim,
2290  const int cb_type, const int ob_type)
2291  : FiniteElementCollection(order + 1)
2292  , dim(dim)
2293  , cb_type(cb_type)
2294  , ob_type(ob_type)
2295 {
2296  int p = order;
2297  MFEM_VERIFY(p >= 0, "RT_FECollection requires order >= 0.");
2298 
2299  int cp_type = BasisType::GetQuadrature1D(cb_type);
2300  int op_type = BasisType::GetQuadrature1D(ob_type);
2301 
2303  {
2304  const char *cb_name = BasisType::Name(cb_type); // this may abort
2305  MFEM_ABORT("unknown closed BasisType: " << cb_name);
2306  }
2308  ob_type != BasisType::IntegratedGLL)
2309  {
2310  const char *ob_name = BasisType::Name(ob_type); // this may abort
2311  MFEM_ABORT("unknown open BasisType: " << ob_name);
2312  }
2313 
2314  InitFaces(p, dim, FiniteElement::INTEGRAL, true);
2315 
2316  if (cb_type == BasisType::GaussLobatto &&
2317  ob_type == BasisType::GaussLegendre)
2318  {
2319  snprintf(rt_name, 32, "RT_%dD_P%d", dim, p);
2320  }
2321  else
2322  {
2323  snprintf(rt_name, 32, "RT@%c%c_%dD_P%d", (int)BasisType::GetChar(cb_type),
2324  (int)BasisType::GetChar(ob_type), dim, p);
2325  }
2326 
2327  const int pp1 = p + 1;
2328  if (dim == 2)
2329  {
2330  // TODO: cb_type, ob_type for triangles
2332  RT_dof[Geometry::TRIANGLE] = p*pp1;
2333 
2335  ob_type);
2336  // two vector components * n_unk_face *
2337  RT_dof[Geometry::SQUARE] = 2*p*pp1;
2338  }
2339  else if (dim == 3)
2340  {
2341  // TODO: cb_type, ob_type for tets
2343  RT_dof[Geometry::TETRAHEDRON] = p*pp1*(p + 2)/2;
2344 
2345  RT_Elements[Geometry::CUBE] = new RT_HexahedronElement(p, cb_type, ob_type);
2346  RT_dof[Geometry::CUBE] = 3*p*pp1*pp1;
2347 
2348  RT_Elements[Geometry::PRISM] = new RT_WedgeElement(p);
2349  RT_dof[Geometry::PRISM] = p*pp1*(3*p + 4)/2;
2350 
2351  RT_Elements[Geometry::PYRAMID] = new RT0PyrFiniteElement(false);
2353  }
2354  else
2355  {
2356  MFEM_ABORT("invalid dim = " << dim);
2357  }
2358 }
2359 
2360 // This is a special protected constructor only used by RT_Trace_FECollection
2361 // and DG_Interface_FECollection
2363  const int map_type, const bool signs,
2364  const int ob_type)
2365  : FiniteElementCollection(p + 1)
2366  , ob_type(ob_type)
2367 {
2370  {
2371  const char *ob_name = BasisType::Name(ob_type); // this may abort
2372  MFEM_ABORT("Invalid open basis type: " << ob_name);
2373  }
2374  InitFaces(p, dim, map_type, signs);
2375 }
2376 
2377 void RT_FECollection::InitFaces(const int p, const int dim_,
2378  const int map_type,
2379  const bool signs)
2380 {
2381  int op_type = BasisType::GetQuadrature1D(ob_type);
2382 
2383  MFEM_VERIFY(Quadrature1D::CheckOpen(op_type) != Quadrature1D::Invalid,
2384  "invalid open point type");
2385 
2386  const int pp1 = p + 1, pp2 = p + 2;
2387 
2388  for (int g = 0; g < Geometry::NumGeom; g++)
2389  {
2390  RT_Elements[g] = NULL;
2391  RT_dof[g] = 0;
2392  }
2393  // Degree of Freedom orderings
2394  for (int i = 0; i < 2; i++)
2395  {
2396  SegDofOrd[i] = NULL;
2397  }
2398  for (int i = 0; i < 6; i++)
2399  {
2400  TriDofOrd[i] = NULL;
2401  }
2402  for (int i = 0; i < 8; i++)
2403  {
2404  QuadDofOrd[i] = NULL;
2405  }
2406 
2407  if (dim_ == 2)
2408  {
2409  L2_SegmentElement *l2_seg = new L2_SegmentElement(p, ob_type);
2410  l2_seg->SetMapType(map_type);
2411  RT_Elements[Geometry::SEGMENT] = l2_seg;
2412  RT_dof[Geometry::SEGMENT] = pp1;
2413 
2414  SegDofOrd[0] = new int[2*pp1];
2415  SegDofOrd[1] = SegDofOrd[0] + pp1;
2416  for (int i = 0; i <= p; i++)
2417  {
2418  SegDofOrd[0][i] = i;
2419  SegDofOrd[1][i] = signs ? (-1 - (p - i)) : (p - i);
2420  }
2421  }
2422  else if (dim_ == 3)
2423  {
2424  L2_TriangleElement *l2_tri = new L2_TriangleElement(p, ob_type);
2425  l2_tri->SetMapType(map_type);
2426  RT_Elements[Geometry::TRIANGLE] = l2_tri;
2427  RT_dof[Geometry::TRIANGLE] = pp1*pp2/2;
2428 
2430  l2_quad->SetMapType(map_type);
2431  RT_Elements[Geometry::SQUARE] = l2_quad;
2432  RT_dof[Geometry::SQUARE] = pp1*pp1;
2433 
2434  int TriDof = RT_dof[Geometry::TRIANGLE];
2435  TriDofOrd[0] = new int[6*TriDof];
2436  for (int i = 1; i < 6; i++)
2437  {
2438  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
2439  }
2440  // see Mesh::GetTriOrientation in mesh/mesh.cpp,
2441  // the constructor of H1_FECollection
2442  for (int j = 0; j <= p; j++)
2443  {
2444  for (int i = 0; i + j <= p; i++)
2445  {
2446  int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
2447  int k = p - j - i;
2448  TriDofOrd[0][o] = o; // (0,1,2)
2449  TriDofOrd[1][o] = -1-(TriDof-((pp2-j)*(pp1-j))/2+k); // (1,0,2)
2450  TriDofOrd[2][o] = TriDof-((pp2-i)*(pp1-i))/2+k; // (2,0,1)
2451  TriDofOrd[3][o] = -1-(TriDof-((pp2-k)*(pp1-k))/2+i); // (2,1,0)
2452  TriDofOrd[4][o] = TriDof-((pp2-k)*(pp1-k))/2+j; // (1,2,0)
2453  TriDofOrd[5][o] = -1-(TriDof-((pp2-i)*(pp1-i))/2+j); // (0,2,1)
2454  if (!signs)
2455  {
2456  for (int kk = 1; kk < 6; kk += 2)
2457  {
2458  TriDofOrd[kk][o] = -1 - TriDofOrd[kk][o];
2459  }
2460  }
2461  }
2462  }
2463 
2464  int QuadDof = RT_dof[Geometry::SQUARE];
2465  QuadDofOrd[0] = new int[8*QuadDof];
2466  for (int i = 1; i < 8; i++)
2467  {
2468  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
2469  }
2470  // see Mesh::GetQuadOrientation in mesh/mesh.cpp
2471  for (int j = 0; j <= p; j++)
2472  {
2473  for (int i = 0; i <= p; i++)
2474  {
2475  int o = i + j*pp1;
2476  QuadDofOrd[0][o] = i + j*pp1; // (0,1,2,3)
2477  QuadDofOrd[1][o] = -1 - (j + i*pp1); // (0,3,2,1)
2478  QuadDofOrd[2][o] = j + (p - i)*pp1; // (1,2,3,0)
2479  QuadDofOrd[3][o] = -1 - ((p - i) + j*pp1); // (1,0,3,2)
2480  QuadDofOrd[4][o] = (p - i) + (p - j)*pp1; // (2,3,0,1)
2481  QuadDofOrd[5][o] = -1 - ((p - j) + (p - i)*pp1); // (2,1,0,3)
2482  QuadDofOrd[6][o] = (p - j) + i*pp1; // (3,0,1,2)
2483  QuadDofOrd[7][o] = -1 - (i + (p - j)*pp1); // (3,2,1,0)
2484  if (!signs)
2485  {
2486  for (int k = 1; k < 8; k += 2)
2487  {
2488  QuadDofOrd[k][o] = -1 - QuadDofOrd[k][o];
2489  }
2490  }
2491  }
2492  }
2493  }
2494 }
2495 
2496 const FiniteElement *
2498 {
2499  if (GeomType != Geometry::PYRAMID || this->GetOrder() == 1)
2500  {
2501  return RT_Elements[GeomType];
2502  }
2503  else
2504  {
2505  MFEM_ABORT("RT Pyramid basis functions are not yet supported "
2506  "for order > 0.");
2507  return NULL;
2508  }
2509 }
2510 
2512  int Or) const
2513 {
2514  if (GeomType == Geometry::SEGMENT)
2515  {
2516  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
2517  }
2518  else if (GeomType == Geometry::TRIANGLE)
2519  {
2520  return TriDofOrd[Or%6];
2521  }
2522  else if (GeomType == Geometry::SQUARE)
2523  {
2524  return QuadDofOrd[Or%8];
2525  }
2526  return NULL;
2527 }
2528 
2530 {
2531  int tr_dim, tr_p;
2532  if (!strncmp(rt_name, "RT_", 3))
2533  {
2534  tr_dim = atoi(rt_name + 3);
2535  tr_p = atoi(rt_name + 7);
2536  }
2537  else // rt_name = RT@.._.D_P*
2538  {
2539  tr_dim = atoi(rt_name + 6);
2540  tr_p = atoi(rt_name + 10);
2541  }
2542  return new RT_Trace_FECollection(tr_p, tr_dim, FiniteElement::INTEGRAL,
2543  ob_type);
2544 }
2545 
2547 {
2548  delete [] SegDofOrd[0];
2549  delete [] TriDofOrd[0];
2550  delete [] QuadDofOrd[0];
2551  for (int g = 0; g < Geometry::NumGeom; g++)
2552  {
2553  delete RT_Elements[g];
2554  }
2555 }
2556 
2557 
2559  const int map_type,
2560  const int ob_type)
2561  : RT_FECollection(p, dim, map_type, true, ob_type)
2562 {
2563  const char *prefix =
2564  (map_type == FiniteElement::INTEGRAL) ? "RT_Trace" : "RT_ValTrace";
2565  char ob_str[3] = { '\0', '\0', '\0' };
2566 
2567  if (ob_type != BasisType::GaussLegendre)
2568  {
2569  ob_str[0] = '@';
2570  ob_str[1] = BasisType::GetChar(ob_type);
2571  }
2572  snprintf(rt_name, 32, "%s%s_%dD_P%d", prefix, ob_str, dim, p);
2573 
2574  MFEM_VERIFY(dim == 2 || dim == 3, "Wrong dimension, dim = " << dim);
2575 }
2576 
2577 
2579  const int map_type,
2580  const int ob_type)
2581  : RT_FECollection(p, dim, map_type, false, ob_type)
2582 {
2583  MFEM_VERIFY(dim == 2 || dim == 3, "Wrong dimension, dim = " << dim);
2584 
2585  const char *prefix =
2586  (map_type == FiniteElement::VALUE) ? "DG_Iface" : "DG_IntIface";
2587  if (ob_type == BasisType::GaussLegendre)
2588  {
2589  snprintf(rt_name, 32, "%s_%dD_P%d", prefix, dim, p);
2590  }
2591  else
2592  {
2593  snprintf(rt_name, 32, "%s@%c_%dD_P%d", prefix,
2594  (int)BasisType::GetChar(ob_type), dim, p);
2595  }
2596 }
2597 
2599  const int cb_type, const int ob_type)
2600  : FiniteElementCollection(dim > 1 ? p : p - 1)
2601  , dim(dim)
2602  , cb_type(cb_type)
2603  , ob_type(ob_type)
2604 {
2605  MFEM_VERIFY(p >= 1, "ND_FECollection requires order >= 1.");
2606  MFEM_VERIFY(dim >= 1 && dim <= 3, "ND_FECollection requires 1 <= dim <= 3.");
2607 
2608  const int pm1 = p - 1, pm2 = p - 2;
2609 
2610  if (cb_type == BasisType::GaussLobatto &&
2611  ob_type == BasisType::GaussLegendre)
2612  {
2613  snprintf(nd_name, 32, "ND_%dD_P%d", dim, p);
2614  }
2615  else
2616  {
2617  snprintf(nd_name, 32, "ND@%c%c_%dD_P%d", (int)BasisType::GetChar(cb_type),
2618  (int)BasisType::GetChar(ob_type), dim, p);
2619  }
2620 
2621  for (int g = 0; g < Geometry::NumGeom; g++)
2622  {
2623  ND_Elements[g] = NULL;
2624  ND_dof[g] = 0;
2625  }
2626  for (int i = 0; i < 2; i++)
2627  {
2628  SegDofOrd[i] = NULL;
2629  }
2630  for (int i = 0; i < 6; i++)
2631  {
2632  TriDofOrd[i] = NULL;
2633  }
2634  for (int i = 0; i < 8; i++)
2635  {
2636  QuadDofOrd[i] = NULL;
2637  }
2638 
2639  int op_type = BasisType::GetQuadrature1D(ob_type);
2640  int cp_type = BasisType::GetQuadrature1D(cb_type);
2641 
2642  // Error checking
2644  ob_type != BasisType::IntegratedGLL)
2645  {
2646  const char *ob_name = BasisType::Name(ob_type);
2647  MFEM_ABORT("Invalid open basis point type: " << ob_name);
2648  }
2650  {
2651  const char *cb_name = BasisType::Name(cb_type);
2652  MFEM_ABORT("Invalid closed basis point type: " << cb_name);
2653  }
2654 
2655  if (dim >= 1)
2656  {
2657  ND_Elements[Geometry::SEGMENT] = new ND_SegmentElement(p, ob_type);
2659 
2660  SegDofOrd[0] = new int[2*p];
2661  SegDofOrd[1] = SegDofOrd[0] + p;
2662  for (int i = 0; i < p; i++)
2663  {
2664  SegDofOrd[0][i] = i;
2665  SegDofOrd[1][i] = -1 - (pm1 - i);
2666  }
2667  }
2668 
2669  if (dim >= 2)
2670  {
2672  ob_type);
2673  ND_dof[Geometry::SQUARE] = 2*p*pm1;
2674 
2675  // TODO: cb_type and ob_type for triangles
2677  ND_dof[Geometry::TRIANGLE] = p*pm1;
2678 
2679  int QuadDof = ND_dof[Geometry::SQUARE];
2680  QuadDofOrd[0] = new int[8*QuadDof];
2681  for (int i = 1; i < 8; i++)
2682  {
2683  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
2684  }
2685  // see Mesh::GetQuadOrientation in mesh/mesh.cpp
2686  for (int j = 0; j < pm1; j++)
2687  {
2688  for (int i = 0; i < p; i++)
2689  {
2690  int d1 = i + j*p; // x-component
2691  int d2 = p*pm1 + j + i*pm1; // y-component
2692  // (0,1,2,3)
2693  QuadDofOrd[0][d1] = d1;
2694  QuadDofOrd[0][d2] = d2;
2695  // (0,3,2,1)
2696  QuadDofOrd[1][d1] = d2;
2697  QuadDofOrd[1][d2] = d1;
2698  // (1,2,3,0)
2699  // QuadDofOrd[2][d1] = p*pm1 + (pm2 - j) + i*pm1;
2700  // QuadDofOrd[2][d2] = -1 - ((pm1 - i) + j*p);
2701  QuadDofOrd[2][d1] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2702  QuadDofOrd[2][d2] = i + (pm2 - j)*p;
2703  // (1,0,3,2)
2704  QuadDofOrd[3][d1] = -1 - ((pm1 - i) + j*p);
2705  QuadDofOrd[3][d2] = p*pm1 + (pm2 - j) + i*pm1;
2706  // (2,3,0,1)
2707  QuadDofOrd[4][d1] = -1 - ((pm1 - i) + (pm2 - j)*p);
2708  QuadDofOrd[4][d2] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2709  // (2,1,0,3)
2710  QuadDofOrd[5][d1] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2711  QuadDofOrd[5][d2] = -1 - ((pm1 - i) + (pm2 - j)*p);
2712  // (3,0,1,2)
2713  // QuadDofOrd[6][d1] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2714  // QuadDofOrd[6][d2] = i + (pm2 - j)*p;
2715  QuadDofOrd[6][d1] = p*pm1 + (pm2 - j) + i*pm1;
2716  QuadDofOrd[6][d2] = -1 - ((pm1 - i) + j*p);
2717  // (3,2,1,0)
2718  QuadDofOrd[7][d1] = i + (pm2 - j)*p;
2719  QuadDofOrd[7][d2] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2720  }
2721  }
2722 
2723  int TriDof = ND_dof[Geometry::TRIANGLE];
2724  TriDofOrd[0] = new int[6*TriDof];
2725  for (int i = 1; i < 6; i++)
2726  {
2727  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
2728  }
2729  // see Mesh::GetTriOrientation in mesh/mesh.cpp,
2730  // the constructor of H1_FECollection
2731  for (int j = 0; j <= pm2; j++)
2732  {
2733  for (int i = 0; i + j <= pm2; i++)
2734  {
2735  int k0 = p*pm1 - (p - j)*(pm1 - j) + 2*i;
2736  int k1 = 2*pm2 - 2*i + ((2*p-3)-j)*j;
2737  int k2 = 2*pm2 - 2*j + ((2*p-3)-i)*i;
2738  int k3 = p*pm1 - 2 - 3*j - i - (i+j)*(i+j);
2739  int k4 = p*pm1 - 2 - 3*i - j - (i+j)*(i+j);
2740  int k5 = p*pm1 - (p - i)*(pm1 - i) + 2*j;
2741 
2742  // (0,1,2)
2743  TriDofOrd[0][k0 ] = k0;
2744  TriDofOrd[0][k0+1] = k0 + 1;
2745  // (1,0,2)
2746  TriDofOrd[1][k0 ] = k1;
2747  TriDofOrd[1][k0+1] = k1 + 1;
2748  // (2,0,1)
2749  TriDofOrd[2][k0 ] = k2;
2750  TriDofOrd[2][k0+1] = k2 + 1;
2751  // (2,1,0)
2752  TriDofOrd[3][k0 ] = k3;
2753  TriDofOrd[3][k0+1] = k3 + 1;
2754  // (1,2,0)
2755  TriDofOrd[4][k0 ] = k4;
2756  TriDofOrd[4][k0+1] = k4 + 1;
2757  // (0,2,1)
2758  TriDofOrd[5][k0 ] = k5;
2759  TriDofOrd[5][k0+1] = k5 + 1;
2760  }
2761  }
2762  }
2763 
2764  if (dim >= 3)
2765  {
2766  ND_Elements[Geometry::CUBE] = new ND_HexahedronElement(p, cb_type, ob_type);
2767  ND_dof[Geometry::CUBE] = 3*p*pm1*pm1;
2768 
2769  // TODO: cb_type and ob_type for tets
2771  ND_dof[Geometry::TETRAHEDRON] = p*pm1*pm2/2;
2772 
2774  ND_dof[Geometry::PRISM] = p*pm1*(3*p-4)/2;
2775 
2778  }
2779 }
2780 
2781 const FiniteElement *
2783 {
2784  if (GeomType != Geometry::PYRAMID || this->GetOrder() == 1)
2785  {
2786  return ND_Elements[GeomType];
2787  }
2788  else
2789  {
2790  MFEM_ABORT("ND Pyramid basis functions are not yet supported "
2791  "for order > 1.");
2792  return NULL;
2793  }
2794 }
2795 
2797  int Or) const
2798 {
2799  if (GeomType == Geometry::SEGMENT)
2800  {
2801  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
2802  }
2803  else if (GeomType == Geometry::TRIANGLE)
2804  {
2805  return TriDofOrd[Or%6];
2806  }
2807  else if (GeomType == Geometry::SQUARE)
2808  {
2809  return QuadDofOrd[Or%8];
2810  }
2811  return NULL;
2812 }
2813 
2815 {
2816  int tr_p, tr_dim, tr_cb_type, tr_ob_type;
2817 
2818  tr_p = ND_dof[Geometry::SEGMENT];
2819  if (nd_name[2] == '_') // ND_
2820  {
2821  tr_dim = atoi(nd_name + 3);
2822  tr_cb_type = BasisType::GaussLobatto;
2823  tr_ob_type = BasisType::GaussLegendre;
2824  }
2825  else // ND@
2826  {
2827  tr_dim = atoi(nd_name + 6);
2828  tr_cb_type = BasisType::GetType(nd_name[3]);
2829  tr_ob_type = BasisType::GetType(nd_name[4]);
2830  }
2831  return new ND_Trace_FECollection(tr_p, tr_dim, tr_cb_type, tr_ob_type);
2832 }
2833 
2835 {
2836  delete [] SegDofOrd[0];
2837  delete [] TriDofOrd[0];
2838  delete [] QuadDofOrd[0];
2839  for (int g = 0; g < Geometry::NumGeom; g++)
2840  {
2841  delete ND_Elements[g];
2842  }
2843 }
2844 
2845 
2847  const int cb_type,
2848  const int ob_type)
2849  : ND_FECollection(p, dim-1, cb_type, ob_type)
2850 {
2851  if (cb_type == BasisType::GaussLobatto &&
2852  ob_type == BasisType::GaussLegendre)
2853  {
2854  snprintf(nd_name, 32, "ND_Trace_%dD_P%d", dim, p);
2855  }
2856  else
2857  {
2858  snprintf(nd_name, 32, "ND_Trace@%c%c_%dD_P%d",
2859  (int)BasisType::GetChar(cb_type),
2860  (int)BasisType::GetChar(ob_type), dim, p);
2861  }
2862 }
2863 
2864 
2866  const int cb_type, const int ob_type)
2868 {
2869  MFEM_VERIFY(p >= 1, "ND_R1D_FECollection requires order >= 1.");
2870  MFEM_VERIFY(dim == 1, "ND_R1D_FECollection requires dim == 1.");
2871 
2872  if (cb_type == BasisType::GaussLobatto &&
2873  ob_type == BasisType::GaussLegendre)
2874  {
2875  snprintf(nd_name, 32, "ND_R1D_%dD_P%d", dim, p);
2876  }
2877  else
2878  {
2879  snprintf(nd_name, 32, "ND_R1D@%c%c_%dD_P%d",
2880  (int)BasisType::GetChar(cb_type),
2881  (int)BasisType::GetChar(ob_type), dim, p);
2882  }
2883 
2884  for (int g = 0; g < Geometry::NumGeom; g++)
2885  {
2886  ND_Elements[g] = NULL;
2887  ND_dof[g] = 0;
2888  }
2889 
2890  int op_type = BasisType::GetQuadrature1D(ob_type);
2891  int cp_type = BasisType::GetQuadrature1D(cb_type);
2892 
2893  // Error checking
2895  {
2896  const char *ob_name = BasisType::Name(ob_type);
2897  MFEM_ABORT("Invalid open basis point type: " << ob_name);
2898  }
2900  {
2901  const char *cb_name = BasisType::Name(cb_type);
2902  MFEM_ABORT("Invalid closed basis point type: " << cb_name);
2903  }
2904 
2906  ND_dof[Geometry::POINT] = 2;
2907 
2909  cb_type,
2910  ob_type);
2911  ND_dof[Geometry::SEGMENT] = 3 * p - 2;
2912 }
2913 
2915  int Or) const
2916 {
2917  return NULL;
2918 }
2919 
2921 {
2922  return NULL;
2923 }
2924 
2926 {
2927  for (int g = 0; g < Geometry::NumGeom; g++)
2928  {
2929  delete ND_Elements[g];
2930  }
2931 }
2932 
2933 
2935  const int cb_type, const int ob_type)
2936  : FiniteElementCollection(p + 1)
2937 {
2938  MFEM_VERIFY(p >= 0, "RT_R1D_FECollection requires order >= 0.");
2939  MFEM_VERIFY(dim == 1, "RT_R1D_FECollection requires dim == 1.");
2940 
2941  if (cb_type == BasisType::GaussLobatto &&
2942  ob_type == BasisType::GaussLegendre)
2943  {
2944  snprintf(rt_name, 32, "RT_R1D_%dD_P%d", dim, p);
2945  }
2946  else
2947  {
2948  snprintf(rt_name, 32, "RT_R1D@%c%c_%dD_P%d",
2949  (int)BasisType::GetChar(cb_type),
2950  (int)BasisType::GetChar(ob_type), dim, p);
2951  }
2952 
2953  for (int g = 0; g < Geometry::NumGeom; g++)
2954  {
2955  RT_Elements[g] = NULL;
2956  RT_dof[g] = 0;
2957  }
2958 
2959  int op_type = BasisType::GetQuadrature1D(ob_type);
2960  int cp_type = BasisType::GetQuadrature1D(cb_type);
2961 
2962  // Error checking
2964  {
2965  const char *ob_name = BasisType::Name(ob_type);
2966  MFEM_ABORT("Invalid open basis point type: " << ob_name);
2967  }
2969  {
2970  const char *cb_name = BasisType::Name(cb_type);
2971  MFEM_ABORT("Invalid closed basis point type: " << cb_name);
2972  }
2973 
2975  RT_dof[Geometry::POINT] = 1;
2976 
2978  cb_type,
2979  ob_type);
2980  RT_dof[Geometry::SEGMENT] = 3 * p + 2;
2981 }
2982 
2984  int Or) const
2985 {
2986  return NULL;
2987 }
2988 
2990 {
2991  MFEM_ABORT("this method is not implemented in RT_R1D_FECollection!");
2992  return NULL;
2993 }
2994 
2996 {
2997  for (int g = 0; g < Geometry::NumGeom; g++)
2998  {
2999  delete RT_Elements[g];
3000  }
3001 }
3002 
3003 
3005  const int cb_type, const int ob_type)
3007 {
3008  MFEM_VERIFY(p >= 1, "ND_R2D_FECollection requires order >= 1.");
3009  MFEM_VERIFY(dim >= 1 && dim <= 2,
3010  "ND_R2D_FECollection requires 1 <= dim <= 2.");
3011 
3012  const int pm1 = p - 1, pm2 = p - 2;
3013 
3014  if (cb_type == BasisType::GaussLobatto &&
3015  ob_type == BasisType::GaussLegendre)
3016  {
3017  snprintf(nd_name, 32, "ND_R2D_%dD_P%d", dim, p);
3018  }
3019  else
3020  {
3021  snprintf(nd_name, 32, "ND_R2D@%c%c_%dD_P%d",
3022  (int)BasisType::GetChar(cb_type),
3023  (int)BasisType::GetChar(ob_type), dim, p);
3024  }
3025 
3026  for (int g = 0; g < Geometry::NumGeom; g++)
3027  {
3028  ND_Elements[g] = NULL;
3029  ND_dof[g] = 0;
3030  }
3031  for (int i = 0; i < 2; i++)
3032  {
3033  SegDofOrd[i] = NULL;
3034  }
3035 
3036  int op_type = BasisType::GetQuadrature1D(ob_type);
3037  int cp_type = BasisType::GetQuadrature1D(cb_type);
3038 
3039  // Error checking
3041  {
3042  const char *ob_name = BasisType::Name(ob_type);
3043  MFEM_ABORT("Invalid open basis point type: " << ob_name);
3044  }
3046  {
3047  const char *cb_name = BasisType::Name(cb_type);
3048  MFEM_ABORT("Invalid closed basis point type: " << cb_name);
3049  }
3050 
3051  ND_dof[Geometry::POINT] = 1;
3052 
3053  if (dim >= 1)
3054  {
3056  cb_type,
3057  ob_type);
3058  ND_dof[Geometry::SEGMENT] = 2 * p - 1;
3059 
3060  SegDofOrd[0] = new int[4 * p - 2];
3061  SegDofOrd[1] = SegDofOrd[0] + 2 * p - 1;
3062  for (int i = 0; i < p; i++)
3063  {
3064  SegDofOrd[0][i] = i;
3065  SegDofOrd[1][i] = -1 - (pm1 - i);
3066  }
3067  for (int i = 0; i < pm1; i++)
3068  {
3069  SegDofOrd[0][p+i] = p + i;
3070  SegDofOrd[1][p+i] = 2 * pm1 - i;
3071  }
3072  }
3073 
3074  if (dim >= 2)
3075  {
3077  cb_type,
3078  ob_type);
3079  ND_dof[Geometry::SQUARE] = 2*p*pm1 + pm1*pm1;
3080 
3081  // TODO: cb_type and ob_type for triangles
3083  ND_dof[Geometry::TRIANGLE] = p*pm1 + (pm1*pm2)/2;
3084  }
3085 }
3086 
3088  int Or) const
3089 {
3090  if (GeomType == Geometry::SEGMENT)
3091  {
3092  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
3093  }
3094  return NULL;
3095 }
3096 
3098 {
3099  int p, dim, cb_type, ob_type;
3100 
3102  if (nd_name[5] == '_') // ND_R2D_
3103  {
3104  dim = atoi(nd_name + 6);
3105  cb_type = BasisType::GaussLobatto;
3106  ob_type = BasisType::GaussLegendre;
3107  }
3108  else // ND_R2D@
3109  {
3110  dim = atoi(nd_name + 9);
3111  cb_type = BasisType::GetType(nd_name[6]);
3112  ob_type = BasisType::GetType(nd_name[7]);
3113  }
3114  return new ND_R2D_Trace_FECollection(p, dim, cb_type, ob_type);
3115 }
3116 
3118 {
3119  delete [] SegDofOrd[0];
3120  for (int g = 0; g < Geometry::NumGeom; g++)
3121  {
3122  delete ND_Elements[g];
3123  }
3124 }
3125 
3126 
3128  const int cb_type,
3129  const int ob_type)
3130  : ND_R2D_FECollection(p, dim-1, cb_type, ob_type)
3131 {
3132  if (cb_type == BasisType::GaussLobatto &&
3133  ob_type == BasisType::GaussLegendre)
3134  {
3135  snprintf(nd_name, 32, "ND_R2D_Trace_%dD_P%d", dim, p);
3136  }
3137  else
3138  {
3139  snprintf(nd_name, 32, "ND_R2D_Trace@%c%c_%dD_P%d",
3140  (int)BasisType::GetChar(cb_type),
3141  (int)BasisType::GetChar(ob_type), dim, p);
3142  }
3143 }
3144 
3145 
3147  const int cb_type, const int ob_type)
3148  : FiniteElementCollection(p + 1),
3149  ob_type(ob_type)
3150 {
3151  MFEM_VERIFY(p >= 0, "RT_R2D_FECollection requires order >= 0.");
3152  MFEM_VERIFY(dim >= 1 && dim <= 2,
3153  "RT_R2D_FECollection requires 1 <= dim <= 2.");
3154 
3155  int cp_type = BasisType::GetQuadrature1D(cb_type);
3156  int op_type = BasisType::GetQuadrature1D(ob_type);
3157 
3159  {
3160  const char *cb_name = BasisType::Name(cb_type); // this may abort
3161  MFEM_ABORT("unknown closed BasisType: " << cb_name);
3162  }
3164  {
3165  const char *ob_name = BasisType::Name(ob_type); // this may abort
3166  MFEM_ABORT("unknown open BasisType: " << ob_name);
3167  }
3168 
3169  InitFaces(p, dim, FiniteElement::INTEGRAL, true);
3170 
3171  if (cb_type == BasisType::GaussLobatto &&
3172  ob_type == BasisType::GaussLegendre)
3173  {
3174  snprintf(rt_name, 32, "RT_R2D_%dD_P%d", dim, p);
3175  }
3176  else
3177  {
3178  snprintf(rt_name, 32, "RT_R2D@%c%c_%dD_P%d",
3179  (int)BasisType::GetChar(cb_type),
3180  (int)BasisType::GetChar(ob_type), dim, p);
3181  }
3182 
3183  const int pp1 = p + 1;
3184  const int pp2 = p + 2;
3185  if (dim == 2)
3186  {
3187  // TODO: cb_type, ob_type for triangles
3189  RT_dof[Geometry::TRIANGLE] = p*pp1 + (pp1 * pp2) / 2;
3190 
3192  cb_type,
3193  ob_type);
3194  // two vector components * n_unk_face *
3195  RT_dof[Geometry::SQUARE] = 2*p*pp1 + pp1*pp1;
3196  }
3197 }
3198 
3199 // This is a special protected constructor only used by RT_Trace_FECollection
3200 // and DG_Interface_FECollection
3202  const int map_type,
3203  const bool signs, const int ob_type)
3204  : ob_type(ob_type)
3205 {
3208  {
3209  const char *ob_name = BasisType::Name(ob_type); // this may abort
3210  MFEM_ABORT("Invalid open basis type: " << ob_name);
3211  }
3212  InitFaces(p, dim, map_type, signs);
3213 }
3214 
3215 void RT_R2D_FECollection::InitFaces(const int p, const int dim,
3216  const int map_type,
3217  const bool signs)
3218 {
3219  int op_type = BasisType::GetQuadrature1D(ob_type);
3220 
3221  MFEM_VERIFY(Quadrature1D::CheckOpen(op_type) != Quadrature1D::Invalid,
3222  "invalid open point type");
3223 
3224  const int pp1 = p + 1;
3225 
3226  for (int g = 0; g < Geometry::NumGeom; g++)
3227  {
3228  RT_Elements[g] = NULL;
3229  RT_dof[g] = 0;
3230  }
3231  // Degree of Freedom orderings
3232  for (int i = 0; i < 2; i++)
3233  {
3234  SegDofOrd[i] = NULL;
3235  }
3236 
3237  if (dim == 2)
3238  {
3239  L2_SegmentElement *l2_seg = new L2_SegmentElement(p, ob_type);
3240  l2_seg->SetMapType(map_type);
3241  RT_Elements[Geometry::SEGMENT] = l2_seg;
3242  RT_dof[Geometry::SEGMENT] = pp1;
3243 
3244  SegDofOrd[0] = new int[2*pp1];
3245  SegDofOrd[1] = SegDofOrd[0] + pp1;
3246  for (int i = 0; i <= p; i++)
3247  {
3248  SegDofOrd[0][i] = i;
3249  SegDofOrd[1][i] = signs ? (-1 - (p - i)) : (p - i);
3250  }
3251  }
3252 }
3253 
3255  int Or) const
3256 {
3257  if (GeomType == Geometry::SEGMENT)
3258  {
3259  return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
3260  }
3261  return NULL;
3262 }
3263 
3265 {
3266  int dim, p;
3267  if (!strncmp(rt_name, "RT_R2D_", 7))
3268  {
3269  dim = atoi(rt_name + 7);
3270  p = atoi(rt_name + 11);
3271  }
3272  else // rt_name = RT_R2D@.._.D_P*
3273  {
3274  dim = atoi(rt_name + 10);
3275  p = atoi(rt_name + 14);
3276  }
3278 }
3279 
3281 {
3282  delete [] SegDofOrd[0];
3283  for (int g = 0; g < Geometry::NumGeom; g++)
3284  {
3285  delete RT_Elements[g];
3286  }
3287 }
3288 
3289 
3291  const int map_type,
3292  const int ob_type)
3293  : RT_R2D_FECollection(p, dim-1, map_type, true, ob_type)
3294 {
3295  const char *prefix =
3296  (map_type == FiniteElement::INTEGRAL) ? "RT_R2D_Trace" : "RT_R2D_ValTrace";
3297  char ob_str[3] = { '\0', '\0', '\0' };
3298 
3299  if (ob_type != BasisType::GaussLegendre)
3300  {
3301  ob_str[0] = '@';
3302  ob_str[1] = BasisType::GetChar(ob_type);
3303  }
3304  snprintf(rt_name, 32, "%s%s_%dD_P%d", prefix, ob_str, dim, p);
3305 
3306  MFEM_VERIFY(dim == 2, "Wrong dimension, dim = " << dim);
3307 }
3308 
3309 
3311 {
3312  snprintf(d_name, 32, "Local_%s", fe_name);
3313 
3314  Local_Element = NULL;
3315 
3316  if (!strcmp(fe_name, "BiCubic2DFiniteElement") ||
3317  !strcmp(fe_name, "Quad_Q3"))
3318  {
3319  GeomType = Geometry::SQUARE;
3320  Local_Element = new BiCubic2DFiniteElement;
3321  }
3322  else if (!strcmp(fe_name, "Nedelec1HexFiniteElement") ||
3323  !strcmp(fe_name, "Hex_ND1"))
3324  {
3325  GeomType = Geometry::CUBE;
3326  Local_Element = new Nedelec1HexFiniteElement;
3327  }
3328  else if (!strncmp(fe_name, "H1_", 3))
3329  {
3330  GeomType = Geometry::SQUARE;
3331  Local_Element = new H1_QuadrilateralElement(atoi(fe_name + 7));
3332  }
3333  else if (!strncmp(fe_name, "H1Pos_", 6))
3334  {
3335  GeomType = Geometry::SQUARE;
3336  Local_Element = new H1Pos_QuadrilateralElement(atoi(fe_name + 10));
3337  }
3338  else if (!strncmp(fe_name, "L2_", 3))
3339  {
3340  GeomType = Geometry::SQUARE;
3341  Local_Element = new L2_QuadrilateralElement(atoi(fe_name + 7));
3342  }
3343  else
3344  {
3345  mfem::err << "Local_FECollection::Local_FECollection : fe_name = "
3346  << fe_name << endl;
3347  mfem_error();
3348  }
3349 }
3350 
3351 
3353  : FiniteElementCollection((Order == VariableOrder) ? 1 : Order)
3354 {
3355  const int order = (Order == VariableOrder) ? 1 : Order;
3356  SegmentFE = new NURBS1DFiniteElement(order);
3357  QuadrilateralFE = new NURBS2DFiniteElement(order);
3358  ParallelepipedFE = new NURBS3DFiniteElement(order);
3359 
3360  SetOrder(Order);
3361 }
3362 
3363 void NURBSFECollection::SetOrder(int Order) const
3364 {
3365  mOrder = Order;
3366  if (Order != VariableOrder)
3367  {
3368  snprintf(name, 16, "NURBS%i", Order);
3369  }
3370  else
3371  {
3372  snprintf(name, 16, "NURBS");
3373  }
3374 }
3375 
3377 {
3378  delete ParallelepipedFE;
3379  delete QuadrilateralFE;
3380  delete SegmentFE;
3381 }
3382 
3383 const FiniteElement *
3385 {
3386  switch (GeomType)
3387  {
3388  case Geometry::SEGMENT: return SegmentFE;
3389  case Geometry::SQUARE: return QuadrilateralFE;
3390  case Geometry::CUBE: return ParallelepipedFE;
3391  default:
3392  mfem_error ("NURBSFECollection: unknown geometry type.");
3393  }
3394  return SegmentFE; // Make some compilers happy
3395 }
3396 
3398 {
3399  mfem_error("NURBSFECollection::DofForGeometry");
3400  return 0; // Make some compilers happy
3401 }
3402 
3404  int Or) const
3405 {
3406  mfem_error("NURBSFECollection::DofOrderForOrientation");
3407  return NULL;
3408 }
3409 
3411 {
3412  MFEM_ABORT("NURBS finite elements can not be statically condensed!");
3413  return NULL;
3414 }
3415 
3416 }
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:869
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:986
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:606
A 0D Nedelec finite element for the boundary of a 1D domain.
Definition: fe_nd.hpp:399
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1955
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1200
Arbitrary order L2 elements in 3D on a wedge.
Definition: fe_l2.hpp:140
Arbitrary order Nedelec elements in 1D on a segment.
Definition: fe_nd.hpp:279
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
Piecewise-cubic discontinuous finite elements in 2D. This class is kept only for backward compatibili...
Definition: fe_coll.hpp:1051
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1192
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:1026
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:349
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1062
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1008
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:422
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1241
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2782
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1157
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1273
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
Definition: fe_base.hpp:44
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:496
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a triangle.
Definition: fe_pos.hpp:312
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:45
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2242
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:3376
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2256
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:2362
Arbitrary order L2 elements in 2D on a triangle.
Definition: fe_l2.hpp:93
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:605
Arbitrary order Nedelec 3D elements in 2D on a square.
Definition: fe_nd.hpp:629
class LinearPyramidFiniteElement PyramidFE
Definition: fe.cpp:44
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2511
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3410
static const int NumGeom
Definition: geom.hpp:42
A 3D 1st order Nedelec element on a pyramid.
Arbitrary order L2 elements in 3D on a cube.
Definition: fe_l2.hpp:68
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2983
virtual ~L2_FECollection()
Definition: fe_coll.cpp:2275
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1862
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3384
Arbitrary order L2 elements in 1D utilizing the Bernstein basis on a segment.
Definition: fe_pos.hpp:258
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:790
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:647
T * GetData()
Returns the data.
Definition: array.hpp:112
Arbitrary order H1 elements in 2D on a square.
Definition: fe_h1.hpp:41
Piecewise-(bi/tri)linear continuous finite elements.
Definition: fe_coll.hpp:663
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1876
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:959
A 0D point finite element.
Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the interface between mesh elements...
Definition: fe_coll.hpp:548
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1416
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
Definition: fe_base.hpp:60
Arbitrary order Nedelec elements in 3D on a tetrahedron.
Definition: fe_nd.hpp:166
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition: fe_coll.hpp:518
Arbitrary order, three component, Raviart-Thomas elements in 1D on a segment.
Definition: fe_rt.hpp:331
virtual ~RT_R2D_FECollection()
Definition: fe_coll.cpp:3280
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2814
virtual void SetMapType(const int map_type_)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.cpp:2418
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1024
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:693
const int base_p
Order as returned by GetOrder().
Definition: fe_coll.hpp:209
Second order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:848
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2377
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2529
Arbitrary order, three component, Nedelec elements in 1D on a segment.
Definition: fe_nd.hpp:420
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:2578
RT_R2D_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3201
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe_base.hpp:91
A 2D bi-cubic element on a square with uniformly spaces nodes.
Arbitrary order Nedelec elements in 2D on a triangle.
Definition: fe_nd.hpp:224
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Definition: fe_base.hpp:110
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:951
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:199
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
Definition: fe_coll.cpp:3352
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1215
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1281
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:657
Arbitrary order Raviart-Thomas 3D elements in 2D on a triangle.
Definition: fe_rt.hpp:465
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:683
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:727
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2546
Arbitrary order Raviart-Thomas 3D elements in 2D on a square.
Definition: fe_rt.hpp:493
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1975
Class for standard nodal finite elements.
Definition: fe_base.hpp:706
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:743
PointFiniteElement PointFE
Definition: point.cpp:30
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2989
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:855
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1917
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:910
Arbitrary order L2 elements in 2D utilizing the Bernstein basis on a square.
Definition: fe_pos.hpp:276
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:569
int HasFaceDofs(Geometry::Type geom, int p) const
Definition: fe_coll.cpp:25
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.hpp:668
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:227
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
Arbitrary order L2 elements in 1D on a segment.
Definition: fe_l2.hpp:21
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, Geometry::Type &fg, int &fo, const int face_info)
Definition: fe_coll.cpp:352
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:348
static const int Dimension[NumGeom]
Definition: geom.hpp:47
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:937
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:322
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:468
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1402
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1392
Piecewise-linear discontinuous finite elements in 2D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:921
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:522
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1551
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:2846
Arbitrary order &quot;H^{1/2}-conforming&quot; trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:280
ND_R2D_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3127
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2598
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:1251
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces...
Definition: fe_coll.hpp:455
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:3087
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:1008
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2834
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2497
Bernstein polynomials.
Definition: fe_base.hpp:32
RT_R1D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2934
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:842
Array< FiniteElementCollection * > var_orders
Definition: fe_coll.hpp:216
Third order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:872
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:994
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1203
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:813
void SubDofOrder(Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:397
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1090
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1296
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:720
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:613
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1942
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1352
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1075
Arbitrary order H1 elements in 3D on a tetrahedron.
Definition: fe_h1.hpp:105
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1314
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:469
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1898
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:897
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:226
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3264
static const char * Name[NumGeom]
Definition: geom.hpp:45
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:797
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a wedge.
Definition: fe_pos.hpp:235
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:774
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:800
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3097
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:562
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1335
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1141
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1032
Arbitrary order L2 elements in 3D on a tetrahedron.
Definition: fe_l2.hpp:118
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:332
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1483
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
Definition: intrules.cpp:854
Arbitrary order Raviart-Thomas elements in 3D on a tetrahedron.
Definition: fe_rt.hpp:216
void SetOrder(int Order) const
Set the order and the name, based on the given Order: either a positive number for fixed order...
Definition: fe_coll.cpp:3363
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:896
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2796
Arbitrary order H1 elements in 3D utilizing the Bernstein basis on a cube.
Definition: fe_pos.hpp:159
ND_R1D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2865
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:776
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition: fe_coll.hpp:1127
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a wedge.
Definition: fe_pos.hpp:349
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:630
Arbitrary order H1 elements in 3D on a cube.
Definition: fe_h1.hpp:62
Serendipity basis (squares / cubes)
Definition: fe_base.hpp:36
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1322
class Linear3DFiniteElement TetrahedronFE
Definition: fe.cpp:36
const Array< int > & GetLexicographicOrdering() const
Get an Array&lt;int&gt; that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe_base.hpp:781
Arbitrary order Raviart-Thomas elements in 2D on a triangle.
Definition: fe_rt.hpp:156
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1449
Linear2DFiniteElement TriangleFE
Definition: fe.cpp:32
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Definition: intrules.cpp:866
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1376
Integrated GLL indicator functions.
Definition: fe_base.hpp:38
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:944
ND_R2D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3004
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:51
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
Definition: fe_base.hpp:103
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1359
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:296
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition: fe_coll.hpp:260
Lowest order Nedelec finite elements in 3D. This class is kept only for backward compatibility, consider using the new ND_FECollection instead.
Definition: fe_coll.hpp:1178
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1054
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:762
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1041
virtual const char * Name() const
Definition: fe_coll.hpp:65
Arbitrary order Nedelec elements in 2D on a square.
Definition: fe_nd.hpp:102
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1514
void InitVarOrder(int p) const
Definition: fe_coll.cpp:304
Arbitrary order H1 serendipity elements in 2D on a quad.
Definition: fe_ser.hpp:21
Arbitrary order Raviart-Thomas elements in 2D on a square.
Definition: fe_rt.hpp:23
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1529
class LinearWedgeFiniteElement WedgeFE
Definition: fe.cpp:40
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:3403
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2920
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:670
Arbitrary order H1 elements in 2D on a triangle.
Definition: fe_h1.hpp:83
int dim
Definition: ex24.cpp:53
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:3215
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:1101
RT_R2D_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3290
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1098
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1110
Arbitrary order L2 elements in 2D on a square.
Definition: fe_l2.hpp:39
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:523
Arbitrary order H1 elements in 1D.
Definition: fe_h1.hpp:21
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1255
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:827
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:563
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:710
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1465
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:2558
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:3254
Arbitrary order H1 elements in 1D utilizing the Bernstein basis.
Definition: fe_pos.hpp:117
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:587
Arbitrary order Nedelec 3D elements in 2D on a triangle.
Definition: fe_nd.hpp:598
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:691
First order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:824
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1231
Arbitrary order Nedelec elements in 3D on a cube.
Definition: fe_nd.hpp:22
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a square.
Definition: fe_pos.hpp:140
A 3D constant element on a pyramid.
Arbitrary order L2 elements in 3D utilizing the Bernstein basis on a cube.
Definition: fe_pos.hpp:294
virtual ~RT_R1D_FECollection()
Definition: fe_coll.cpp:2995
A linear element defined on a square pyramid.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:883
BiLinear2DFiniteElement QuadrilateralFE
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:415
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1075
Arbitrary order &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:395
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
Definition: fe_coll.hpp:557
virtual ~ND_R1D_FECollection()
Definition: fe_coll.cpp:2925
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
Second order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1229
A 3D 0th order Raviert-Thomas element on a pyramid.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1126
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3397
TriLinear3DFiniteElement HexahedronFE
Definition: hexahedron.cpp:52
Arbitrary order H1 elements in 3D on a wedge.
Definition: fe_h1.hpp:130
Piecewise-quadratic discontinuous finite elements in 2D. This class is kept only for backward compati...
Definition: fe_coll.hpp:986
Arbitrary order 3D &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh e...
Definition: fe_coll.hpp:597
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:972
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:925
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:3310
Arbitrary order Raviart-Thomas elements in 3D on a cube.
Definition: fe_rt.hpp:92
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1178
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:968
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1434
An arbitrary order 3D NURBS element on a cube.
Definition: fe_nurbs.hpp:115
An arbitrary order 1D NURBS element on a segment.
Definition: fe_nurbs.hpp:63
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1165
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:423
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:495
An arbitrary order 2D NURBS element on a square.
Definition: fe_nurbs.hpp:83
Arbitrary order H1 elements in 2D utilizing the Bernstein basis on a triangle.
Definition: fe_pos.hpp:178
A 3D 1st order Nedelec element on a cube.
virtual ~ND_R2D_FECollection()
Definition: fe_coll.cpp:3117
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1501
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:288
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2914