17 #define snprintf _snprintf_s
35 mfem_error (
"FiniteElementCollection::HasFaceDofs:"
36 " unknown geometry type.");
43 MFEM_ABORT(
"this method is not implemented in this derived class!");
51 if (!strcmp(name,
"Linear"))
55 else if (!strcmp(name,
"Quadratic"))
59 else if (!strcmp(name,
"QuadraticPos"))
63 else if (!strcmp(name,
"Cubic"))
67 else if (!strcmp(name,
"Const3D"))
71 else if (!strcmp(name,
"Const2D"))
75 else if (!strcmp(name,
"LinearDiscont2D"))
79 else if (!strcmp(name,
"GaussLinearDiscont2D"))
83 else if (!strcmp(name,
"P1OnQuad"))
87 else if (!strcmp(name,
"QuadraticDiscont2D"))
91 else if (!strcmp(name,
"QuadraticPosDiscont2D"))
95 else if (!strcmp(name,
"GaussQuadraticDiscont2D"))
99 else if (!strcmp(name,
"CubicDiscont2D"))
103 else if (!strcmp(name,
"LinearDiscont3D"))
107 else if (!strcmp(name,
"QuadraticDiscont3D"))
111 else if (!strcmp(name,
"LinearNonConf3D"))
115 else if (!strcmp(name,
"CrouzeixRaviart"))
119 else if (!strcmp(name,
"ND1_3D"))
123 else if (!strcmp(name,
"RT0_2D"))
127 else if (!strcmp(name,
"RT1_2D"))
131 else if (!strcmp(name,
"RT2_2D"))
135 else if (!strcmp(name,
"RT0_3D"))
139 else if (!strcmp(name,
"RT1_3D"))
143 else if (!strncmp(name,
"H1_Trace_", 9))
147 else if (!strncmp(name,
"H1_Trace@", 9))
152 else if (!strncmp(name,
"H1_", 3))
156 else if (!strncmp(name,
"H1Pos_Trace_", 12))
161 else if (!strncmp(name,
"H1Pos_", 6))
165 else if (!strncmp(name,
"H1Ser_", 6))
169 else if (!strncmp(name,
"H1@", 3))
174 else if (!strncmp(name,
"L2_T", 4))
177 else if (!strncmp(name,
"L2_", 3))
181 else if (!strncmp(name,
"L2Int_T", 7))
186 else if (!strncmp(name,
"L2Int_", 6))
192 else if (!strncmp(name,
"RT_Trace_", 9))
196 else if (!strncmp(name,
"RT_ValTrace_", 12))
201 else if (!strncmp(name,
"RT_Trace@", 9))
207 else if (!strncmp(name,
"RT_ValTrace@", 12))
213 else if (!strncmp(name,
"DG_Iface_", 9))
217 else if (!strncmp(name,
"DG_Iface@", 9))
223 else if (!strncmp(name,
"DG_IntIface_", 12))
228 else if (!strncmp(name,
"DG_IntIface@", 12))
234 else if (!strncmp(name,
"RT_", 3))
238 else if (!strncmp(name,
"RT@", 3))
244 else if (!strncmp(name,
"ND_Trace_", 9))
248 else if (!strncmp(name,
"ND_Trace@", 9))
254 else if (!strncmp(name,
"ND_", 3))
258 else if (!strncmp(name,
"ND@", 3))
264 else if (!strncmp(name,
"Local_", 6))
268 else if (!strncmp(name,
"NURBS", 5))
283 MFEM_ABORT(
"unknown FiniteElementCollection: " << name);
285 MFEM_VERIFY(!strcmp(fec->
Name(), name),
"input name: \"" << name
286 <<
"\" does not match the created collection name: \""
287 << fec->
Name() <<
'"');
292 template <Geometry::Type geom>
297 nv = g_consts::NumVert;
298 ne = g_consts::NumEdges;
301 template <Geometry::Type geom,
typename v_t>
303 GetEdge(
int &nv, v_t &v,
int &ne,
int &e,
int &eo,
const int edge_info)
308 nv = e_consts::NumVert;
312 MFEM_ASSERT(0 <= e && e < g_consts::NumEdges,
"");
313 MFEM_ASSERT(0 <= eo && eo < e_consts::NumOrient,
"");
314 v[0] = g_consts::Edges[e][0];
315 v[1] = g_consts::Edges[e][1];
316 v[0] = e_consts::Orient[eo][v[0]];
317 v[1] = e_consts::Orient[eo][v[1]];
321 typename v_t,
typename e_t,
typename eo_t>
323 GetFace(
int &nv, v_t &v,
int &ne, e_t &e, eo_t &eo,
324 int &nf,
int &f,
Geometry::Type &fg,
int &fo,
const int face_info)
329 nv = f_consts::NumVert;
334 MFEM_ASSERT(0 <= f && f < g_consts::NumFaces,
"");
335 MFEM_ASSERT(0 <= fo && fo < f_consts::NumOrient,
"");
336 for (
int i = 0; i < f_consts::NumVert; i++)
338 v[i] = f_consts::Orient[fo][i];
339 v[i] = g_consts::FaceVert[f][v[i]];
341 ne = f_consts::NumEdges;
342 for (
int i = 0; i < f_consts::NumEdges; i++)
344 int v0 = v[f_consts::Edges[i][0]];
345 int v1 = v[f_consts::Edges[i][1]];
347 if (v0 > v1) { swap(v0, v1); eor = 1; }
348 for (
int j = g_consts::VertToVert::I[v0];
true; j++)
350 MFEM_ASSERT(j < g_consts::VertToVert::I[v0+1],
351 "internal error, edge not found");
352 if (v1 == g_consts::VertToVert::J[j][0])
354 int en = g_consts::VertToVert::J[j][1];
374 "invalid Geom = " << Geom);
376 "invalid SDim = " << SDim <<
382 const int off = nvd*(Info/64);
384 for (
int i = 0; i < nvd; i++)
391 int v[4], e[4], eo[4], f[1], fo[1];
392 int av = 0, nv = 0, ae = 0, ne = 0, nf = 0;
399 GetNVE<Geometry::SEGMENT>(av, ae);
400 GetEdge<Geometry::SEGMENT>(nv, v, ne, e[0], eo[0], Info);
406 GetNVE<Geometry::TRIANGLE>(av, ae);
410 GetEdge<Geometry::TRIANGLE>(nv, v, ne, e[0], eo[0], Info);
413 GetFace<Geometry::TRIANGLE,Geometry::TRIANGLE>(
414 nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
424 GetNVE<Geometry::SQUARE>(av, ae);
428 GetEdge<Geometry::SQUARE>(nv, v, ne, e[0], eo[0], Info);
431 GetFace<Geometry::SQUARE,Geometry::SQUARE>(
432 nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
442 GetNVE<Geometry::TETRAHEDRON>(av, ae);
446 GetEdge<Geometry::TETRAHEDRON>(nv, v, ne, e[0], eo[0], Info);
449 GetFace<Geometry::TETRAHEDRON,Geometry::TRIANGLE>(
450 nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
460 GetNVE<Geometry::CUBE>(av, ae);
464 GetEdge<Geometry::CUBE>(nv, v, ne, e[0], eo[0], Info);
467 GetFace<Geometry::CUBE,Geometry::SQUARE>(
468 nv, v, ne, e, eo, nf, f[0], fg[0], fo[0], Info);
477 MFEM_ABORT(
"invalid Geom = " << Geom);
484 for (
int i = 0; i < nv; i++)
486 for (
int j = 0; j < nvd; j++)
488 dofs[i*nvd+j] = v[i]*nvd+j;
491 int l_off = nv*nvd, g_off = av*nvd;
496 for (
int i = 0; i < ne; i++)
500 for (
int j = 0; j < ned; j++)
502 dofs[l_off+i*ned+j] =
504 g_off+e[i]*ned+ed[j] :
505 -1-(g_off+e[i]*ned+(-1-ed[j]));
515 const int nfd = DofForGeometry(fg[0]);
517 for (
int i = 0; i < nf; i++)
519 const int *fd = DofOrderForOrientation(fg[i], fo[i]);
520 for (
int j = 0; j < nfd; j++)
522 dofs[l_off+i*nfd+j] =
524 g_off+f[i]*nfd+fd[j] :
525 -1-(g_off+f[i]*nfd+(-1-fd[j]));
536 ", SDim = " << SDim <<
" is not supported");
552 mfem_error (
"LinearFECollection: unknown geometry type.");
569 mfem_error (
"LinearFECollection: unknown geometry type.");
594 mfem_error (
"QuadraticFECollection: unknown geometry type.");
611 mfem_error (
"QuadraticFECollection: unknown geometry type.");
619 static int indexes[] = { 0 };
634 mfem_error (
"QuadraticPosFECollection: unknown geometry type.");
647 mfem_error (
"QuadraticPosFECollection: unknown geometry type.");
655 static int indexes[] = { 0 };
674 mfem_error (
"CubicFECollection: unknown geometry type.");
691 mfem_error (
"CubicFECollection: unknown geometry type.");
701 static int ind_pos[] = { 0, 1 };
702 static int ind_neg[] = { 1, 0 };
712 static int indexes[] = { 0 };
718 static int sq_ind[8][4] = {{0, 1, 2, 3}, {0, 2, 1, 3},
719 {2, 0, 3, 1}, {1, 0, 3, 2},
720 {3, 2, 1, 0}, {3, 1, 2, 0},
721 {1, 3, 0, 2}, {2, 3, 0, 1}
740 mfem_error (
"CrouzeixRaviartFECollection: unknown geometry type.");
754 mfem_error (
"CrouzeixRaviartFECollection: unknown geometry type.");
762 static int indexes[] = { 0 };
777 mfem_error (
"RT0_2DFECollection: unknown geometry type.");
791 mfem_error (
"RT0_2DFECollection: unknown geometry type.");
799 static int ind_pos[] = { 0 };
800 static int ind_neg[] = { -1 };
819 mfem_error (
"RT1_2DFECollection: unknown geometry type.");
833 mfem_error (
"RT1_2DFECollection: unknown geometry type.");
841 static int ind_pos[] = { 0, 1 };
842 static int ind_neg[] = { -2, -1 };
860 mfem_error (
"RT2_2DFECollection: unknown geometry type.");
874 mfem_error (
"RT2_2DFECollection: unknown geometry type.");
882 static int ind_pos[] = { 0, 1, 2 };
883 static int ind_neg[] = { -3, -2, -1 };
901 mfem_error (
"Const2DFECollection: unknown geometry type.");
915 mfem_error (
"Const2DFECollection: unknown geometry type.");
936 mfem_error (
"LinearDiscont2DFECollection: unknown geometry type.");
950 mfem_error (
"LinearDiscont2DFECollection: unknown geometry type.");
971 mfem_error (
"GaussLinearDiscont2DFECollection:"
972 " unknown geometry type.");
987 mfem_error (
"GaussLinearDiscont2DFECollection:"
988 " unknown geometry type.");
1005 mfem_error (
"P1OnQuadFECollection: unknown geometry type.");
1018 mfem_error (
"P1OnQuadFECollection: unknown geometry type.");
1039 mfem_error (
"QuadraticDiscont2DFECollection: unknown geometry type.");
1054 mfem_error (
"QuadraticDiscont2DFECollection: unknown geometry type.");
1074 mfem_error (
"QuadraticPosDiscont2DFECollection: unknown geometry type.");
1088 mfem_error (
"QuadraticPosDiscont2DFECollection: unknown geometry type.");
1104 mfem_error (
"GaussQuadraticDiscont2DFECollection:"
1105 " unknown geometry type.");
1120 mfem_error (
"GaussQuadraticDiscont2DFECollection:"
1121 " unknown geometry type.");
1142 mfem_error (
"CubicDiscont2DFECollection: unknown geometry type.");
1156 mfem_error (
"CubicDiscont2DFECollection: unknown geometry type.");
1179 mfem_error (
"LinearNonConf3DFECollection: unknown geometry type.");
1195 mfem_error (
"LinearNonConf3DFECollection: unknown geometry type.");
1203 static int indexes[] = { 0 };
1218 mfem_error (
"Const3DFECollection: unknown geometry type.");
1235 mfem_error (
"Const3DFECollection: unknown geometry type.");
1256 mfem_error (
"LinearDiscont3DFECollection: unknown geometry type.");
1272 mfem_error (
"LinearDiscont3DFECollection: unknown geometry type.");
1293 mfem_error (
"QuadraticDiscont3DFECollection: unknown geometry type.");
1310 mfem_error (
"QuadraticDiscont3DFECollection: unknown geometry type.");
1334 mfem_error (
"RefinedLinearFECollection: unknown geometry type.");
1350 mfem_error (
"RefinedLinearFECollection: unknown geometry type.");
1358 static int indexes[] = { 0 };
1372 mfem_error (
"ND1_3DFECollection: unknown geometry type.");
1388 mfem_error (
"ND1_3DFECollection: unknown geometry type.");
1396 static int ind_pos[] = { 0 };
1397 static int ind_neg[] = { -1 };
1417 mfem_error (
"RT0_3DFECollection: unknown geometry type.");
1433 mfem_error (
"RT0_3DFECollection: unknown geometry type.");
1441 static int ind_pos[] = { 0 };
1442 static int ind_neg[] = { -1 };
1464 mfem_error (
"RT1_3DFECollection: unknown geometry type.");
1479 mfem_error (
"RT1_3DFECollection: unknown geometry type.");
1489 static int sq_ind[8][4] =
1491 {0, 1, 2, 3}, {-1, -3, -2, -4},
1492 {2, 0, 3, 1}, {-2, -1, -4, -3},
1493 {3, 2, 1, 0}, {-4, -2, -3, -1},
1494 {1, 3, 0, 2}, {-3, -4, -1, -2}
1508 MFEM_VERIFY(p >= 1,
"H1_FECollection requires order >= 1.");
1509 MFEM_VERIFY(dim >= 0 && dim <= 3,
"H1_FECollection requires 0 <= dim <= 3.");
1511 const int pm1 = p - 1, pm2 = pm1 - 1, pm3 = pm2 - 1;
1519 snprintf(h1_name, 32,
"H1_%dD_P%d", dim, p);
1524 snprintf(h1_name, 32,
"H1Pos_%dD_P%d", dim, p);
1529 snprintf(h1_name, 32,
"H1Ser_%dD_P%d", dim, p);
1538 snprintf(h1_name, 32,
"H1@%c_%dD_P%d",
1546 H1_Elements[g] = NULL;
1548 for (
int i = 0; i < 2; i++)
1550 SegDofOrd[i] = NULL;
1552 for (
int i = 0; i < 6; i++)
1554 TriDofOrd[i] = NULL;
1556 for (
int i = 0; i < 8; i++)
1558 QuadDofOrd[i] = NULL;
1576 SegDofOrd[0] =
new int[2*pm1];
1577 SegDofOrd[1] = SegDofOrd[0] + pm1;
1578 for (
int i = 0; i < pm1; i++)
1580 SegDofOrd[0][i] = i;
1581 SegDofOrd[1][i] = pm2 - i;
1614 TriDofOrd[0] =
new int[6*TriDof];
1615 for (
int i = 1; i < 6; i++)
1617 TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1620 for (
int j = 0; j < pm2; j++)
1622 for (
int i = 0; i + j < pm2; i++)
1624 int o = TriDof - ((pm1 - j)*(pm2 - j))/2 + i;
1625 int k = pm3 - j - i;
1626 TriDofOrd[0][o] = o;
1627 TriDofOrd[1][o] = TriDof - ((pm1-j)*(pm2-j))/2 + k;
1628 TriDofOrd[2][o] = TriDof - ((pm1-i)*(pm2-i))/2 + k;
1629 TriDofOrd[3][o] = TriDof - ((pm1-k)*(pm2-k))/2 + i;
1630 TriDofOrd[4][o] = TriDof - ((pm1-k)*(pm2-k))/2 + j;
1631 TriDofOrd[5][o] = TriDof - ((pm1-i)*(pm2-i))/2 + j;
1635 QuadDofOrd[0] =
new int[8*QuadDof];
1636 for (
int i = 1; i < 8; i++)
1638 QuadDofOrd[i] = QuadDofOrd[i-1] + QuadDof;
1658 const int pm4 = pm3 -1;
1660 for (
int j = 0; j < pm3; j++)
1662 for (
int i = 0; i < pm3; i++)
1665 QuadDofOrd[0][o] = i + j*pm3;
1666 QuadDofOrd[1][o] = j + i*pm3;
1667 QuadDofOrd[2][o] = j + (pm4 - i)*pm3;
1668 QuadDofOrd[3][o] = (pm4 - i) + j*pm3;
1669 QuadDofOrd[4][o] = (pm4 - i) + (pm4 - j)*pm3;
1670 QuadDofOrd[5][o] = (pm4 - j) + (pm4 - i)*pm3;
1671 QuadDofOrd[6][o] = (pm4 - j) + i*pm3;
1672 QuadDofOrd[7][o] = i + (pm4 - j)*pm3;
1680 for (
int j = 0; j < pm1; j++)
1682 for (
int i = 0; i < pm1; i++)
1685 QuadDofOrd[0][o] = i + j*pm1;
1686 QuadDofOrd[1][o] = j + i*pm1;
1687 QuadDofOrd[2][o] = j + (pm2 - i)*pm1;
1688 QuadDofOrd[3][o] = (pm2 - i) + j*pm1;
1689 QuadDofOrd[4][o] = (pm2 - i) + (pm2 - j)*pm1;
1690 QuadDofOrd[5][o] = (pm2 - j) + (pm2 - i)*pm1;
1691 QuadDofOrd[6][o] = (pm2 - j) + i*pm1;
1692 QuadDofOrd[7][o] = i + (pm2 - j)*pm1;
1724 return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
1728 return TriDofOrd[Or%6];
1732 return QuadDofOrd[Or%8];
1741 if (!strncmp(h1_name,
"H1_", 3))
1743 dim = atoi(h1_name + 3);
1745 else if (!strncmp(h1_name,
"H1Pos_", 6))
1747 dim = atoi(h1_name + 6);
1749 else if (!strncmp(h1_name,
"H1@", 3))
1751 dim = atoi(h1_name + 5);
1758 const int *dof_map = NULL;
1766 ->GetDofMap().GetData();
1769 MFEM_ABORT(
"Geometry type " <<
Geometry::Name[GeomType] <<
" is not "
1779 delete [] SegDofOrd[0];
1780 delete [] TriDofOrd[0];
1781 delete [] QuadDofOrd[0];
1784 delete H1_Elements[g];
1795 snprintf(
h1_name, 32,
"H1_Trace_%dD_P%d", dim, p);
1799 snprintf(
h1_name, 32,
"H1Pos_Trace_%dD_P%d", dim, p);
1803 snprintf(
h1_name, 32,
"H1_Trace@%c_%dD_P%d",
1812 MFEM_VERIFY(p >= 0,
"L2_FECollection requires order >= 0.");
1815 const char *prefix = NULL;
1821 MFEM_ABORT(
"invalid map_type: " << map_type);
1826 snprintf(d_name, 32,
"%s_%dD_P%d", prefix, dim, p);
1829 snprintf(d_name, 32,
"%s_T%d_%dD_P%d", prefix, btype, dim, p);
1834 L2_Elements[g] = NULL;
1835 Tr_Elements[g] = NULL;
1837 for (
int i = 0; i < 2; i++)
1839 SegDofOrd[i] = NULL;
1841 for (
int i = 0; i < 6; i++)
1843 TriDofOrd[i] = NULL;
1866 const int pp1 = p + 1;
1867 SegDofOrd[0] =
new int[2*pp1];
1868 SegDofOrd[1] = SegDofOrd[0] + pp1;
1869 for (
int i = 0; i <= p; i++)
1871 SegDofOrd[0][i] = i;
1872 SegDofOrd[1][i] = p - i;
1900 TriDofOrd[0] =
new int[6*TriDof];
1901 for (
int i = 1; i < 6; i++)
1903 TriDofOrd[i] = TriDofOrd[i-1] + TriDof;
1905 const int pp1 = p + 1, pp2 = pp1 + 1;
1906 for (
int j = 0; j <= p; j++)
1908 for (
int i = 0; i + j <= p; i++)
1910 int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
1912 TriDofOrd[0][o] = o;
1913 TriDofOrd[1][o] = TriDof - ((pp2-j)*(pp1-j))/2 + k;
1914 TriDofOrd[2][o] = TriDof - ((pp2-i)*(pp1-i))/2 + k;
1915 TriDofOrd[3][o] = TriDof - ((pp2-k)*(pp1-k))/2 + i;
1916 TriDofOrd[4][o] = TriDof - ((pp2-k)*(pp1-k))/2 + j;
1917 TriDofOrd[5][o] = TriDof - ((pp2-i)*(pp1-i))/2 + j;
1921 OtherDofOrd =
new int[QuadDof];
1922 for (
int j = 0; j < QuadDof; j++)
1960 const int MaxDof = std::max(TetDof, std::max(PriDof, HexDof));
1961 OtherDofOrd =
new int[MaxDof];
1962 for (
int j = 0; j < MaxDof; j++)
1969 mfem::err <<
"L2_FECollection::L2_FECollection : dim = "
1981 return (Or > 0) ? SegDofOrd[0] : SegDofOrd[1];
1984 return TriDofOrd[Or%6];
1987 return (Or == 0) ? OtherDofOrd : NULL;
1993 delete [] OtherDofOrd;
1994 delete [] SegDofOrd[0];
1995 delete [] TriDofOrd[0];
1998 delete L2_Elements[i];
1999 delete Tr_Elements[i];
2005 const int cb_type,
const int ob_type)
2008 MFEM_VERIFY(p >= 0,
"RT_FECollection requires order >= 0.");
2016 MFEM_ABORT(
"unknown closed BasisType: " << cb_name);
2021 MFEM_ABORT(
"unknown open BasisType: " << ob_name);
2029 snprintf(
rt_name, 32,
"RT_%dD_P%d", dim, p);
2037 const int pp1 = p + 1;
2060 MFEM_ABORT(
"invalid dim = " << dim);
2067 const bool signs,
const int ob_type)
2074 MFEM_ABORT(
"Invalid open basis type: " << ob_name);
2085 "invalid open point type");
2087 const int pp1 = p + 1, pp2 = p + 2;
2095 for (
int i = 0; i < 2; i++)
2099 for (
int i = 0; i < 6; i++)
2103 for (
int i = 0; i < 8; i++)
2117 for (
int i = 0; i <= p; i++)
2120 SegDofOrd[1][i] = signs ? (-1 - (p - i)) : (p - i);
2137 for (
int i = 1; i < 6; i++)
2143 for (
int j = 0; j <= p; j++)
2145 for (
int i = 0; i + j <= p; i++)
2147 int o = TriDof - ((pp2 - j)*(pp1 - j))/2 + i;
2150 TriDofOrd[1][o] = -1-(TriDof-((pp2-j)*(pp1-j))/2+k);
2151 TriDofOrd[2][o] = TriDof-((pp2-i)*(pp1-i))/2+k;
2152 TriDofOrd[3][o] = -1-(TriDof-((pp2-k)*(pp1-k))/2+i);
2153 TriDofOrd[4][o] = TriDof-((pp2-k)*(pp1-k))/2+j;
2154 TriDofOrd[5][o] = -1-(TriDof-((pp2-i)*(pp1-i))/2+j);
2157 for (
int k = 1; k < 6; k += 2)
2159 TriDofOrd[k][o] = -1 - TriDofOrd[k][o];
2167 for (
int i = 1; i < 8; i++)
2172 for (
int j = 0; j <= p; j++)
2174 for (
int i = 0; i <= p; i++)
2181 QuadDofOrd[4][o] = (p - i) + (p - j)*pp1;
2182 QuadDofOrd[5][o] = -1 - ((p - j) + (p - i)*pp1);
2183 QuadDofOrd[6][o] = (p - j) + i*pp1;
2184 QuadDofOrd[7][o] = -1 - (i + (p - j)*pp1);
2187 for (
int k = 1; k < 8; k += 2)
2189 QuadDofOrd[k][o] = -1 - QuadDofOrd[k][o];
2218 if (!strncmp(
rt_name,
"RT_", 3))
2248 const char *prefix =
2250 char ob_str[3] = {
'\0',
'\0',
'\0' };
2257 snprintf(
rt_name, 32,
"%s%s_%dD_P%d", prefix, ob_str, dim, p);
2259 MFEM_VERIFY(dim == 2 || dim == 3,
"Wrong dimension, dim = " << dim);
2268 MFEM_VERIFY(dim == 2 || dim == 3,
"Wrong dimension, dim = " << dim);
2270 const char *prefix =
2274 snprintf(
rt_name, 32,
"%s_%dD_P%d", prefix, dim, p);
2278 snprintf(
rt_name, 32,
"%s@%c_%dD_P%d", prefix,
2284 const int cb_type,
const int ob_type)
2286 MFEM_VERIFY(p >= 1,
"ND_FECollection requires order >= 1.");
2287 MFEM_VERIFY(dim >= 1 && dim <= 3,
"ND_FECollection requires 1 <= dim <= 3.");
2289 const int pm1 = p - 1, pm2 = p - 2;
2294 snprintf(
nd_name, 32,
"ND_%dD_P%d", dim, p);
2307 for (
int i = 0; i < 2; i++)
2311 for (
int i = 0; i < 6; i++)
2315 for (
int i = 0; i < 8; i++)
2327 MFEM_ABORT(
"Invalid open basis point type: " << ob_name);
2332 MFEM_ABORT(
"Invalid closed basis point type: " << cb_name);
2342 for (
int i = 0; i < p; i++)
2361 for (
int i = 1; i < 8; i++)
2366 for (
int j = 0; j < pm1; j++)
2368 for (
int i = 0; i < p; i++)
2371 int d2 = p*pm1 + j + i*pm1;
2381 QuadDofOrd[2][d1] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2385 QuadDofOrd[3][d2] = p*pm1 + (pm2 - j) + i*pm1;
2387 QuadDofOrd[4][d1] = -1 - ((pm1 - i) + (pm2 - j)*p);
2388 QuadDofOrd[4][d2] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2390 QuadDofOrd[5][d1] = -1 - (p*pm1 + (pm2 - j) + (pm1 - i)*pm1);
2391 QuadDofOrd[5][d2] = -1 - ((pm1 - i) + (pm2 - j)*p);
2395 QuadDofOrd[6][d1] = p*pm1 + (pm2 - j) + i*pm1;
2396 QuadDofOrd[6][d2] = -1 - ((pm1 - i) + j*p);
2398 QuadDofOrd[7][d1] = i + (pm2 - j)*p;
2399 QuadDofOrd[7][d2] = -1 - (p*pm1 + j + (pm1 - i)*pm1);
2405 for (
int i = 1; i < 6; i++)
2411 for (
int j = 0; j <= pm2; j++)
2413 for (
int i = 0; i + j <= pm2; i++)
2415 int k1 = p*pm1 - (p - j)*(pm1 - j) + 2*i;
2416 int k2 = p*pm1 - (p - i)*(pm1 - i) + 2*j;
2451 if (Or != 0 && Or != 5)
2453 MFEM_ABORT(
"triangle face orientation " << Or <<
" is not supported! "
2454 "Use Mesh::ReorientTetMesh to fix it.");
2467 int p,
dim, cb_type, ob_type;
2505 snprintf(
nd_name, 32,
"ND_Trace_%dD_P%d", dim, p);
2509 snprintf(
nd_name, 32,
"ND_Trace@%c%c_%dD_P%d",
2518 snprintf(d_name, 32,
"Local_%s", fe_name);
2520 Local_Element = NULL;
2522 if (!strcmp(fe_name,
"BiCubic2DFiniteElement") ||
2523 !strcmp(fe_name,
"Quad_Q3"))
2528 else if (!strcmp(fe_name,
"Nedelec1HexFiniteElement") ||
2529 !strcmp(fe_name,
"Hex_ND1"))
2534 else if (!strncmp(fe_name,
"H1_", 3))
2539 else if (!strncmp(fe_name,
"H1Pos_", 6))
2544 else if (!strncmp(fe_name,
"L2_", 3))
2551 mfem::err <<
"Local_FECollection::Local_FECollection : fe_name = "
2573 snprintf(name, 16,
"NURBS%i", Order);
2577 snprintf(name, 16,
"NURBS");
2583 delete ParallelepipedFE;
2584 delete QuadrilateralFE;
2597 mfem_error (
"NURBSFECollection: unknown geometry type.");
2604 mfem_error(
"NURBSFECollection::DofForGeometry");
2611 mfem_error(
"NURBSFECollection::DofOrderForOrientation");
2617 MFEM_ABORT(
"NURBS finite elements can not be statically condensed!");
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Abstract class for Finite Elements.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
int Size() const
Logical size of the array.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
For scalar fields; preserves point values.
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
int RT_dof[Geometry::NumGeom]
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
FiniteElement * ND_Elements[Geometry::NumGeom]
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
virtual FiniteElementCollection * GetTraceCollection() const
virtual ~NURBSFECollection()
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
RT_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
FiniteElementCollection * GetTraceCollection() const
virtual ~L2_FECollection()
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Piecewise-(bi)linear continuous finite elements.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
FiniteElementCollection * GetTraceCollection() const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
FiniteElementCollection * GetTraceCollection() const
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE, const int ob_type=BasisType::GaussLegendre)
int HasFaceDofs(Geometry::Type GeomType) const
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual ~RT_FECollection()
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Piecewise-(bi)cubic continuous finite elements.
PointFiniteElement PointFE
virtual int DofForGeometry(Geometry::Type GeomType) const
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
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)
FiniteElement * RT_Elements[Geometry::NumGeom]
static const int Dimension[NumGeom]
virtual int DofForGeometry(Geometry::Type GeomType) const
static void GetNVE(int &nv, int &ne)
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
For scalar fields; preserves volume integrals.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
ND_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
TriLinear3DFiniteElement HexahedronFE
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Discontinuous collection defined locally by a given finite element.
Version of QuadraticDiscont2DFECollection with positive basis functions.
virtual ~ND_FECollection()
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
void SubDofOrder(Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual int DofForGeometry(Geometry::Type GeomType) const
Version of QuadraticFECollection with positive basis functions.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual ~H1_FECollection()
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
FiniteElementCollection * GetTraceCollection() const
static const char * Name[NumGeom]
Piecewise-linear nonconforming finite elements in 3D.
Crouzeix-Raviart nonconforming elements in 2D.
class H1_WedgeElement WedgeFE
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
void SetOrder(int Order) const
Set the order and the name, based on the given Order: either a positive number for fixed order...
virtual int DofForGeometry(Geometry::Type GeomType) const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual int DofForGeometry(Geometry::Type GeomType) const
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
class Linear3DFiniteElement TetrahedronFE
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Linear2DFiniteElement TriangleFE
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Serendipity basis (squares / cubes)
virtual int DofForGeometry(Geometry::Type GeomType) const
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
static char GetChar(int b_type)
Check and convert a BasisType constant to a char basis identifier.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const char * Name() const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
virtual int DofForGeometry(Geometry::Type GeomType) const
Piecewise-(bi)quadratic continuous finite elements.
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Arbitrary order H(curl)-conforming Nedelec finite elements.
virtual int DofForGeometry(Geometry::Type GeomType) const
Arbitrary order H1-conforming (continuous) finite elements.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual int DofForGeometry(Geometry::Type GeomType) const
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Local_FECollection(const char *fe_name)
virtual int DofForGeometry(Geometry::Type GeomType) const
BiLinear2DFiniteElement QuadrilateralFE
Linear (P1) finite elements on quadrilaterals.
Linear1DFiniteElement SegmentFE
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
int ND_dof[Geometry::NumGeom]
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Arbitrary order "L2-conforming" discontinuous finite elements.