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