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