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