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