MFEM  v3.2
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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #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 
25 int FiniteElementCollection::HasFaceDofs(int GeomType) const
26 {
27  switch (GeomType)
28  {
29  case Geometry::TETRAHEDRON: return DofForGeometry (Geometry::TRIANGLE);
30  case Geometry::CUBE: return DofForGeometry (Geometry::SQUARE);
31  default:
32  mfem_error ("FiniteElementCollection::HasFaceDofs:"
33  " unknown geometry type.");
34  }
35  return 0;
36 }
37 
39 {
40  MFEM_ABORT("this method is not implemented in this derived class!");
41  return NULL;
42 }
43 
45 {
46  FiniteElementCollection *fec = NULL;
47 
48  if (!strcmp(name, "Linear"))
49  {
50  fec = new LinearFECollection;
51  }
52  else if (!strcmp(name, "Quadratic"))
53  {
54  fec = new QuadraticFECollection;
55  }
56  else if (!strcmp(name, "QuadraticPos"))
57  {
58  fec = new QuadraticPosFECollection;
59  }
60  else if (!strcmp(name, "Cubic"))
61  {
62  fec = new CubicFECollection;
63  }
64  else if (!strcmp(name, "Const3D"))
65  {
66  fec = new Const3DFECollection;
67  }
68  else if (!strcmp(name, "Const2D"))
69  {
70  fec = new Const2DFECollection;
71  }
72  else if (!strcmp(name, "LinearDiscont2D"))
73  {
75  }
76  else if (!strcmp(name, "GaussLinearDiscont2D"))
77  {
79  }
80  else if (!strcmp(name, "P1OnQuad"))
81  {
82  fec = new P1OnQuadFECollection;
83  }
84  else if (!strcmp(name, "QuadraticDiscont2D"))
85  {
87  }
88  else if (!strcmp(name, "QuadraticPosDiscont2D"))
89  {
91  }
92  else if (!strcmp(name, "GaussQuadraticDiscont2D"))
93  {
95  }
96  else if (!strcmp(name, "CubicDiscont2D"))
97  {
99  }
100  else if (!strcmp(name, "LinearDiscont3D"))
101  {
102  fec = new LinearDiscont3DFECollection;
103  }
104  else if (!strcmp(name, "QuadraticDiscont3D"))
105  {
107  }
108  else if (!strcmp(name, "LinearNonConf3D"))
109  {
110  fec = new LinearNonConf3DFECollection;
111  }
112  else if (!strcmp(name, "CrouzeixRaviart"))
113  {
114  fec = new CrouzeixRaviartFECollection;
115  }
116  else if (!strcmp(name, "ND1_3D"))
117  {
118  fec = new ND1_3DFECollection;
119  }
120  else if (!strcmp(name, "RT0_2D"))
121  {
122  fec = new RT0_2DFECollection;
123  }
124  else if (!strcmp(name, "RT1_2D"))
125  {
126  fec = new RT1_2DFECollection;
127  }
128  else if (!strcmp(name, "RT2_2D"))
129  {
130  fec = new RT2_2DFECollection;
131  }
132  else if (!strcmp(name, "RT0_3D"))
133  {
134  fec = new RT0_3DFECollection;
135  }
136  else if (!strcmp(name, "RT1_3D"))
137  {
138  fec = new RT1_3DFECollection;
139  }
140  else if (!strncmp(name, "H1_Trace_", 9))
141  {
142  fec = new H1_Trace_FECollection(atoi(name + 13), atoi(name + 9));
143  }
144  else if (!strncmp(name, "H1_", 3))
145  {
146  fec = new H1_FECollection(atoi(name + 7), atoi(name + 3));
147  }
148  else if (!strncmp(name, "H1Pos_Trace_", 12))
149  {
150  fec = new H1_Trace_FECollection(atoi(name + 16), atoi(name + 12), 1);
151  }
152  else if (!strncmp(name, "H1Pos_", 6))
153  {
154  fec = new H1Pos_FECollection(atoi(name + 10), atoi(name + 6));
155  }
156  else if (!strncmp(name, "L2_T", 4))
157  fec = new L2_FECollection(atoi(name + 10), atoi(name + 6),
158  atoi(name + 4));
159  else if (!strncmp(name, "L2_", 3))
160  {
161  fec = new L2_FECollection(atoi(name + 7), atoi(name + 3));
162  }
163  else if (!strncmp(name, "RT_Trace_", 9))
164  {
165  fec = new RT_Trace_FECollection(atoi(name + 13), atoi(name + 9));
166  }
167  else if (!strncmp(name, "RT_ValTrace_", 12))
168  {
169  fec = new RT_Trace_FECollection(atoi(name + 16), atoi(name + 12),
171  }
172  else if (!strncmp(name, "DG_Iface_", 9))
173  {
174  fec = new DG_Interface_FECollection(atoi(name + 13), atoi(name + 9));
175  }
176  else if (!strncmp(name, "DG_IntIface_", 12))
177  {
178  fec = new DG_Interface_FECollection(atoi(name + 16), atoi(name + 12),
180  }
181  else if (!strncmp(name, "RT_", 3))
182  {
183  fec = new RT_FECollection(atoi(name + 7), atoi(name + 3));
184  }
185  else if (!strncmp(name, "ND_Trace_", 9))
186  {
187  fec = new ND_Trace_FECollection(atoi(name + 13), atoi(name + 9));
188  }
189  else if (!strncmp(name, "ND_", 3))
190  {
191  fec = new ND_FECollection(atoi(name + 7), atoi(name + 3));
192  }
193  else if (!strncmp(name, "Local_", 6))
194  {
195  fec = new Local_FECollection(name + 6);
196  }
197  else if (!strncmp(name, "NURBS", 5))
198  {
199  fec = new NURBSFECollection(atoi(name + 5));
200  }
201  else
202  mfem_error ("FiniteElementCollection::New : "
203  "Unknown FiniteElementCollection!");
204 
205  return fec;
206 }
207 
208 template <Geometry::Type geom>
209 inline void FiniteElementCollection::GetNVE(int &nv, int &ne)
210 {
211  typedef typename Geometry::Constants<geom> g_consts;
212 
213  nv = g_consts::NumVert;
214  ne = g_consts::NumEdges;
215 }
216 
217 template <Geometry::Type geom, typename v_t>
218 inline void FiniteElementCollection::
219 GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
220 {
221  typedef typename Geometry::Constants<Geometry::SEGMENT> e_consts;
222  typedef typename Geometry::Constants<geom> g_consts;
223 
224  nv = e_consts::NumVert;
225  ne = 1;
226  e = edge_info/64;
227  eo = edge_info%64;
228  MFEM_ASSERT(0 <= e && e < g_consts::NumEdges, "");
229  MFEM_ASSERT(0 <= eo && eo < e_consts::NumOrient, "");
230  v[0] = g_consts::Edges[e][0];
231  v[1] = g_consts::Edges[e][1];
232  v[0] = e_consts::Orient[eo][v[0]];
233  v[1] = e_consts::Orient[eo][v[1]];
234 }
235 
236 template <Geometry::Type geom, Geometry::Type f_geom,
237  typename v_t, typename e_t, typename eo_t>
238 inline void FiniteElementCollection::
239 GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo,
240  int &nf, int &f, int &fg, int &fo, const int face_info)
241 {
242  typedef typename Geometry::Constants< geom> g_consts;
243  typedef typename Geometry::Constants<f_geom> f_consts;
244 
245  nv = f_consts::NumVert;
246  nf = 1;
247  f = face_info/64;
248  fg = f_geom;
249  fo = face_info%64;
250  MFEM_ASSERT(0 <= f && f < g_consts::NumFaces, "");
251  MFEM_ASSERT(0 <= fo && fo < f_consts::NumOrient, "");
252  for (int i = 0; i < f_consts::NumVert; i++)
253  {
254  v[i] = f_consts::Orient[fo][i];
255  v[i] = g_consts::FaceVert[f][v[i]];
256  }
257  ne = f_consts::NumEdges;
258  for (int i = 0; i < f_consts::NumEdges; i++)
259  {
260  int v0 = v[f_consts::Edges[i][0]];
261  int v1 = v[f_consts::Edges[i][1]];
262  int eor = 0;
263  if (v0 > v1) { swap(v0, v1); eor = 1; }
264  for (int j = g_consts::VertToVert::I[v0]; true; j++)
265  {
266  MFEM_ASSERT(j < g_consts::VertToVert::I[v0+1],
267  "internal error, edge not found");
268  if (v1 == g_consts::VertToVert::J[j][0])
269  {
270  int en = g_consts::VertToVert::J[j][1];
271  if (en < 0)
272  {
273  en = -1-en;
274  eor = 1-eor;
275  }
276  e[i] = en;
277  eo[i] = eor;
278  break;
279  }
280  }
281  }
282 }
283 
284 void FiniteElementCollection::SubDofOrder(int Geom, int SDim, int Info,
285  Array<int> &dofs) const
286 {
287  // Info = 64 * SubIndex + SubOrientation
288  MFEM_ASSERT(0 <= Geom && Geom < Geometry::NumGeom,
289  "invalid Geom = " << Geom);
290  const int Dim = Geometry::Dimension[Geom];
291  MFEM_ASSERT(0 <= SDim && SDim <= Dim, "invalid SDim = " << SDim
292  << " for Geom = " << Geometry::Name[Geom]);
293 
294  const int nvd = DofForGeometry(Geometry::POINT);
295  if (SDim == 0) // vertex
296  {
297  const int off = nvd*(Info/64);
298  dofs.SetSize(nvd);
299  for (int i = 0; i < nvd; i++)
300  {
301  dofs[i] = off + i;
302  }
303  }
304  else
305  {
306  int v[4], e[4], eo[4], f[1], fg[1], fo[1];
307  int av = 0, nv = 0, ae = 0, ne = 0, nf = 0;
308 
309  switch (Geom)
310  {
311  case Geometry::SEGMENT:
312  {
313  GetNVE<Geometry::SEGMENT>(av, ae);
314  GetEdge<Geometry::SEGMENT>(nv, v, ne, e[0], eo[0], Info);
315  break;
316  }
317 
318  case Geometry::TRIANGLE:
319  {
320  GetNVE<Geometry::TRIANGLE>(av, ae);
321  switch (SDim)
322  {
323  case 1:
324  GetEdge<Geometry::TRIANGLE>(nv, v, ne, e[0], eo[0], Info);
325  break;
326  case 2:
327  GetFace<Geometry::TRIANGLE,Geometry::TRIANGLE>(
328  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
329  break;
330  default:
331  goto not_supp;
332  }
333  break;
334  }
335 
336  case Geometry::SQUARE:
337  {
338  GetNVE<Geometry::SQUARE>(av, ae);
339  switch (SDim)
340  {
341  case 1:
342  GetEdge<Geometry::SQUARE>(nv, v, ne, e[0], eo[0], Info);
343  break;
344  case 2:
345  GetFace<Geometry::SQUARE,Geometry::SQUARE>(
346  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
347  break;
348  default:
349  goto not_supp;
350  }
351  break;
352  }
353 
355  {
356  GetNVE<Geometry::TETRAHEDRON>(av, ae);
357  switch (SDim)
358  {
359  case 1:
360  GetEdge<Geometry::TETRAHEDRON>(nv, v, ne, e[0], eo[0], Info);
361  break;
362  case 2:
363  GetFace<Geometry::TETRAHEDRON,Geometry::TRIANGLE>(
364  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
365  break;
366  default:
367  goto not_supp;
368  }
369  break;
370  }
371 
372  case Geometry::CUBE:
373  {
374  GetNVE<Geometry::CUBE>(av, ae);
375  switch (SDim)
376  {
377  case 1:
378  GetEdge<Geometry::CUBE>(nv, v, ne, e[0], eo[0], Info);
379  break;
380  case 2:
381  GetFace<Geometry::CUBE,Geometry::SQUARE>(
382  nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
383  break;
384  default:
385  goto not_supp;
386  }
387  break;
388  }
389 
390  default:
391  MFEM_ABORT("invalid Geom = " << Geom);
392  }
393 
394  int ned = (ne > 0) ? DofForGeometry(Geometry::SEGMENT) : 0;
395 
396  // add vertex dofs
397  dofs.SetSize(nv*nvd+ne*ned);
398  for (int i = 0; i < nv; i++)
399  {
400  for (int j = 0; j < nvd; j++)
401  {
402  dofs[i*nvd+j] = v[i]*nvd+j;
403  }
404  }
405  int l_off = nv*nvd, g_off = av*nvd;
406 
407  // add edge dofs
408  if (ned > 0)
409  {
410  for (int i = 0; i < ne; i++)
411  {
412  const int *ed = DofOrderForOrientation(Geometry::SEGMENT,
413  eo[i] ? -1 : 1);
414  for (int j = 0; j < ned; j++)
415  {
416  dofs[l_off+i*ned+j] =
417  ed[j] >= 0 ?
418  g_off+e[i]*ned+ed[j] :
419  -1-(g_off+e[i]*ned+(-1-ed[j]));
420  }
421  }
422  l_off += ne*ned;
423  g_off += ae*ned;
424  }
425 
426  // add face dofs
427  if (nf > 0)
428  {
429  const int nfd = DofForGeometry(fg[0]); // assume same face geometry
430  dofs.SetSize(dofs.Size()+nf*nfd);
431  for (int i = 0; i < nf; i++)
432  {
433  const int *fd = DofOrderForOrientation(fg[i], fo[i]);
434  for (int j = 0; j < nfd; j++)
435  {
436  dofs[l_off+i*nfd+j] =
437  fd[j] >= 0 ?
438  g_off+f[i]*nfd+fd[j] :
439  -1-(g_off+f[i]*nfd+(-1-fd[j]));
440  }
441  }
442  }
443 
444  // add volume dofs ...
445  }
446  return;
447 
448 not_supp:
449  MFEM_ABORT("Geom = " << Geometry::Name[Geom] <<
450  ", SDim = " << SDim << " is not supported");
451 }
452 
453 const FiniteElement *
455 {
456  switch (GeomType)
457  {
458  case Geometry::POINT: return &PointFE;
459  case Geometry::SEGMENT: return &SegmentFE;
460  case Geometry::TRIANGLE: return &TriangleFE;
461  case Geometry::SQUARE: return &QuadrilateralFE;
462  case Geometry::TETRAHEDRON: return &TetrahedronFE;
463  case Geometry::CUBE: return &ParallelepipedFE;
464  default:
465  mfem_error ("LinearFECollection: unknown geometry type.");
466  }
467  return &SegmentFE; // Make some compilers happy
468 }
469 
470 int LinearFECollection::DofForGeometry(int GeomType) const
471 {
472  switch (GeomType)
473  {
474  case Geometry::POINT: return 1;
475  case Geometry::SEGMENT: return 0;
476  case Geometry::TRIANGLE: return 0;
477  case Geometry::SQUARE: return 0;
478  case Geometry::TETRAHEDRON: return 0;
479  case Geometry::CUBE: return 0;
480  default:
481  mfem_error ("LinearFECollection: unknown geometry type.");
482  }
483  return 0; // Make some compilers happy
484 }
485 
486 int * LinearFECollection::DofOrderForOrientation(int GeomType, int Or) const
487 {
488  return NULL;
489 }
490 
491 
492 const FiniteElement *
494 {
495  switch (GeomType)
496  {
497  case Geometry::POINT: return &PointFE;
498  case Geometry::SEGMENT: return &SegmentFE;
499  case Geometry::TRIANGLE: return &TriangleFE;
500  case Geometry::SQUARE: return &QuadrilateralFE;
501  case Geometry::TETRAHEDRON: return &TetrahedronFE;
502  case Geometry::CUBE: return &ParallelepipedFE;
503  default:
504  mfem_error ("QuadraticFECollection: unknown geometry type.");
505  }
506  return &SegmentFE; // Make some compilers happy
507 }
508 
510 {
511  switch (GeomType)
512  {
513  case Geometry::POINT: return 1;
514  case Geometry::SEGMENT: return 1;
515  case Geometry::TRIANGLE: return 0;
516  case Geometry::SQUARE: return 1;
517  case Geometry::TETRAHEDRON: return 0;
518  case Geometry::CUBE: return 1;
519  default:
520  mfem_error ("QuadraticFECollection: unknown geometry type.");
521  }
522  return 0; // Make some compilers happy
523 }
524 
525 int * QuadraticFECollection::DofOrderForOrientation(int GeomType, int Or) const
526 {
527  static int indexes[] = { 0 };
528 
529  return indexes;
530 }
531 
532 
533 const FiniteElement *
535 {
536  switch (GeomType)
537  {
538  case Geometry::SEGMENT: return &SegmentFE;
539  case Geometry::SQUARE: return &QuadrilateralFE;
540  default:
541  mfem_error ("QuadraticPosFECollection: unknown geometry type.");
542  }
543  return NULL; // Make some compilers happy
544 }
545 
547 {
548  switch (GeomType)
549  {
550  case Geometry::POINT: return 1;
551  case Geometry::SEGMENT: return 1;
552  case Geometry::SQUARE: return 1;
553  default:
554  mfem_error ("QuadraticPosFECollection: unknown geometry type.");
555  }
556  return 0; // Make some compilers happy
557 }
558 
560 const
561 {
562  static int indexes[] = { 0 };
563 
564  return indexes;
565 }
566 
567 
568 const FiniteElement *
570 {
571  switch (GeomType)
572  {
573  case Geometry::POINT: return &PointFE;
574  case Geometry::SEGMENT: return &SegmentFE;
575  case Geometry::TRIANGLE: return &TriangleFE;
576  case Geometry::SQUARE: return &QuadrilateralFE;
577  case Geometry::TETRAHEDRON: return &TetrahedronFE;
578  case Geometry::CUBE: return &ParallelepipedFE;
579  default:
580  mfem_error ("CubicFECollection: unknown geometry type.");
581  }
582  return &SegmentFE; // Make some compilers happy
583 }
584 
585 int CubicFECollection::DofForGeometry(int GeomType) const
586 {
587  switch (GeomType)
588  {
589  case Geometry::POINT: return 1;
590  case Geometry::SEGMENT: return 2;
591  case Geometry::TRIANGLE: return 1;
592  case Geometry::SQUARE: return 4;
593  case Geometry::TETRAHEDRON: return 0;
594  case Geometry::CUBE: return 8;
595  default:
596  mfem_error ("CubicFECollection: unknown geometry type.");
597  }
598  return 0; // Make some compilers happy
599 }
600 
601 int * CubicFECollection::DofOrderForOrientation(int GeomType, int Or) const
602 {
603  if (GeomType == Geometry::SEGMENT)
604  {
605  static int ind_pos[] = { 0, 1 };
606  static int ind_neg[] = { 1, 0 };
607 
608  if (Or < 0)
609  {
610  return ind_neg;
611  }
612  return ind_pos;
613  }
614  else if (GeomType == Geometry::TRIANGLE)
615  {
616  static int indexes[] = { 0 };
617 
618  return indexes;
619  }
620  else if (GeomType == Geometry::SQUARE)
621  {
622  static int sq_ind[8][4] = {{0, 1, 2, 3}, {0, 2, 1, 3},
623  {2, 0, 3, 1}, {1, 0, 3, 2},
624  {3, 2, 1, 0}, {3, 1, 2, 0},
625  {1, 3, 0, 2}, {2, 3, 0, 1}
626  };
627  return sq_ind[Or];
628  }
629 
630  return NULL;
631 }
632 
633 
634 const FiniteElement *
636 {
637  switch (GeomType)
638  {
639  case Geometry::SEGMENT: return &SegmentFE;
640  case Geometry::TRIANGLE: return &TriangleFE;
641  case Geometry::SQUARE: return &QuadrilateralFE;
642  default:
643  mfem_error ("CrouzeixRaviartFECollection: unknown geometry type.");
644  }
645  return &SegmentFE; // Make some compilers happy
646 }
647 
649 {
650  switch (GeomType)
651  {
652  case Geometry::POINT: return 0;
653  case Geometry::SEGMENT: return 1;
654  case Geometry::TRIANGLE: return 0;
655  case Geometry::SQUARE: return 0;
656  default:
657  mfem_error ("CrouzeixRaviartFECollection: unknown geometry type.");
658  }
659  return 0; // Make some compilers happy
660 }
661 
663 const
664 {
665  static int indexes[] = { 0 };
666 
667  return indexes;
668 }
669 
670 
671 const FiniteElement *
673 {
674  switch (GeomType)
675  {
676  case Geometry::SEGMENT: return &SegmentFE;
677  case Geometry::TRIANGLE: return &TriangleFE;
678  case Geometry::SQUARE: return &QuadrilateralFE;
679  default:
680  mfem_error ("RT0_2DFECollection: unknown geometry type.");
681  }
682  return &SegmentFE; // Make some compilers happy
683 }
684 
685 int RT0_2DFECollection::DofForGeometry(int GeomType) const
686 {
687  switch (GeomType)
688  {
689  case Geometry::POINT: return 0;
690  case Geometry::SEGMENT: return 1;
691  case Geometry::TRIANGLE: return 0;
692  case Geometry::SQUARE: return 0;
693  default:
694  mfem_error ("RT0_2DFECollection: unknown geometry type.");
695  }
696  return 0; // Make some compilers happy
697 }
698 
700 const
701 {
702  static int ind_pos[] = { 0 };
703  static int ind_neg[] = { -1 };
704 
705  if (Or > 0)
706  {
707  return ind_pos;
708  }
709  return ind_neg;
710 }
711 
712 
713 const FiniteElement *
715 {
716  switch (GeomType)
717  {
718  case Geometry::SEGMENT: return &SegmentFE;
719  case Geometry::TRIANGLE: return &TriangleFE;
720  case Geometry::SQUARE: return &QuadrilateralFE;
721  default:
722  mfem_error ("RT1_2DFECollection: unknown geometry type.");
723  }
724  return &SegmentFE; // Make some compilers happy
725 }
726 
727 int RT1_2DFECollection::DofForGeometry(int GeomType) const
728 {
729  switch (GeomType)
730  {
731  case Geometry::POINT: return 0;
732  case Geometry::SEGMENT: return 2;
733  case Geometry::TRIANGLE: return 2;
734  case Geometry::SQUARE: return 4;
735  default:
736  mfem_error ("RT1_2DFECollection: unknown geometry type.");
737  }
738  return 0; // Make some compilers happy
739 }
740 
742 const
743 {
744  static int ind_pos[] = { 0, 1 };
745  static int ind_neg[] = { -2, -1 };
746 
747  if (Or > 0)
748  {
749  return ind_pos;
750  }
751  return ind_neg;
752 }
753 
754 const FiniteElement *
756 {
757  switch (GeomType)
758  {
759  case Geometry::SEGMENT: return &SegmentFE;
760  case Geometry::TRIANGLE: return &TriangleFE;
761  case Geometry::SQUARE: return &QuadrilateralFE;
762  default:
763  mfem_error ("RT2_2DFECollection: unknown geometry type.");
764  }
765  return &SegmentFE; // Make some compilers happy
766 }
767 
768 int RT2_2DFECollection::DofForGeometry(int GeomType) const
769 {
770  switch (GeomType)
771  {
772  case Geometry::POINT: return 0;
773  case Geometry::SEGMENT: return 3;
774  case Geometry::TRIANGLE: return 6;
775  case Geometry::SQUARE: return 12;
776  default:
777  mfem_error ("RT2_2DFECollection: unknown geometry type.");
778  }
779  return 0; // Make some compilers happy
780 }
781 
783 const
784 {
785  static int ind_pos[] = { 0, 1, 2 };
786  static int ind_neg[] = { -3, -2, -1 };
787 
788  if (Or > 0)
789  {
790  return ind_pos;
791  }
792  return ind_neg;
793 }
794 
795 
796 const FiniteElement *
798 {
799  switch (GeomType)
800  {
801  case Geometry::TRIANGLE: return &TriangleFE;
802  case Geometry::SQUARE: return &QuadrilateralFE;
803  default:
804  mfem_error ("Const2DFECollection: unknown geometry type.");
805  }
806  return &TriangleFE; // Make some compilers happy
807 }
808 
809 int Const2DFECollection::DofForGeometry(int GeomType) const
810 {
811  switch (GeomType)
812  {
813  case Geometry::POINT: return 0;
814  case Geometry::SEGMENT: return 0;
815  case Geometry::TRIANGLE: return 1;
816  case Geometry::SQUARE: return 1;
817  default:
818  mfem_error ("Const2DFECollection: unknown geometry type.");
819  }
820  return 0; // Make some compilers happy
821 }
822 
824 const
825 {
826  return NULL;
827 }
828 
829 
830 const FiniteElement *
832 {
833  switch (GeomType)
834  {
835  case Geometry::TRIANGLE: return &TriangleFE;
836  case Geometry::SQUARE: return &QuadrilateralFE;
837  default:
838  mfem_error ("LinearDiscont2DFECollection: unknown geometry type.");
839  }
840  return &TriangleFE; // Make some compilers happy
841 }
842 
844 {
845  switch (GeomType)
846  {
847  case Geometry::POINT: return 0;
848  case Geometry::SEGMENT: return 0;
849  case Geometry::TRIANGLE: return 3;
850  case Geometry::SQUARE: return 4;
851  default:
852  mfem_error ("LinearDiscont2DFECollection: unknown geometry type.");
853  }
854  return 0; // Make some compilers happy
855 }
856 
858 const
859 {
860  return NULL;
861 }
862 
863 
864 const FiniteElement *
866 {
867  switch (GeomType)
868  {
869  case Geometry::TRIANGLE: return &TriangleFE;
870  case Geometry::SQUARE: return &QuadrilateralFE;
871  default:
872  mfem_error ("GaussLinearDiscont2DFECollection:"
873  " unknown geometry type.");
874  }
875  return &TriangleFE; // Make some compilers happy
876 }
877 
879 {
880  switch (GeomType)
881  {
882  case Geometry::POINT: return 0;
883  case Geometry::SEGMENT: return 0;
884  case Geometry::TRIANGLE: return 3;
885  case Geometry::SQUARE: return 4;
886  default:
887  mfem_error ("GaussLinearDiscont2DFECollection:"
888  " unknown geometry type.");
889  }
890  return 0; // Make some compilers happy
891 }
892 
894  int GeomType, int Or) const
895 {
896  return NULL;
897 }
898 
899 
900 const FiniteElement *
902 {
903  if (GeomType != Geometry::SQUARE)
904  {
905  mfem_error ("P1OnQuadFECollection: unknown geometry type.");
906  }
907  return &QuadrilateralFE;
908 }
909 
911 {
912  switch (GeomType)
913  {
914  case Geometry::POINT: return 0;
915  case Geometry::SEGMENT: return 0;
916  case Geometry::SQUARE: return 3;
917  default:
918  mfem_error ("P1OnQuadFECollection: unknown geometry type.");
919  }
920  return 0; // Make some compilers happy
921 }
922 
924  int GeomType, int Or) const
925 {
926  return NULL;
927 }
928 
929 
930 const FiniteElement *
932 {
933  switch (GeomType)
934  {
935  case Geometry::TRIANGLE: return &TriangleFE;
936  case Geometry::SQUARE: return &QuadrilateralFE;
937  default:
938  mfem_error ("QuadraticDiscont2DFECollection: unknown geometry type.");
939  }
940  return &TriangleFE; // Make some compilers happy
941 }
942 
944 {
945  switch (GeomType)
946  {
947  case Geometry::POINT: return 0;
948  case Geometry::SEGMENT: return 0;
949  case Geometry::TRIANGLE: return 6;
950  case Geometry::SQUARE: return 9;
951  default:
952  mfem_error ("QuadraticDiscont2DFECollection: unknown geometry type.");
953  }
954  return 0; // Make some compilers happy
955 }
956 
958  int GeomType, int Or) const
959 {
960  return NULL;
961 }
962 
963 
964 const FiniteElement *
966 {
967  switch (GeomType)
968  {
969  case Geometry::SQUARE: return &QuadrilateralFE;
970  default:
971  mfem_error ("QuadraticPosDiscont2DFECollection: unknown geometry type.");
972  }
973  return NULL; // Make some compilers happy
974 }
975 
977 {
978  switch (GeomType)
979  {
980  case Geometry::POINT: return 0;
981  case Geometry::SEGMENT: return 0;
982  case Geometry::SQUARE: return 9;
983  default:
984  mfem_error ("QuadraticPosDiscont2DFECollection: unknown geometry type.");
985  }
986  return 0; // Make some compilers happy
987 }
988 
989 
990 const FiniteElement *
992 const
993 {
994  switch (GeomType)
995  {
996  case Geometry::TRIANGLE: return &TriangleFE;
997  case Geometry::SQUARE: return &QuadrilateralFE;
998  default:
999  mfem_error ("GaussQuadraticDiscont2DFECollection:"
1000  " unknown geometry type.");
1001  }
1002  return &QuadrilateralFE; // Make some compilers happy
1003 }
1004 
1006 {
1007  switch (GeomType)
1008  {
1009  case Geometry::POINT: return 0;
1010  case Geometry::SEGMENT: return 0;
1011  case Geometry::TRIANGLE: return 6;
1012  case Geometry::SQUARE: return 9;
1013  default:
1014  mfem_error ("GaussQuadraticDiscont2DFECollection:"
1015  " unknown geometry type.");
1016  }
1017  return 0; // Make some compilers happy
1018 }
1019 
1021  int GeomType, int Or) const
1022 {
1023  return NULL;
1024 }
1025 
1026 
1027 const FiniteElement *
1029 {
1030  switch (GeomType)
1031  {
1032  case Geometry::TRIANGLE: return &TriangleFE;
1033  case Geometry::SQUARE: return &QuadrilateralFE;
1034  default:
1035  mfem_error ("CubicDiscont2DFECollection: unknown geometry type.");
1036  }
1037  return &TriangleFE; // Make some compilers happy
1038 }
1039 
1041 {
1042  switch (GeomType)
1043  {
1044  case Geometry::POINT: return 0;
1045  case Geometry::SEGMENT: return 0;
1046  case Geometry::TRIANGLE: return 10;
1047  case Geometry::SQUARE: return 16;
1048  default:
1049  mfem_error ("CubicDiscont2DFECollection: unknown geometry type.");
1050  }
1051  return 0; // Make some compilers happy
1052 }
1053 
1055 const
1056 {
1057  return NULL;
1058 }
1059 
1060 
1061 const FiniteElement *
1063 {
1064  switch (GeomType)
1065  {
1066  case Geometry::TRIANGLE: return &TriangleFE;
1067  case Geometry::SQUARE: return &QuadrilateralFE;
1068  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1069  case Geometry::CUBE: return &ParallelepipedFE;
1070  default:
1071  mfem_error ("LinearNonConf3DFECollection: unknown geometry type.");
1072  }
1073  return &TriangleFE; // Make some compilers happy
1074 }
1075 
1077 {
1078  switch (GeomType)
1079  {
1080  case Geometry::POINT: return 0;
1081  case Geometry::SEGMENT: return 0;
1082  case Geometry::TRIANGLE: return 1;
1083  case Geometry::SQUARE: return 1;
1084  case Geometry::TETRAHEDRON: return 0;
1085  case Geometry::CUBE: return 0;
1086  default:
1087  mfem_error ("LinearNonConf3DFECollection: unknown geometry type.");
1088  }
1089  return 0; // Make some compilers happy
1090 }
1091 
1093 const
1094 {
1095  static int indexes[] = { 0 };
1096 
1097  return indexes;
1098 }
1099 
1100 
1101 const FiniteElement *
1103 {
1104  switch (GeomType)
1105  {
1106  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1107  case Geometry::CUBE: return &ParallelepipedFE;
1108  default:
1109  mfem_error ("Const3DFECollection: unknown geometry type.");
1110  }
1111  return &TetrahedronFE; // Make some compilers happy
1112 }
1113 
1115 {
1116  switch (GeomType)
1117  {
1118  case Geometry::POINT: return 0;
1119  case Geometry::SEGMENT: return 0;
1120  case Geometry::TRIANGLE: return 0;
1121  case Geometry::TETRAHEDRON: return 1;
1122  case Geometry::SQUARE: return 0;
1123  case Geometry::CUBE: return 1;
1124  default:
1125  mfem_error ("Const3DFECollection: unknown geometry type.");
1126  }
1127  return 0; // Make some compilers happy
1128 }
1129 
1131 const
1132 {
1133  return NULL;
1134 }
1135 
1136 
1137 const FiniteElement *
1139 {
1140  switch (GeomType)
1141  {
1142  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1143  case Geometry::CUBE: return &ParallelepipedFE;
1144  default:
1145  mfem_error ("LinearDiscont3DFECollection: unknown geometry type.");
1146  }
1147  return &TetrahedronFE; // Make some compilers happy
1148 }
1149 
1151 {
1152  switch (GeomType)
1153  {
1154  case Geometry::POINT: return 0;
1155  case Geometry::SEGMENT: return 0;
1156  case Geometry::TRIANGLE: return 0;
1157  case Geometry::SQUARE: return 0;
1158  case Geometry::TETRAHEDRON: return 4;
1159  case Geometry::CUBE: return 8;
1160  default:
1161  mfem_error ("LinearDiscont3DFECollection: unknown geometry type.");
1162  }
1163  return 0; // Make some compilers happy
1164 }
1165 
1167 const
1168 {
1169  return NULL;
1170 }
1171 
1172 
1173 const FiniteElement *
1175 {
1176  switch (GeomType)
1177  {
1178  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1179  case Geometry::CUBE: return &ParallelepipedFE;
1180  default:
1181  mfem_error ("QuadraticDiscont3DFECollection: unknown geometry type.");
1182  }
1183  return &TetrahedronFE; // Make some compilers happy
1184 }
1185 
1187 {
1188  switch (GeomType)
1189  {
1190  case Geometry::POINT: return 0;
1191  case Geometry::SEGMENT: return 0;
1192  case Geometry::TRIANGLE: return 0;
1193  case Geometry::SQUARE: return 0;
1194  case Geometry::TETRAHEDRON: return 10;
1195  case Geometry::CUBE: return 27;
1196  default:
1197  mfem_error ("QuadraticDiscont3DFECollection: unknown geometry type.");
1198  }
1199  return 0; // Make some compilers happy
1200 }
1201 
1203  int GeomType, int Or) const
1204 {
1205  return NULL;
1206 }
1207 
1208 const FiniteElement *
1210 {
1211  switch (GeomType)
1212  {
1213  case Geometry::POINT: return &PointFE;
1214  case Geometry::SEGMENT: return &SegmentFE;
1215  case Geometry::TRIANGLE: return &TriangleFE;
1216  case Geometry::SQUARE: return &QuadrilateralFE;
1217  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1218  case Geometry::CUBE: return &ParallelepipedFE;
1219  default:
1220  mfem_error ("RefinedLinearFECollection: unknown geometry type.");
1221  }
1222  return &SegmentFE; // Make some compilers happy
1223 }
1224 
1226 {
1227  switch (GeomType)
1228  {
1229  case Geometry::POINT: return 1;
1230  case Geometry::SEGMENT: return 1;
1231  case Geometry::TRIANGLE: return 0;
1232  case Geometry::SQUARE: return 1;
1233  case Geometry::TETRAHEDRON: return 0;
1234  case Geometry::CUBE: return 1;
1235  default:
1236  mfem_error ("RefinedLinearFECollection: unknown geometry type.");
1237  }
1238  return 0; // Make some compilers happy
1239 }
1240 
1242  int Or) const
1243 {
1244  static int indexes[] = { 0 };
1245 
1246  return indexes;
1247 }
1248 
1249 
1250 const FiniteElement *
1252 {
1253  switch (GeomType)
1254  {
1255  case Geometry::CUBE: return &HexahedronFE;
1256  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1257  default:
1258  mfem_error ("ND1_3DFECollection: unknown geometry type.");
1259  }
1260  return &HexahedronFE; // Make some compilers happy
1261 }
1262 
1263 int ND1_3DFECollection::DofForGeometry(int GeomType) const
1264 {
1265  switch (GeomType)
1266  {
1267  case Geometry::POINT: return 0;
1268  case Geometry::SEGMENT: return 1;
1269  case Geometry::TRIANGLE: return 0;
1270  case Geometry::SQUARE: return 0;
1271  case Geometry::TETRAHEDRON: return 0;
1272  case Geometry::CUBE: return 0;
1273  default:
1274  mfem_error ("ND1_3DFECollection: unknown geometry type.");
1275  }
1276  return 0; // Make some compilers happy
1277 }
1278 
1280 const
1281 {
1282  static int ind_pos[] = { 0 };
1283  static int ind_neg[] = { -1 };
1284 
1285  if (Or > 0)
1286  {
1287  return ind_pos;
1288  }
1289  return ind_neg;
1290 }
1291 
1292 
1293 const FiniteElement *
1295 {
1296  switch (GeomType)
1297  {
1298  case Geometry::TRIANGLE: return &TriangleFE;
1299  case Geometry::SQUARE: return &QuadrilateralFE;
1300  case Geometry::CUBE: return &HexahedronFE;
1301  case Geometry::TETRAHEDRON: return &TetrahedronFE;
1302  default:
1303  mfem_error ("RT0_3DFECollection: unknown geometry type.");
1304  }
1305  return &HexahedronFE; // Make some compilers happy
1306 }
1307 
1308 int RT0_3DFECollection::DofForGeometry(int GeomType) const
1309 {
1310  switch (GeomType)
1311  {
1312  case Geometry::POINT: return 0;
1313  case Geometry::SEGMENT: return 0;
1314  case Geometry::TRIANGLE: return 1;
1315  case Geometry::SQUARE: return 1;
1316  case Geometry::TETRAHEDRON: return 0;
1317  case Geometry::CUBE: return 0;
1318  default:
1319  mfem_error ("RT0_3DFECollection: unknown geometry type.");
1320  }
1321  return 0; // Make some compilers happy
1322 }
1323 
1325 const
1326 {
1327  static int ind_pos[] = { 0 };
1328  static int ind_neg[] = { -1 };
1329 
1330  if ((GeomType == Geometry::TRIANGLE) || (GeomType == Geometry::SQUARE))
1331  {
1332  if (Or % 2 == 0)
1333  {
1334  return ind_pos;
1335  }
1336  return ind_neg;
1337  }
1338  return NULL;
1339 }
1340 
1341 const FiniteElement *
1343 {
1344  switch (GeomType)
1345  {
1346  case Geometry::TRIANGLE: return &TriangleFE;
1347  case Geometry::SQUARE: return &QuadrilateralFE;
1348  case Geometry::CUBE: return &HexahedronFE;
1349  default:
1350  mfem_error ("RT1_3DFECollection: unknown geometry type.");
1351  }
1352  return &HexahedronFE; // Make some compilers happy
1353 }
1354 
1355 int RT1_3DFECollection::DofForGeometry(int GeomType) const
1356 {
1357  switch (GeomType)
1358  {
1359  case Geometry::POINT: return 0;
1360  case Geometry::SEGMENT: return 0;
1361  case Geometry::TRIANGLE: return 2;
1362  case Geometry::SQUARE: return 4;
1363  case Geometry::CUBE: return 12;
1364  default:
1365  mfem_error ("RT1_3DFECollection: unknown geometry type.");
1366  }
1367  return 0; // Make some compilers happy
1368 }
1369 
1371 const
1372 {
1373  if (GeomType == Geometry::SQUARE)
1374  {
1375  static int sq_ind[8][4] =
1376  {
1377  {0, 1, 2, 3}, {-1, -3, -2, -4},
1378  {2, 0, 3, 1}, {-2, -1, -4, -3},
1379  {3, 2, 1, 0}, {-4, -2, -3, -1},
1380  {1, 3, 0, 2}, {-3, -4, -1, -2}
1381  };
1382 
1383  return sq_ind[Or];
1384  }
1385  else
1386  {
1387  return NULL;
1388  }
1389 }
1390 
1391 
1392 H1_FECollection::H1_FECollection(const int p, const int dim, const int type)
1393 {
1394  const int pm1 = p - 1, pm2 = pm1 - 1, pm3 = pm2 - 1;
1395 
1396  m_type = (BasisType)type;
1397  if (type == GaussLobatto)
1398  {
1399  snprintf(h1_name, 32, "H1_%dD_P%d", dim, p);
1400  }
1401  else
1402  {
1403  snprintf(h1_name, 32, "H1Pos_%dD_P%d", dim, p);
1404  }
1405 
1406  for (int g = 0; g < Geometry::NumGeom; g++)
1407  {
1408  H1_dof[g] = 0;
1409  H1_Elements[g] = NULL;
1410  }
1411  for (int i = 0; i < 2; i++)
1412  {
1413  SegDofOrd[i] = NULL;
1414  }
1415  for (int i = 0; i < 6; i++)
1416  {
1417  TriDofOrd[i] = NULL;
1418  }
1419  for (int i = 0; i < 8; i++)
1420  {
1421  QuadDofOrd[i] = NULL;
1422  }
1423 
1424  H1_dof[Geometry::POINT] = 1;
1425  H1_Elements[Geometry::POINT] = new PointFiniteElement;
1426 
1427  if (dim >= 1)
1428  {
1429  H1_dof[Geometry::SEGMENT] = pm1;
1430  if (type == GaussLobatto)
1431  {
1432  H1_Elements[Geometry::SEGMENT] = new H1_SegmentElement(p);
1433  }
1434  else
1435  {
1436  H1_Elements[Geometry::SEGMENT] = new H1Pos_SegmentElement(p);
1437  }
1438 
1439  SegDofOrd[0] = new int[2*pm1];
1440  SegDofOrd[1] = SegDofOrd[0] + pm1;
1441  for (int i = 0; i < pm1; i++)
1442  {
1443  SegDofOrd[0][i] = i;
1444  SegDofOrd[1][i] = pm2 - i;
1445  }
1446  }
1447 
1448  if (dim >= 2)
1449  {
1450  H1_dof[Geometry::TRIANGLE] = (pm1*pm2)/2;
1451  H1_dof[Geometry::SQUARE] = pm1*pm1;
1452  if (type == GaussLobatto)
1453  {
1454  H1_Elements[Geometry::TRIANGLE] = new H1_TriangleElement(p);
1455  H1_Elements[Geometry::SQUARE] = new H1_QuadrilateralElement(p);
1456  }
1457  else
1458  {
1459  H1_Elements[Geometry::TRIANGLE] = new H1Pos_TriangleElement(p);
1460  H1_Elements[Geometry::SQUARE] = new H1Pos_QuadrilateralElement(p);
1461  }
1462 
1463  const int &TriDof = H1_dof[Geometry::TRIANGLE];
1464  const int &QuadDof = H1_dof[Geometry::SQUARE];
1465  TriDofOrd[0] = new int[6*TriDof];
1466  for (int i = 1; i < 6; i++)
1467  {
1468  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1469  }
1470  // see Mesh::GetTriOrientation in mesh/mesh.cpp
1471  for (int j = 0; j < pm2; j++)
1472  for (int i = 0; i + j < pm2; i++)
1473  {
1474  int o = TriDof - ((pm1 - j)*(pm2 - j))/2 + i;
1475  int k = pm3 - j - i;
1476  TriDofOrd[0][o] = o; // (0,1,2)
1477  TriDofOrd[1][o] = TriDof - ((pm1-j)*(pm2-j))/2 + k; // (1,0,2)
1478  TriDofOrd[2][o] = TriDof - ((pm1-i)*(pm2-i))/2 + k; // (2,0,1)
1479  TriDofOrd[3][o] = TriDof - ((pm1-k)*(pm2-k))/2 + i; // (2,1,0)
1480  TriDofOrd[4][o] = TriDof - ((pm1-k)*(pm2-k))/2 + j; // (1,2,0)
1481  TriDofOrd[5][o] = TriDof - ((pm1-i)*(pm2-i))/2 + j; // (0,2,1)
1482  }
1483 
1484  QuadDofOrd[0] = new int[8*QuadDof];
1485  for (int i = 1; i < 8; i++)
1486  {
1487  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
1488  }
1489  // see Mesh::GetQuadOrientation in mesh/mesh.cpp
1490  for (int j = 0; j < pm1; j++)
1491  for (int i = 0; i < pm1; i++)
1492  {
1493  int o = i + j*pm1;
1494  QuadDofOrd[0][o] = i + j*pm1; // (0,1,2,3)
1495  QuadDofOrd[1][o] = j + i*pm1; // (0,3,2,1)
1496  QuadDofOrd[2][o] = j + (pm2 - i)*pm1; // (1,2,3,0)
1497  QuadDofOrd[3][o] = (pm2 - i) + j*pm1; // (1,0,3,2)
1498  QuadDofOrd[4][o] = (pm2 - i) + (pm2 - j)*pm1; // (2,3,0,1)
1499  QuadDofOrd[5][o] = (pm2 - j) + (pm2 - i)*pm1; // (2,1,0,3)
1500  QuadDofOrd[6][o] = (pm2 - j) + i*pm1; // (3,0,1,2)
1501  QuadDofOrd[7][o] = i + (pm2 - j)*pm1; // (3,2,1,0)
1502  }
1503 
1504  if (dim >= 3)
1505  {
1506  H1_dof[Geometry::TETRAHEDRON] = (TriDof*pm3)/3;
1507  H1_dof[Geometry::CUBE] = QuadDof*pm1;
1508  if (type == GaussLobatto)
1509  {
1510  H1_Elements[Geometry::TETRAHEDRON] = new H1_TetrahedronElement(p);
1511  H1_Elements[Geometry::CUBE] = new H1_HexahedronElement(p);
1512  }
1513  else
1514  {
1515  H1_Elements[Geometry::TETRAHEDRON] = new H1Pos_TetrahedronElement(p);
1516  H1_Elements[Geometry::CUBE] = new H1Pos_HexahedronElement(p);
1517  }
1518  }
1519  }
1520 }
1521 
1522 int *H1_FECollection::DofOrderForOrientation(int GeomType, int Or) const
1523 {
1524  if (GeomType == Geometry::SEGMENT)
1525  {
1526  if (Or > 0)
1527  {
1528  return SegDofOrd[0];
1529  }
1530  return SegDofOrd[1];
1531  }
1532  else if (GeomType == Geometry::TRIANGLE)
1533  {
1534  return TriDofOrd[Or%6];
1535  }
1536  else if (GeomType == Geometry::SQUARE)
1537  {
1538  return QuadDofOrd[Or%8];
1539  }
1540  return NULL;
1541 }
1542 
1544 {
1545  int p = H1_dof[Geometry::SEGMENT] + 1;
1546  if (!strncmp(h1_name, "H1_", 3))
1547  {
1548  return new H1_Trace_FECollection(p, atoi(h1_name + 3));
1549  }
1550  else if (!strncmp(h1_name, "H1Pos_", 6))
1551  {
1552  return new H1_Trace_FECollection(p, atoi(h1_name + 6), 1);
1553  }
1554  return NULL;
1555 }
1556 
1558 {
1559  delete [] SegDofOrd[0];
1560  delete [] TriDofOrd[0];
1561  delete [] QuadDofOrd[0];
1562  for (int g = 0; g < Geometry::NumGeom; g++)
1563  {
1564  delete H1_Elements[g];
1565  }
1566 }
1567 
1568 
1570  const int type)
1571  : H1_FECollection(p, dim-1, type)
1572 {
1573  if (type == 0)
1574  {
1575  snprintf(h1_name, 32, "H1_Trace_%dD_P%d", dim, p);
1576  }
1577  else
1578  {
1579  snprintf(h1_name, 32, "H1Pos_Trace_%dD_P%d", dim, p);
1580  }
1581 }
1582 
1583 
1584 L2_FECollection::L2_FECollection(const int p, const int dim, const int type)
1585 {
1586  m_type = (BasisType)type;
1587  if (type == 0)
1588  {
1589  snprintf(d_name, 32, "L2_%dD_P%d", dim, p);
1590  }
1591  else
1592  {
1593  snprintf(d_name, 32, "L2_T%d_%dD_P%d", type, dim, p);
1594  }
1595 
1596  for (int g = 0; g < Geometry::NumGeom; g++)
1597  {
1598  L2_Elements[g] = NULL;
1599  Tr_Elements[g] = NULL;
1600  }
1601  for (int i = 0; i < 2; i++)
1602  {
1603  SegDofOrd[i] = NULL;
1604  }
1605  for (int i = 0; i < 6; i++)
1606  {
1607  TriDofOrd[i] = NULL;
1608  }
1609 
1610  if (dim == 1)
1611  {
1612  if (type == 0 || type == 1)
1613  {
1614  L2_Elements[Geometry::SEGMENT] = new L2_SegmentElement(p, type);
1615  }
1616  else
1617  {
1618  L2_Elements[Geometry::SEGMENT] = new L2Pos_SegmentElement(p);
1619  }
1620 
1621  Tr_Elements[Geometry::POINT] = new PointFiniteElement;
1622 
1623  const int pp1 = p + 1;
1624  SegDofOrd[0] = new int[2*pp1];
1625  SegDofOrd[1] = SegDofOrd[0] + pp1;
1626  for (int i = 0; i <= p; i++)
1627  {
1628  SegDofOrd[0][i] = i;
1629  SegDofOrd[1][i] = p - i;
1630  }
1631  }
1632  else if (dim == 2)
1633  {
1634  if (type == 0 || type == 1)
1635  {
1636  L2_Elements[Geometry::TRIANGLE] = new L2_TriangleElement(p, type);
1637  L2_Elements[Geometry::SQUARE] = new L2_QuadrilateralElement(p, type);
1638  }
1639  else
1640  {
1641  L2_Elements[Geometry::TRIANGLE] = new L2Pos_TriangleElement(p);
1642  L2_Elements[Geometry::SQUARE] = new L2Pos_QuadrilateralElement(p);
1643  }
1644  Tr_Elements[Geometry::SEGMENT] = new L2_SegmentElement(p, 0);
1645 
1646  const int TriDof = L2_Elements[Geometry::TRIANGLE]->GetDof();
1647  TriDofOrd[0] = new int[6*TriDof];
1648  for (int i = 1; i < 6; i++)
1649  {
1650  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1651  }
1652  const int pp1 = p + 1, pp2 = pp1 + 1;
1653  for (int j = 0; j <= p; j++)
1654  for (int i = 0; i + j <= p; i++)
1655  {
1656  int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
1657  int k = p - j - i;
1658  TriDofOrd[0][o] = o; // (0,1,2)
1659  TriDofOrd[1][o] = TriDof - ((pp2-j)*(pp1-j))/2 + k; // (1,0,2)
1660  TriDofOrd[2][o] = TriDof - ((pp2-i)*(pp1-i))/2 + k; // (2,0,1)
1661  TriDofOrd[3][o] = TriDof - ((pp2-k)*(pp1-k))/2 + i; // (2,1,0)
1662  TriDofOrd[4][o] = TriDof - ((pp2-k)*(pp1-k))/2 + j; // (1,2,0)
1663  TriDofOrd[5][o] = TriDof - ((pp2-i)*(pp1-i))/2 + j; // (0,2,1)
1664  }
1665  }
1666  else if (dim == 3)
1667  {
1668  if (type == 0 || type == 1)
1669  {
1670  L2_Elements[Geometry::TETRAHEDRON] =
1671  new L2_TetrahedronElement(p, type);
1672  L2_Elements[Geometry::CUBE] = new L2_HexahedronElement(p, type);
1673  }
1674  else
1675  {
1676  L2_Elements[Geometry::TETRAHEDRON] = new L2Pos_TetrahedronElement(p);
1677  L2_Elements[Geometry::CUBE] = new L2Pos_HexahedronElement(p);
1678  }
1679  Tr_Elements[Geometry::TRIANGLE] = new L2_TriangleElement(p, 0);
1680  Tr_Elements[Geometry::SQUARE] = new L2_QuadrilateralElement(p, 0);
1681  }
1682  else
1683  {
1684  cerr << "L2_FECollection::L2_FECollection : dim = "
1685  << dim << endl;
1686  mfem_error();
1687  }
1688 }
1689 
1690 int *L2_FECollection::DofOrderForOrientation(int GeomType, int Or) const
1691 {
1692  if (GeomType == Geometry::SEGMENT)
1693  {
1694  if (Or > 0)
1695  {
1696  return SegDofOrd[0];
1697  }
1698  return SegDofOrd[1];
1699  }
1700  else if (GeomType == Geometry::TRIANGLE)
1701  {
1702  return TriDofOrd[Or%6];
1703  }
1704  return NULL;
1705 }
1706 
1708 {
1709  delete [] SegDofOrd[0];
1710  delete [] TriDofOrd[0];
1711  for (int i = 0; i < Geometry::NumGeom; i++)
1712  {
1713  delete L2_Elements[i];
1714  delete Tr_Elements[i];
1715  }
1716 }
1717 
1718 
1719 RT_FECollection::RT_FECollection(const int p, const int dim)
1720 {
1721  InitFaces(p, dim, FiniteElement::INTEGRAL, true);
1722 
1723  snprintf(rt_name, 32, "RT_%dD_P%d", dim, p);
1724 
1725  const int pp1 = p + 1;
1726  if (dim == 2)
1727  {
1729  RT_dof[Geometry::TRIANGLE] = p*pp1;
1730 
1732  RT_dof[Geometry::SQUARE] = 2*p*pp1;
1733  }
1734  else if (dim == 3)
1735  {
1737  RT_dof[Geometry::TETRAHEDRON] = p*pp1*(p + 2)/2;
1738 
1740  RT_dof[Geometry::CUBE] = 3*p*pp1*pp1;
1741  }
1742  else
1743  {
1744  cerr << "RT_FECollection::RT_FECollection : dim = " << dim << endl;
1745  mfem_error();
1746  }
1747 }
1748 
1749 void RT_FECollection::InitFaces(const int p, const int dim, const int map_type,
1750  const bool signs)
1751 {
1752  const int pp1 = p + 1, pp2 = p + 2;
1753 
1754  for (int g = 0; g < Geometry::NumGeom; g++)
1755  {
1756  RT_Elements[g] = NULL;
1757  RT_dof[g] = 0;
1758  }
1759  for (int i = 0; i < 2; i++)
1760  {
1761  SegDofOrd[i] = NULL;
1762  }
1763  for (int i = 0; i < 6; i++)
1764  {
1765  TriDofOrd[i] = NULL;
1766  }
1767  for (int i = 0; i < 8; i++)
1768  {
1769  QuadDofOrd[i] = NULL;
1770  }
1771 
1772  if (dim == 2)
1773  {
1774  L2_SegmentElement *l2_seg = new L2_SegmentElement(p);
1775  l2_seg->SetMapType(map_type);
1776  RT_Elements[Geometry::SEGMENT] = l2_seg;
1777  RT_dof[Geometry::SEGMENT] = pp1;
1778 
1779  SegDofOrd[0] = new int[2*pp1];
1780  SegDofOrd[1] = SegDofOrd[0] + pp1;
1781  for (int i = 0; i <= p; i++)
1782  {
1783  SegDofOrd[0][i] = i;
1784  SegDofOrd[1][i] = signs ? (-1 - (p - i)) : (p - i);
1785  }
1786  }
1787  else if (dim == 3)
1788  {
1789  L2_TriangleElement *l2_tri = new L2_TriangleElement(p);
1790  l2_tri->SetMapType(map_type);
1791  RT_Elements[Geometry::TRIANGLE] = l2_tri;
1792  RT_dof[Geometry::TRIANGLE] = pp1*pp2/2;
1793 
1795  l2_quad->SetMapType(map_type);
1796  RT_Elements[Geometry::SQUARE] = l2_quad;
1797  RT_dof[Geometry::SQUARE] = pp1*pp1;
1798 
1799  int TriDof = RT_dof[Geometry::TRIANGLE];
1800  TriDofOrd[0] = new int[6*TriDof];
1801  for (int i = 1; i < 6; i++)
1802  {
1803  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1804  }
1805  // see Mesh::GetTriOrientation in mesh/mesh.cpp,
1806  // the constructor of H1_FECollection
1807  for (int j = 0; j <= p; j++)
1808  for (int i = 0; i + j <= p; i++)
1809  {
1810  int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
1811  int k = p - j - i;
1812  TriDofOrd[0][o] = o; // (0,1,2)
1813  TriDofOrd[1][o] = -1-(TriDof-((pp2-j)*(pp1-j))/2+k); // (1,0,2)
1814  TriDofOrd[2][o] = TriDof-((pp2-i)*(pp1-i))/2+k; // (2,0,1)
1815  TriDofOrd[3][o] = -1-(TriDof-((pp2-k)*(pp1-k))/2+i); // (2,1,0)
1816  TriDofOrd[4][o] = TriDof-((pp2-k)*(pp1-k))/2+j; // (1,2,0)
1817  TriDofOrd[5][o] = -1-(TriDof-((pp2-i)*(pp1-i))/2+j); // (0,2,1)
1818  if (!signs)
1819  {
1820  for (int k = 1; k < 6; k += 2)
1821  {
1822  TriDofOrd[k][o] = -1 - TriDofOrd[k][o];
1823  }
1824  }
1825  }
1826 
1827  int QuadDof = RT_dof[Geometry::SQUARE];
1828  QuadDofOrd[0] = new int[8*QuadDof];
1829  for (int i = 1; i < 8; i++)
1830  {
1831  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
1832  }
1833  // see Mesh::GetQuadOrientation in mesh/mesh.cpp
1834  for (int j = 0; j <= p; j++)
1835  for (int i = 0; i <= p; i++)
1836  {
1837  int o = i + j*pp1;
1838  QuadDofOrd[0][o] = i + j*pp1; // (0,1,2,3)
1839  QuadDofOrd[1][o] = -1 - (j + i*pp1); // (0,3,2,1)
1840  QuadDofOrd[2][o] = j + (p - i)*pp1; // (1,2,3,0)
1841  QuadDofOrd[3][o] = -1 - ((p - i) + j*pp1); // (1,0,3,2)
1842  QuadDofOrd[4][o] = (p - i) + (p - j)*pp1; // (2,3,0,1)
1843  QuadDofOrd[5][o] = -1 - ((p - j) + (p - i)*pp1); // (2,1,0,3)
1844  QuadDofOrd[6][o] = (p - j) + i*pp1; // (3,0,1,2)
1845  QuadDofOrd[7][o] = -1 - (i + (p - j)*pp1); // (3,2,1,0)
1846  if (!signs)
1847  {
1848  for (int k = 1; k < 8; k += 2)
1849  {
1850  QuadDofOrd[k][o] = -1 - QuadDofOrd[k][o];
1851  }
1852  }
1853  }
1854  }
1855 }
1856 
1857 int *RT_FECollection::DofOrderForOrientation(int GeomType, int Or) const
1858 {
1859  if (GeomType == Geometry::SEGMENT)
1860  {
1861  if (Or > 0)
1862  {
1863  return SegDofOrd[0];
1864  }
1865  return SegDofOrd[1];
1866  }
1867  else if (GeomType == Geometry::TRIANGLE)
1868  {
1869  return TriDofOrd[Or%6];
1870  }
1871  else if (GeomType == Geometry::SQUARE)
1872  {
1873  return QuadDofOrd[Or%8];
1874  }
1875  return NULL;
1876 }
1877 
1879 {
1880  return new RT_Trace_FECollection(atoi(rt_name + 7), atoi(rt_name + 3));
1881 }
1882 
1884 {
1885  delete [] SegDofOrd[0];
1886  delete [] TriDofOrd[0];
1887  delete [] QuadDofOrd[0];
1888  for (int g = 0; g < Geometry::NumGeom; g++)
1889  {
1890  delete RT_Elements[g];
1891  }
1892 }
1893 
1895  const int map_type)
1896  : RT_FECollection(p, dim, map_type, true)
1897 {
1898  if (map_type == FiniteElement::INTEGRAL)
1899  {
1900  snprintf(rt_name, 32, "RT_Trace_%dD_P%d", dim, p);
1901  }
1902  else
1903  {
1904  snprintf(rt_name, 32, "RT_ValTrace_%dD_P%d", dim, p);
1905  }
1906 
1907  MFEM_VERIFY(dim == 2 || dim == 3, "Wrong dimension, dim = " << dim);
1908 }
1909 
1911  const int map_type)
1912  : RT_FECollection(p, dim, map_type, false)
1913 {
1914  if (map_type == FiniteElement::VALUE)
1915  {
1916  snprintf(rt_name, 32, "DG_Iface_%dD_P%d", dim, p);
1917  }
1918  else
1919  {
1920  snprintf(rt_name, 32, "DG_IntIface_%dD_P%d", dim, p);
1921  }
1922 
1923  MFEM_VERIFY(dim == 2 || dim == 3, "Wrong dimension, dim = " << dim);
1924 }
1925 
1926 ND_FECollection::ND_FECollection(const int p, const int dim)
1927 {
1928  const int pm1 = p - 1, pm2 = p - 2;
1929 
1930  snprintf(nd_name, 32, "ND_%dD_P%d", dim, p);
1931 
1932  for (int g = 0; g < Geometry::NumGeom; g++)
1933  {
1934  ND_Elements[g] = NULL;
1935  ND_dof[g] = 0;
1936  }
1937  for (int i = 0; i < 2; i++)
1938  {
1939  SegDofOrd[i] = NULL;
1940  }
1941  for (int i = 0; i < 6; i++)
1942  {
1943  TriDofOrd[i] = NULL;
1944  }
1945  for (int i = 0; i < 8; i++)
1946  {
1947  QuadDofOrd[i] = NULL;
1948  }
1949 
1950  if (dim >= 1)
1951  {
1954 
1955  SegDofOrd[0] = new int[2*p];
1956  SegDofOrd[1] = SegDofOrd[0] + p;
1957  for (int i = 0; i < p; i++)
1958  {
1959  SegDofOrd[0][i] = i;
1960  SegDofOrd[1][i] = -1 - (pm1 - i);
1961  }
1962  }
1963 
1964  if (dim >= 2)
1965  {
1967  ND_dof[Geometry::SQUARE] = 2*p*pm1;
1968 
1970  ND_dof[Geometry::TRIANGLE] = p*pm1;
1971 
1972  int QuadDof = ND_dof[Geometry::SQUARE];
1973  QuadDofOrd[0] = new int[8*QuadDof];
1974  for (int i = 1; i < 8; i++)
1975  {
1976  QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
1977  }
1978  // see Mesh::GetQuadOrientation in mesh/mesh.cpp
1979  for (int j = 0; j < pm1; j++)
1980  {
1981  for (int i = 0; i < p; i++)
1982  {
1983  int d1 = i + j*p; // x-component
1984  int d2 = p*pm1 + j + i*pm1; // y-component
1985  // (0,1,2,3)
1986  QuadDofOrd[0][d1] = d1;
1987  QuadDofOrd[0][d2] = d2;
1988  // (0,3,2,1)
1989  QuadDofOrd[1][d1] = d2;
1990  QuadDofOrd[1][d2] = d1;
1991  // (1,2,3,0)
1992  // QuadDofOrd[2][d1] = p*pm1 + (pm2 - j) + i*pm1;
1993  // QuadDofOrd[2][d2] = -1 - ((pm1 - i) + j*p);
1994  QuadDofOrd[2][d1] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
1995  QuadDofOrd[2][d2] = i + (pm2 - j)*p;
1996  // (1,0,3,2)
1997  QuadDofOrd[3][d1] = -1 - ((pm1 - i) + j*p);
1998  QuadDofOrd[3][d2] = p*pm1 + (pm2 - j) + i*pm1;
1999  // (2,3,0,1)
2000  QuadDofOrd[4][d1] = -1 - ((pm1 - i) + (pm2 - j)*p);
2001  QuadDofOrd[4][d2] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2002  // (2,1,0,3)
2003  QuadDofOrd[5][d1] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2004  QuadDofOrd[5][d2] = -1 - ((pm1 - i) + (pm2 - j)*p);
2005  // (3,0,1,2)
2006  // QuadDofOrd[6][d1] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2007  // QuadDofOrd[6][d2] = i + (pm2 - j)*p;
2008  QuadDofOrd[6][d1] = p*pm1 + (pm2 - j) + i*pm1;
2009  QuadDofOrd[6][d2] = -1 - ((pm1 - i) + j*p);
2010  // (3,2,1,0)
2011  QuadDofOrd[7][d1] = i + (pm2 - j)*p;
2012  QuadDofOrd[7][d2] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2013  }
2014  }
2015 
2016  int TriDof = ND_dof[Geometry::TRIANGLE];
2017  TriDofOrd[0] = new int[6*TriDof];
2018  for (int i = 1; i < 6; i++)
2019  {
2020  TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
2021  }
2022  // see Mesh::GetTriOrientation in mesh/mesh.cpp,
2023  // the constructor of H1_FECollection
2024  for (int j = 0; j <= pm2; j++)
2025  {
2026  for (int i = 0; i + j <= pm2; i++)
2027  {
2028  int k1 = p*pm1 - (p - j)*(pm1 - j) + 2*i;
2029  int k2 = p*pm1 - (p - i)*(pm1 - i) + 2*j;
2030  // (0,1,2)
2031  TriDofOrd[0][k1 ] = k1;
2032  TriDofOrd[0][k1+1] = k1 + 1;
2033  // (0,2,1)
2034  TriDofOrd[5][k1 ] = k2 + 1;
2035  TriDofOrd[5][k1+1] = k2;
2036 
2037  // The other orientations can not be supported with the current
2038  // interface. The method Mesh::ReorientTetMesh will ensure that
2039  // only orientations 0 and 5 are generated.
2040  }
2041  }
2042  }
2043 
2044  if (dim >= 3)
2045  {
2047  ND_dof[Geometry::CUBE] = 3*p*pm1*pm1;
2048 
2050  ND_dof[Geometry::TETRAHEDRON] = p*pm1*pm2/2;
2051  }
2052 }
2053 
2054 int *ND_FECollection::DofOrderForOrientation(int GeomType, int Or) const
2055 {
2056  if (GeomType == Geometry::SEGMENT)
2057  {
2058  if (Or > 0)
2059  {
2060  return SegDofOrd[0];
2061  }
2062  return SegDofOrd[1];
2063  }
2064  else if (GeomType == Geometry::TRIANGLE)
2065  {
2066  if (Or != 0 && Or != 5)
2067  {
2068  MFEM_ABORT("ND_FECollection::DofOrderForOrientation: "
2069  "triangle face orientation " << Or << " is not supported! "
2070  "Use Mesh::ReorientTetMesh to fix it.");
2071  }
2072  return TriDofOrd[Or%6];
2073  }
2074  else if (GeomType == Geometry::SQUARE)
2075  {
2076  return QuadDofOrd[Or%8];
2077  }
2078  return NULL;
2079 }
2080 
2082 {
2084  atoi(nd_name + 3));
2085 }
2086 
2088 {
2089  delete [] SegDofOrd[0];
2090  delete [] TriDofOrd[0];
2091  delete [] QuadDofOrd[0];
2092  for (int g = 0; g < Geometry::NumGeom; g++)
2093  {
2094  delete ND_Elements[g];
2095  }
2096 }
2097 
2098 
2100  : ND_FECollection(p, dim-1)
2101 {
2102  snprintf(nd_name, 32, "ND_Trace_%dD_P%d", dim, p);
2103 }
2104 
2105 
2107 {
2108  snprintf(d_name, 32, "Local_%s", fe_name);
2109 
2110  Local_Element = NULL;
2111 
2112  if (!strcmp(fe_name, "BiCubic2DFiniteElement") ||
2113  !strcmp(fe_name, "Quad_Q3"))
2114  {
2115  GeomType = Geometry::SQUARE;
2116  Local_Element = new BiCubic2DFiniteElement;
2117  }
2118  else if (!strcmp(fe_name, "Nedelec1HexFiniteElement") ||
2119  !strcmp(fe_name, "Hex_ND1"))
2120  {
2121  GeomType = Geometry::CUBE;
2122  Local_Element = new Nedelec1HexFiniteElement;
2123  }
2124  else if (!strncmp(fe_name, "H1_", 3))
2125  {
2126  GeomType = Geometry::SQUARE;
2127  Local_Element = new H1_QuadrilateralElement(atoi(fe_name + 7));
2128  }
2129  else if (!strncmp(fe_name, "L2_", 3))
2130  {
2131  GeomType = Geometry::SQUARE;
2132  Local_Element = new L2_QuadrilateralElement(atoi(fe_name + 7));
2133  }
2134  else
2135  {
2136  cerr << "Local_FECollection::Local_FECollection : fe_name = "
2137  << fe_name << endl;
2138  mfem_error();
2139  }
2140 }
2141 
2142 
2143 void NURBSFECollection::Allocate(int Order)
2144 {
2145  SegmentFE = new NURBS1DFiniteElement(Order);
2146  QuadrilateralFE = new NURBS2DFiniteElement(Order);
2147  ParallelepipedFE = new NURBS3DFiniteElement(Order);
2148 
2149  snprintf(name, 16, "NURBS%i", Order);
2150 }
2151 
2152 void NURBSFECollection::Deallocate()
2153 {
2154  delete ParallelepipedFE;
2155  delete QuadrilateralFE;
2156  delete SegmentFE;
2157 }
2158 
2159 const FiniteElement *
2161 {
2162  switch (GeomType)
2163  {
2164  case Geometry::SEGMENT: return SegmentFE;
2165  case Geometry::SQUARE: return QuadrilateralFE;
2166  case Geometry::CUBE: return ParallelepipedFE;
2167  default:
2168  mfem_error ("NURBSFECollection: unknown geometry type.");
2169  }
2170  return SegmentFE; // Make some compilers happy
2171 }
2172 
2173 int NURBSFECollection::DofForGeometry(int GeomType) const
2174 {
2175  mfem_error("NURBSFECollection::DofForGeometry");
2176  return 0; // Make some compilers happy
2177 }
2178 
2179 int *NURBSFECollection::DofOrderForOrientation(int GeomType, int Or) const
2180 {
2181  mfem_error("NURBSFECollection::DofOrderForOrientation");
2182  return NULL;
2183 }
2184 
2186 {
2187  MFEM_ABORT("NURBS finite elements can not be statically condensed!");
2188  return NULL;
2189 }
2190 
2191 }
Abstract class for Finite Elements.
Definition: fe.hpp:44
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:263
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1174
int Size() const
Logical size of the array.
Definition: array.hpp:109
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:616
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:184
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:235
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1279
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:38
H1_Trace_FECollection(const int p, const int dim, const int type=GaussLobatto)
Definition: fe_coll.cpp:1569
void SubDofOrder(int Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:284
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:546
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2185
static const int NumGeom
Definition: geom.hpp:34
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2179
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1342
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1102
const Geometry::Type geom
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1020
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1707
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:635
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1263
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:865
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:923
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:901
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:831
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:306
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1370
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1202
ND_FECollection(const int p, const int dim)
Definition: fe_coll.cpp:1926
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:741
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1209
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1308
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2081
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1054
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1040
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:823
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:1749
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1878
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1241
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1138
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:768
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:991
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:714
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:470
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:493
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:755
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:893
virtual ~RT_FECollection()
Definition: fe_coll.cpp:1883
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:809
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1857
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1324
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:373
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:857
PointFiniteElement PointFE
Definition: point.cpp:30
int dim
Definition: ex3.cpp:47
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1355
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:183
static const int Dimension[NumGeom]
Definition: geom.hpp:38
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:209
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:559
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:843
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:486
L2_FECollection(const int p, const int dim, const int type=GaussLegendre)
Definition: fe_coll.cpp:1584
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:910
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1114
ND_Trace_FECollection(const int p, const int dim)
Definition: fe_coll.cpp:2099
TriLinear3DFiniteElement HexahedronFE
Definition: hexahedron.cpp:52
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:809
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:534
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:600
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2087
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1130
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:672
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL)
Definition: fe_coll.cpp:1894
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:943
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:797
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:353
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1062
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1225
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1557
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:454
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1294
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1543
static const char * Name[NumGeom]
Definition: geom.hpp:36
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:417
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:601
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1910
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:397
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:179
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:662
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:648
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:219
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:727
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1150
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:957
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
void mfem_error(const char *msg)
Definition: error.cpp:23
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:331
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:965
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:525
Linear3DFiniteElement TetrahedronFE
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:878
Linear2DFiniteElement TriangleFE
Definition: triangle.cpp:198
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1186
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1028
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1005
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:544
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:2160
static FiniteElementCollection * New(const char *name)
Definition: fe_coll.cpp:44
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1166
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1251
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:509
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:782
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2054
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1522
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:699
RT_FECollection(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.hpp:192
H1_FECollection(const int p, const int dim=3, const int type=GaussLobatto)
Definition: fe_coll.cpp:1392
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:329
int HasFaceDofs(int GeomType) const
Definition: fe_coll.cpp:25
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:231
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, int &fg, int &fo, const int face_info)
Definition: fe_coll.cpp:239
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1092
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1690
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:77
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:2173
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2106
BiLinear2DFiniteElement QuadrilateralFE
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:976
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1076
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:565
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:585
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:236
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:685
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:569
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:130
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:931