16 #include "../coefficient.hpp" 23 const double RT_QuadrilateralElement::nk[8] =
24 { 0., -1., 1., 0., 0., 1., -1., 0. };
32 cp(
poly1d.ClosedPoints(
p + 1, cb_type))
39 const int dof2 =
dof/2;
41 #ifndef MFEM_THREAD_SAFE 52 for (
int i = 0; i <=
p; i++)
54 dof_map[1*dof2 + i + 0*(
p + 1)] = o++;
56 for (
int i = 0; i <=
p; i++)
58 dof_map[0*dof2 + (
p + 1) + i*(
p + 2)] = o++;
60 for (
int i = 0; i <=
p; i++)
62 dof_map[1*dof2 + (
p - i) + (
p + 1)*(
p + 1)] = o++;
64 for (
int i = 0; i <=
p; i++)
66 dof_map[0*dof2 + 0 + (
p - i)*(
p + 2)] = o++;
70 for (
int j = 0; j <=
p; j++)
71 for (
int i = 1; i <=
p; i++)
73 dof_map[0*dof2 + i + j*(
p + 2)] = o++;
75 for (
int j = 1; j <=
p; j++)
76 for (
int i = 0; i <=
p; i++)
78 dof_map[1*dof2 + i + j*(
p + 1)] = o++;
83 for (
int j = 0; j <=
p; j++)
84 for (
int i = 0; i <=
p/2; i++)
86 int idx = 0*dof2 + i + j*(
p + 2);
90 for (
int j =
p/2 + 1; j <=
p; j++)
92 int idx = 0*dof2 + (
p/2 + 1) + j*(
p + 2);
96 for (
int j = 0; j <=
p/2; j++)
97 for (
int i = 0; i <=
p; i++)
99 int idx = 1*dof2 + i + j*(
p + 1);
103 for (
int i = 0; i <=
p/2; i++)
105 int idx = 1*dof2 + i + (
p/2 + 1)*(
p + 1);
110 for (
int j = 0; j <=
p; j++)
111 for (
int i = 0; i <=
p + 1; i++)
125 for (
int j = 0; j <=
p + 1; j++)
126 for (
int i = 0; i <=
p; i++)
145 const int pp1 =
order;
147 #ifdef MFEM_THREAD_SAFE 148 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
153 #ifdef MFEM_THREAD_SAFE 154 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
171 for (
int j = 0; j < pp1; j++)
172 for (
int i = 0; i <= pp1; i++)
177 idx = -1 - idx,
s = -1;
183 shape(idx,0) =
s*shape_cx(i)*shape_oy(j);
186 for (
int j = 0; j <= pp1; j++)
187 for (
int i = 0; i < pp1; i++)
192 idx = -1 - idx,
s = -1;
199 shape(idx,1) =
s*shape_ox(i)*shape_cy(j);
206 const int pp1 =
order;
208 #ifdef MFEM_THREAD_SAFE 209 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
210 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
228 for (
int j = 0; j < pp1; j++)
229 for (
int i = 0; i <= pp1; i++)
234 idx = -1 - idx,
s = -1;
240 divshape(idx) =
s*dshape_cx(i)*shape_oy(j);
242 for (
int j = 0; j <= pp1; j++)
243 for (
int i = 0; i < pp1; i++)
248 idx = -1 - idx,
s = -1;
254 divshape(idx) =
s*shape_ox(i)*dshape_cy(j);
272 for (
int c = 0; c < 2; c++)
276 for (
int j = 0; j < jm; j++)
277 for (
int i = 0; i < im; i++)
280 if (idx < 0) { idx = -1 - idx; }
281 int ic = (c == 0) ? j : i;
282 const double h = cp[ic+1] - cp[ic];
284 for (
int k = 0; k < nqpt; k++)
287 if (c == 0) { ip2d.
Set2(cp[i], cp[j] + (h*ip1d.
x)); }
288 else { ip2d.
Set2(cp[i] + (h*ip1d.
x), cp[j]); }
290 vc.
Eval(xk, Trans, ip2d);
293 nk + dof2nk[idx]*
dim);
305 const int pp1 =
p + 1;
306 const int n_face_dofs =
p;
308 std::vector<int> offsets;
309 std::vector<int> strides = {(face_id == 0 || face_id == 2) ? 1 : pp1};
312 case 0: offsets = {
p*pp1};
break;
313 case 1: offsets = {pp1 - 1};
break;
314 case 2: offsets = {
p*pp1 +
p*(pp1 - 1)};
break;
315 case 3: offsets = {0};
break;
318 std::vector<int> n_dofs(
dim - 1,
p);
319 internal::FillFaceMap(n_face_dofs, offsets, strides, n_dofs, face_map);
323 const double RT_HexahedronElement::nk[18] =
324 { 0.,0.,-1., 0.,-1.,0., 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,0.,1. };
332 cp(
poly1d.ClosedPoints(
p + 1, cb_type))
339 const int dof3 =
dof/3;
341 #ifndef MFEM_THREAD_SAFE 355 for (
int j = 0; j <=
p; j++)
356 for (
int i = 0; i <=
p; i++)
358 dof_map[2*dof3 + i + ((
p - j) + 0*(
p + 1))*(
p + 1)] = o++;
360 for (
int j = 0; j <=
p; j++)
361 for (
int i = 0; i <=
p; i++)
363 dof_map[1*dof3 + i + (0 + j*(
p + 2))*(
p + 1)] = o++;
365 for (
int j = 0; j <=
p; j++)
366 for (
int i = 0; i <=
p; i++)
368 dof_map[0*dof3 + (
p + 1) + (i + j*(
p + 1))*(
p + 2)] = o++;
370 for (
int j = 0; j <=
p; j++)
371 for (
int i = 0; i <=
p; i++)
373 dof_map[1*dof3 + (
p - i) + ((
p + 1) + j*(
p + 2))*(
p + 1)] = o++;
375 for (
int j = 0; j <=
p; j++)
376 for (
int i = 0; i <=
p; i++)
378 dof_map[0*dof3 + 0 + ((
p - i) + j*(
p + 1))*(
p + 2)] = o++;
380 for (
int j = 0; j <=
p; j++)
381 for (
int i = 0; i <=
p; i++)
383 dof_map[2*dof3 + i + (j + (
p + 1)*(
p + 1))*(
p + 1)] = o++;
388 for (
int k = 0; k <=
p; k++)
389 for (
int j = 0; j <=
p; j++)
390 for (
int i = 1; i <=
p; i++)
392 dof_map[0*dof3 + i + (j + k*(
p + 1))*(
p + 2)] = o++;
395 for (
int k = 0; k <=
p; k++)
396 for (
int j = 1; j <=
p; j++)
397 for (
int i = 0; i <=
p; i++)
399 dof_map[1*dof3 + i + (j + k*(
p + 2))*(
p + 1)] = o++;
402 for (
int k = 1; k <=
p; k++)
403 for (
int j = 0; j <=
p; j++)
404 for (
int i = 0; i <=
p; i++)
406 dof_map[2*dof3 + i + (j + k*(
p + 1))*(
p + 1)] = o++;
414 for (
int k = 0; k <=
p; k++)
415 for (
int j = 0; j <=
p; j++)
416 for (
int i = 0; i <=
p/2; i++)
418 int idx = 0*dof3 + i + (j + k*(
p + 1))*(
p + 2);
422 for (
int k = 0; k <=
p; k++)
423 for (
int j = 0; j <=
p/2; j++)
424 for (
int i = 0; i <=
p; i++)
426 int idx = 1*dof3 + i + (j + k*(
p + 2))*(
p + 1);
430 for (
int k = 0; k <=
p/2; k++)
431 for (
int j = 0; j <=
p; j++)
432 for (
int i = 0; i <=
p; i++)
434 int idx = 2*dof3 + i + (j + k*(
p + 1))*(
p + 1);
440 for (
int k = 0; k <=
p; k++)
441 for (
int j = 0; j <=
p; j++)
442 for (
int i = 0; i <=
p + 1; i++)
457 for (
int k = 0; k <=
p; k++)
458 for (
int j = 0; j <=
p + 1; j++)
459 for (
int i = 0; i <=
p; i++)
474 for (
int k = 0; k <=
p + 1; k++)
475 for (
int j = 0; j <=
p; j++)
476 for (
int i = 0; i <=
p; i++)
495 const int pp1 =
order;
497 #ifdef MFEM_THREAD_SAFE 498 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
499 Vector shape_cz(pp1 + 1), shape_oz(pp1);
504 #ifdef MFEM_THREAD_SAFE 505 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
527 for (
int k = 0; k < pp1; k++)
528 for (
int j = 0; j < pp1; j++)
529 for (
int i = 0; i <= pp1; i++)
534 idx = -1 - idx,
s = -1;
540 shape(idx,0) =
s*shape_cx(i)*shape_oy(j)*shape_oz(k);
545 for (
int k = 0; k < pp1; k++)
546 for (
int j = 0; j <= pp1; j++)
547 for (
int i = 0; i < pp1; i++)
552 idx = -1 - idx,
s = -1;
559 shape(idx,1) =
s*shape_ox(i)*shape_cy(j)*shape_oz(k);
563 for (
int k = 0; k <= pp1; k++)
564 for (
int j = 0; j < pp1; j++)
565 for (
int i = 0; i < pp1; i++)
570 idx = -1 - idx,
s = -1;
578 shape(idx,2) =
s*shape_ox(i)*shape_oy(j)*shape_cz(k);
585 const int pp1 =
order;
587 #ifdef MFEM_THREAD_SAFE 588 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
589 Vector shape_cz(pp1 + 1), shape_oz(pp1);
590 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
612 for (
int k = 0; k < pp1; k++)
613 for (
int j = 0; j < pp1; j++)
614 for (
int i = 0; i <= pp1; i++)
619 idx = -1 - idx,
s = -1;
625 divshape(idx) =
s*dshape_cx(i)*shape_oy(j)*shape_oz(k);
628 for (
int k = 0; k < pp1; k++)
629 for (
int j = 0; j <= pp1; j++)
630 for (
int i = 0; i < pp1; i++)
635 idx = -1 - idx,
s = -1;
641 divshape(idx) =
s*shape_ox(i)*dshape_cy(j)*shape_oz(k);
644 for (
int k = 0; k <= pp1; k++)
645 for (
int j = 0; j < pp1; j++)
646 for (
int i = 0; i < pp1; i++)
651 idx = -1 - idx,
s = -1;
657 divshape(idx) =
s*shape_ox(i)*shape_oy(j)*dshape_cz(k);
675 for (
int c = 0; c < 3; c++)
680 for (
int k = 0; k < km; k++)
681 for (
int j = 0; j < jm; j++)
682 for (
int i = 0; i < im; i++)
685 if (idx < 0) { idx = -1 - idx; }
687 if (c == 0) { ic1 = j; ic2 = k; }
688 else if (c == 1) { ic1 = i; ic2 = k; }
689 else { ic1 = i; ic2 = j; }
690 const double h1 = cp[ic1+1] - cp[ic1];
691 const double h2 = cp[ic2+1] - cp[ic2];
693 for (
int q = 0; q < nqpt; q++)
696 if (c == 0) { ip3d.
Set3(cp[i], cp[j] + h1*ip2d.
x, cp[k] + h2*ip2d.
y); }
697 else if (c == 1) { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j], cp[k] + h2*ip2d.
y); }
698 else { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j] + h2*ip2d.
y, cp[k]); }
700 vc.
Eval(xq, Trans, ip3d);
706 dofs(idx) = val*h1*h2;
715 const int pp1 =
p + 1;
716 int n_face_dofs =
p*
p;
717 std::vector<int> strides, offsets;
718 const int n_dof_per_dim =
p*
p*pp1;
719 const auto f = internal::GetFaceNormal3D(face_id);
720 const int face_normal =
f.first, level =
f.second;
721 if (face_normal == 0)
723 offsets = {level ? pp1 - 1 : 0};
724 strides = {pp1,
p*pp1};
726 else if (face_normal == 1)
728 offsets = {n_dof_per_dim + (level ?
p*(pp1 - 1) : 0)};
729 strides = {1,
p*pp1};
731 else if (face_normal == 2)
733 offsets = {2*n_dof_per_dim + (level ?
p*
p*(pp1 - 1) : 0)};
736 std::vector<int> n_dofs = {
p,
p};
737 internal::FillFaceMap(n_face_dofs, offsets, strides, n_dofs, face_map);
741 const double RT_TriangleElement::nk[6] =
742 { 0., -1., 1., 1., -1., 0. };
744 const double RT_TriangleElement::c = 1./3.;
754 #ifndef MFEM_THREAD_SAFE 764 Vector shape_x(
p + 1), shape_y(
p + 1), shape_l(
p + 1);
769 for (
int i = 0; i <=
p; i++)
774 for (
int i = 0; i <=
p; i++)
779 for (
int i = 0; i <=
p; i++)
786 for (
int j = 0; j <
p; j++)
787 for (
int i = 0; i + j <
p; i++)
789 double w = iop[i] + iop[j] + iop[
p-1-i-j];
797 for (
int k = 0; k <
dof; k++)
803 const double *n_k = nk + 2*dof2nk[k];
806 for (
int j = 0; j <=
p; j++)
807 for (
int i = 0; i + j <=
p; i++)
809 double s = shape_x(i)*shape_y(j)*shape_l(
p-i-j);
810 T(o++, k) =
s*n_k[0];
811 T(o++, k) =
s*n_k[1];
813 for (
int i = 0; i <=
p; i++)
815 double s = shape_x(i)*shape_y(
p-i);
816 T(o++, k) =
s*((ip.
x - c)*n_k[0] + (ip.
y - c)*n_k[1]);
829 #ifdef MFEM_THREAD_SAFE 830 Vector shape_x(
p + 1), shape_y(
p + 1), shape_l(
p + 1);
839 for (
int j = 0; j <=
p; j++)
840 for (
int i = 0; i + j <=
p; i++)
842 double s = shape_x(i)*shape_y(j)*shape_l(
p-i-j);
843 u(o,0) =
s;
u(o,1) = 0; o++;
844 u(o,0) = 0;
u(o,1) =
s; o++;
846 for (
int i = 0; i <=
p; i++)
848 double s = shape_x(i)*shape_y(
p-i);
849 u(o,0) = (ip.
x - c)*
s;
850 u(o,1) = (ip.
y - c)*
s;
862 #ifdef MFEM_THREAD_SAFE 863 Vector shape_x(
p + 1), shape_y(
p + 1), shape_l(
p + 1);
864 Vector dshape_x(
p + 1), dshape_y(
p + 1), dshape_l(
p + 1);
873 for (
int j = 0; j <=
p; j++)
874 for (
int i = 0; i + j <=
p; i++)
877 divu(o++) = (dshape_x(i)*shape_l(k) -
878 shape_x(i)*dshape_l(k))*shape_y(j);
879 divu(o++) = (dshape_y(j)*shape_l(k) -
880 shape_y(j)*dshape_l(k))*shape_x(i);
882 for (
int i = 0; i <=
p; i++)
885 divu(o++) = ((shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j) +
886 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i));
889 Ti.
Mult(divu, divshape);
893 const double RT_TetrahedronElement::nk[12] =
894 { 1,1,1, -1,0,0, 0,-1,0, 0,0,-1 };
897 const double RT_TetrahedronElement::c = 1./4.;
907 #ifndef MFEM_THREAD_SAFE 919 Vector shape_x(
p + 1), shape_y(
p + 1), shape_z(
p + 1), shape_l(
p + 1);
925 for (
int j = 0; j <=
p; j++)
926 for (
int i = 0; i + j <=
p; i++)
928 double w = bop[i] + bop[j] + bop[
p-i-j];
932 for (
int j = 0; j <=
p; j++)
933 for (
int i = 0; i + j <=
p; i++)
935 double w = bop[i] + bop[j] + bop[
p-i-j];
939 for (
int j = 0; j <=
p; j++)
940 for (
int i = 0; i + j <=
p; i++)
942 double w = bop[i] + bop[j] + bop[
p-i-j];
946 for (
int j = 0; j <=
p; j++)
947 for (
int i = 0; i + j <=
p; i++)
949 double w = bop[i] + bop[j] + bop[
p-i-j];
955 for (
int k = 0; k <
p; k++)
956 for (
int j = 0; j + k <
p; j++)
957 for (
int i = 0; i + j + k <
p; i++)
959 double w = iop[i] + iop[j] + iop[k] + iop[
p-1-i-j-k];
969 for (
int m = 0; m <
dof; m++)
976 const double *nm = nk + 3*dof2nk[m];
979 for (
int k = 0; k <=
p; k++)
980 for (
int j = 0; j + k <=
p; j++)
981 for (
int i = 0; i + j + k <=
p; i++)
983 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(
p-i-j-k);
984 T(o++, m) =
s * nm[0];
985 T(o++, m) =
s * nm[1];
986 T(o++, m) =
s * nm[2];
988 for (
int j = 0; j <=
p; j++)
989 for (
int i = 0; i + j <=
p; i++)
991 double s = shape_x(i)*shape_y(j)*shape_z(
p-i-j);
992 T(o++, m) =
s*((ip.
x - c)*nm[0] + (ip.
y - c)*nm[1] +
1006 #ifdef MFEM_THREAD_SAFE 1007 Vector shape_x(
p + 1), shape_y(
p + 1), shape_z(
p + 1), shape_l(
p + 1);
1017 for (
int k = 0; k <=
p; k++)
1018 for (
int j = 0; j + k <=
p; j++)
1019 for (
int i = 0; i + j + k <=
p; i++)
1021 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(
p-i-j-k);
1022 u(o,0) =
s;
u(o,1) = 0;
u(o,2) = 0; o++;
1023 u(o,0) = 0;
u(o,1) =
s;
u(o,2) = 0; o++;
1024 u(o,0) = 0;
u(o,1) = 0;
u(o,2) =
s; o++;
1026 for (
int j = 0; j <=
p; j++)
1027 for (
int i = 0; i + j <=
p; i++)
1029 double s = shape_x(i)*shape_y(j)*shape_z(
p-i-j);
1030 u(o,0) = (ip.
x - c)*
s;
u(o,1) = (ip.
y - c)*
s;
u(o,2) = (ip.
z - c)*
s;
1042 #ifdef MFEM_THREAD_SAFE 1043 Vector shape_x(
p + 1), shape_y(
p + 1), shape_z(
p + 1), shape_l(
p + 1);
1044 Vector dshape_x(
p + 1), dshape_y(
p + 1), dshape_z(
p + 1), dshape_l(
p + 1);
1054 for (
int k = 0; k <=
p; k++)
1055 for (
int j = 0; j + k <=
p; j++)
1056 for (
int i = 0; i + j + k <=
p; i++)
1058 int l =
p - i - j - k;
1059 divu(o++) = (dshape_x(i)*shape_l(l) -
1060 shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
1061 divu(o++) = (dshape_y(j)*shape_l(l) -
1062 shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
1063 divu(o++) = (dshape_z(k)*shape_l(l) -
1064 shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
1066 for (
int j = 0; j <=
p; j++)
1067 for (
int i = 0; i + j <=
p; i++)
1071 (shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j)*shape_z(k) +
1072 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i)*shape_z(k) +
1073 (shape_z(k) + (ip.
z - c)*dshape_z(k))*shape_x(i)*shape_y(j);
1076 Ti.
Mult(divu, divshape);
1079 const double RT_WedgeElement::nk[15] =
1080 { 0,0,-1, 0,0,1, 0,-1,0, 1,1,0, -1,0,0};
1084 (
p + 2) * ((
p + 1) * (
p + 2)) / 2 +
1085 (
p + 1) * (
p + 1) * (
p + 3),
p + 1,
1095 MFEM_ASSERT(L2TriangleFE.
GetDof() * H1SegmentFE.
GetDof() +
1097 "Mismatch in number of degrees of freedom " 1098 "when building RT_WedgeElement!");
1100 const int pm1 =
p - 1;
1102 #ifndef MFEM_THREAD_SAFE 1120 for (
int j = 0; j <=
p; j++)
1121 for (
int i = 0; i + j <=
p; i++)
1123 l = j + i * (2 *
p + 3 - i) / 2;
1124 t_dof[o] = l; s_dof[o] = 0; dof2nk[o] = 0;
1131 for (
int j = 0; j <=
p; j++)
1132 for (
int i = 0; i + j <=
p; i++)
1134 t_dof[o] = l; s_dof[o] = 1; dof2nk[o] = 1; l++;
1140 for (
int j = 0; j <=
p; j++)
1141 for (
int i = 0; i <=
p; i++)
1143 t_dof[o] = i; s_dof[o] = j; dof2nk[o] = 2;
1149 for (
int j = 0; j <=
p; j++)
1150 for (
int i = 0; i <=
p; i++)
1152 t_dof[o] =
p + 1 + i; s_dof[o] = j; dof2nk[o] = 3;
1158 for (
int j = 0; j <=
p; j++)
1159 for (
int i = 0; i <=
p; i++)
1161 t_dof[o] = 2 *
p + 2 + i; s_dof[o] = j; dof2nk[o] = 4;
1168 for (
int k = 0; k < L2SegmentFE.
GetDof(); k++)
1171 for (
int j = 0; j <= pm1; j++)
1172 for (
int i = 0; i + j <= pm1; i++)
1174 t_dof[o] = 3 * (
p + 1) + 2 * l; s_dof[o] = k;
1180 t_dof[o] = 3 * (
p + 1) + 2 * l + 1; s_dof[o] = k;
1188 for (
int k = 2; k < H1SegmentFE.
GetDof(); k++)
1190 for (l = 0; l < L2TriangleFE.
GetDof(); l++)
1192 t_dof[o] = l; s_dof[o] = k; dof2nk[o] = 1;
1203 #ifdef MFEM_THREAD_SAFE 1217 for (
int i=0; i<
dof; i++)
1219 if ( dof2nk[i] >= 2 )
1221 shape(i, 0) = trt_shape(t_dof[i], 0) * sl2_shape[s_dof[i]];
1222 shape(i, 1) = trt_shape(t_dof[i], 1) * sl2_shape[s_dof[i]];
1227 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1230 shape(i, 2) =
s * tl2_shape[t_dof[i]] * sh1_shape(s_dof[i]);
1238 #ifdef MFEM_THREAD_SAFE 1253 for (
int i=0; i<
dof; i++)
1255 if ( dof2nk[i] >= 2 )
1257 divshape(i) = trt_dshape(t_dof[i]) * sl2_shape(s_dof[i]);
1261 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1262 divshape(i) =
s * tl2_shape(t_dof[i]) * sh1_dshape(s_dof[i], 0);
1267 const double RT_R1D_SegmentElement::nk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1275 cbasis1d(
poly1d.GetBasis(
p + 1, VerifyClosed(cb_type))),
1276 obasis1d(
poly1d.GetBasis(
p, VerifyOpen(ob_type)))
1284 #ifndef MFEM_THREAD_SAFE 1296 dof_map[0] = o; dof2nk[o++] = 0;
1300 dof_map[
p+1] = o; dof2nk[o++] = 0;
1304 for (
int i = 1; i <=
p; i++)
1307 dof_map[i] = o; dof2nk[o++] = 0;
1310 for (
int i = 0; i <=
p; i++)
1313 dof_map[
p+i+2] = o; dof2nk[o++] = 1;
1316 for (
int i = 0; i <=
p; i++)
1319 dof_map[2*
p+3+i] = o; dof2nk[o++] = 2;
1328 #ifdef MFEM_THREAD_SAFE 1329 Vector shape_cx(
p + 1), shape_ox(
p);
1332 cbasis1d.
Eval(ip.
x, shape_cx);
1333 obasis1d.
Eval(ip.
x, shape_ox);
1337 for (
int i = 0; i <=
p; i++)
1339 int idx = dof_map[o++];
1340 shape(idx,0) = shape_cx(i);
1345 for (
int i = 0; i <
p; i++)
1347 int idx = dof_map[o++];
1349 shape(idx,1) = shape_ox(i);
1353 for (
int i = 0; i <
p; i++)
1355 int idx = dof_map[o++];
1358 shape(idx,2) = shape_ox(i);
1367 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1368 "RT_R1D_SegmentElement cannot be embedded in " 1369 "2 or 3 dimensional spaces");
1370 for (
int i=0; i<
dof; i++)
1372 shape(i, 0) *= J(0,0);
1374 shape *= (1.0 / Trans.
Weight());
1382 #ifdef MFEM_THREAD_SAFE 1387 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
1391 for (
int i = 0; i <=
p; i++)
1393 int idx = dof_map[o++];
1394 divshape(idx) = dshape_cx(i);
1397 for (
int i = 0; i <
p; i++)
1399 int idx = dof_map[o++];
1403 for (
int i = 0; i <
p; i++)
1405 int idx = dof_map[o++];
1418 double * nk_ptr =
const_cast<double*
>(nk);
1420 for (
int k = 0; k <
dof; k++)
1426 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1427 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1430 Trans.
Weight() * vk3(1) * n3(1) +
1431 Trans.
Weight() * vk3(2) * n3(2);
1444 double * nk_ptr =
const_cast<double*
>(nk);
1447 for (
int k = 0; k <
dof; k++)
1451 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1452 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1459 vk[1] = n3[1] * Trans.
Weight();
1460 vk[2] = n3[2] * Trans.
Weight();
1463 double w = 1.0/Trans.
Weight();
1464 for (
int d = 0; d < 1; d++)
1470 for (
int j = 0; j < shape.Size(); j++)
1472 double s = shape(j);
1473 if (fabs(
s) < 1e-12)
1479 for (
int d = 0; d <
vdim; d++)
1481 I(k,j+d*shape.Size()) =
s*vk[d];
1491 double * nk_ptr =
const_cast<double*
>(nk);
1494 for (
int k = 0; k <
dof; k++)
1498 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1499 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1511 I(k, j) +=
vshape(j, 0) * vk[0];
1529 double * nk_ptr =
const_cast<double*
>(nk);
1532 for (
int k = 0; k <
dof; k++)
1535 curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
1536 for (
int j = 0; j < curl_k.Size(); j++)
1538 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1543 const double RT_R2D_SegmentElement::nk[2] = { 0.,1.};
1550 obasis1d(
poly1d.GetBasis(
p, VerifyOpen(ob_type)))
1557 #ifndef MFEM_THREAD_SAFE 1566 for (
int i = 0; i <=
p; i++)
1569 dof_map[i] = o; dof2nk[o++] = 0;
1578 #ifdef MFEM_THREAD_SAFE 1582 obasis1d.
Eval(ip.
x, shape_ox);
1586 for (
int i = 0; i <=
p; i++)
1588 int idx = dof_map[o++];
1589 shape(idx,0) = shape_ox(i);
1599 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1600 "RT_R2D_SegmentElement cannot be embedded in " 1601 "2 or 3 dimensional spaces");
1602 for (
int i=0; i<
dof; i++)
1604 shape(i, 0) *= J(0,0);
1606 shape *= (1.0 / Trans.
Weight());
1624 double * nk_ptr =
const_cast<double*
>(nk);
1631 for (
int k = 0; k <
dof; k++)
1633 Vector n2(&nk_ptr[dof2nk[k] * 2], 2);
1651 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1657 const double *nk_fe)
1673 MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
1674 "RT_R2D_FiniteElement cannot be embedded in " 1675 "3 dimensional spaces");
1676 for (
int i=0; i<
dof; i++)
1678 double sx = shape(i, 0);
1679 double sy = shape(i, 1);
1680 shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
1681 shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
1683 shape *= (1.0 / Trans.
Weight());
1696 double * nk_ptr =
const_cast<double*
>(
nk);
1703 for (
int k = 0; k <
dof; k++)
1717 for (
int i = 0; i <
dim; i++)
1719 Ikj +=
vshape(j, i) * vk[i];
1722 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1734 #ifdef MFEM_THREAD_SAFE 1738 double * nk_ptr =
const_cast<double*
>(
nk);
1742 const double weight = Trans.
Weight();
1743 for (
int j = 0; j <
dof; j++)
1755 for (
int k = 0; k <
dof; k++)
1758 for (
int d = 0; d <
dim; d++)
1760 R_jk +=
vshape(k,d)*pt_data[d];
1762 R_jk +=
vshape(k,2) * n3(2);
1783 double * nk_ptr =
const_cast<double*
>(
nk);
1785 for (
int k = 0; k <
dof; k++)
1795 Trans.
Weight() * vk3(2) * n3(2);
1808 double * nk_ptr =
const_cast<double*
>(
nk);
1811 for (
int k = 0; k <
dof; k++)
1823 vk[2] = n3[2] * Trans.
Weight();
1826 double w = 1.0/Trans.
Weight();
1827 for (
int d = 0; d < 2; d++)
1833 for (
int j = 0; j < shape.Size(); j++)
1835 double s = shape(j);
1836 if (fabs(
s) < 1e-12)
1842 for (
int d = 0; d <
vdim; d++)
1844 I(k,j+d*shape.Size()) =
s*vk[d];
1854 double * nk_ptr =
const_cast<double*
>(
nk);
1857 for (
int k = 0; k <
dof; k++)
1874 for (
int i=0; i<2; i++)
1876 I(k, j) +=
vshape(j, i) * vk[i];
1894 double * nk_ptr =
const_cast<double*
>(
nk);
1897 for (
int k = 0; k <
dof; k++)
1900 curl_shape.Mult(nk_ptr +
dof2nk[k] * 3, curl_k);
1901 for (
int j = 0; j < curl_k.Size(); j++)
1903 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1908 const double RT_R2D_TriangleElement::nk_t[12] =
1909 { 0.,-1.,0., 1.,1.,0., -1.,0.,0., 0.,0.,1. };
1918 #ifndef MFEM_THREAD_SAFE 1929 for (
int e=0; e<3; e++)
1932 for (
int i=0; i<=
p; i++)
1939 for (
int j = 0; j <
p; j++)
1940 for (
int i = 0; i + j <
p; i++)
1947 for (
int j = 0; j <=
p; j++)
1948 for (
int i = 0; i + j <=
p; i++)
1953 MFEM_VERIFY(r == RT_FE.
GetDof(),
1954 "RT_R2D_Triangle incorrect number of RT dofs.");
1955 MFEM_VERIFY(l == L2_FE.
GetDof(),
1956 "RT_R2D_Triangle incorrect number of L2 dofs.");
1957 MFEM_VERIFY(o ==
GetDof(),
1958 "RT_R2D_Triangle incorrect number of dofs.");
1963 for (
int i=0; i<
dof; i++)
1982 #ifdef MFEM_THREAD_SAFE 1990 for (
int i=0; i<
dof; i++)
1995 shape(i, 0) = rt_shape(idx, 0);
1996 shape(i, 1) = rt_shape(idx, 1);
2003 shape(i, 2) = l2_shape(-idx-1);
2011 #ifdef MFEM_THREAD_SAFE 2017 for (
int i=0; i<
dof; i++)
2022 div_shape(i) = rt_dshape(idx);
2031 const double RT_R2D_QuadrilateralElement::nk_q[15] =
2032 { 0., -1., 0., 1., 0., 0., 0., 1., 0., -1., 0., 0., 0., 0., 1. };
2038 cbasis1d(
poly1d.GetBasis(
p + 1, VerifyClosed(cb_type))),
2039 obasis1d(
poly1d.GetBasis(
p, VerifyOpen(ob_type)))
2043 const int dofx = (
p + 1)*(
p + 2);
2044 const int dofy = (
p + 1)*(
p + 2);
2045 const int dofxy = dofx + dofy;
2047 #ifndef MFEM_THREAD_SAFE 2058 for (
int i = 0; i <=
p; i++)
2060 dof_map[dofx + i + 0*(
p + 1)] = o++;
2062 for (
int i = 0; i <=
p; i++)
2066 for (
int i = 0; i <=
p; i++)
2068 dof_map[dofx + (
p - i) + (
p + 1)*(
p + 1)] = o++;
2070 for (
int i = 0; i <=
p; i++)
2076 for (
int j = 0; j <=
p; j++)
2077 for (
int i = 1; i <=
p; i++)
2081 for (
int j = 1; j <=
p; j++)
2082 for (
int i = 0; i <=
p; i++)
2084 dof_map[dofx + i + j*(
p + 1)] = o++;
2086 for (
int j = 0; j <=
p; j++)
2087 for (
int i = 0; i <=
p; i++)
2089 dof_map[dofxy + i + j*(
p + 1)] = o++;
2094 for (
int j = 0; j <=
p; j++)
2095 for (
int i = 0; i <=
p/2; i++)
2097 int idx = i + j*(
p + 2);
2101 for (
int j =
p/2 + 1; j <=
p; j++)
2103 int idx = (
p/2 + 1) + j*(
p + 2);
2107 for (
int j = 0; j <=
p/2; j++)
2108 for (
int i = 0; i <=
p; i++)
2110 int idx = dofx + i + j*(
p + 1);
2114 for (
int i = 0; i <=
p/2; i++)
2116 int idx = dofx + i + (
p/2 + 1)*(
p + 1);
2121 for (
int j = 0; j <=
p; j++)
2122 for (
int i = 0; i <=
p + 1; i++)
2136 for (
int j = 0; j <=
p + 1; j++)
2137 for (
int i = 0; i <=
p; i++)
2151 for (
int j = 0; j <=
p; j++)
2152 for (
int i = 0; i <=
p; i++)
2163 const int pp1 =
order;
2165 #ifdef MFEM_THREAD_SAFE 2166 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2169 cbasis1d.
Eval(ip.
x, shape_cx);
2170 obasis1d.
Eval(ip.
x, shape_ox);
2171 cbasis1d.
Eval(ip.
y, shape_cy);
2172 obasis1d.
Eval(ip.
y, shape_oy);
2175 for (
int j = 0; j < pp1; j++)
2176 for (
int i = 0; i <= pp1; i++)
2181 idx = -1 - idx,
s = -1;
2187 shape(idx,0) =
s*shape_cx(i)*shape_oy(j);
2191 for (
int j = 0; j <= pp1; j++)
2192 for (
int i = 0; i < pp1; i++)
2197 idx = -1 - idx,
s = -1;
2204 shape(idx,1) =
s*shape_ox(i)*shape_cy(j);
2207 for (
int j = 0; j < pp1; j++)
2208 for (
int i = 0; i < pp1; i++)
2213 shape(idx,2) = shape_ox(i)*shape_oy(j);
2220 const int pp1 =
order;
2222 #ifdef MFEM_THREAD_SAFE 2223 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2224 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
2227 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
2228 obasis1d.
Eval(ip.
x, shape_ox);
2229 cbasis1d.
Eval(ip.
y, shape_cy, dshape_cy);
2230 obasis1d.
Eval(ip.
y, shape_oy);
2233 for (
int j = 0; j < pp1; j++)
2234 for (
int i = 0; i <= pp1; i++)
2239 idx = -1 - idx,
s = -1;
2245 divshape(idx) =
s*dshape_cx(i)*shape_oy(j);
2247 for (
int j = 0; j <= pp1; j++)
2248 for (
int i = 0; i < pp1; i++)
2253 idx = -1 - idx,
s = -1;
2259 divshape(idx) =
s*shape_ox(i)*dshape_cy(j);
2261 for (
int j = 0; j < pp1; j++)
2262 for (
int i = 0; i < pp1; i++)
Abstract class for all finite elements.
void Set(const double *p, const int dim)
int GetNPoints() const
Returns the number of the points in the integration rule.
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Class for an integration rule - an Array of IntegrationPoint.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Base class for vector Coefficients that optionally depend on time and space.
void SetRow(int r, const double *row)
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
void SetSize(int s)
Resize the vector to size s.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
int dim
Dimension of reference space.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Data type dense matrix using column-major storage.
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &div_shape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Poly_1D::Basis & obasis1d
void Factor()
Factor the current DenseMatrix, *a.
Geometry::Type geom_type
Geometry::Type of the reference element.
Intermediate class for finite elements whose basis functions return vector values.
virtual void GetFaceMap(const int face_id, Array< int > &face_map) const
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in physical space at the point ...
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
RT_WedgeElement(const int p)
const double * ClosedPoints(const int p, const int btype=BasisType::GaussLobatto)
Get coordinates of a closed (GaussLegendre) set of points if degree p.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
int vdim
Vector dimension of vector-valued basis functions.
static void CalcBasis(const int p, const double x, double *u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
std::function< double(const Vector &)> f(double mass_coeff)
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
virtual void GetFaceMap(const int face_id, Array< int > &face_map) const
Return the mapping from lexicographically ordered face DOFs to lexicographically ordered element DOFs...
RT_R1D_SegmentElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R1D_SegmentElement of order p and closed and open BasisType cb_type and ob_type...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
RT_R2D_FiniteElement(int p, Geometry::Type G, int Do, const double *nk_fe)
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
void Set2(const double x1, const double x2)
int GetVDim()
Returns dimension of the vector.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
const double * OpenPoints(const int p, const int btype=BasisType::GaussLegendre)
Get coordinates of an open (GaussLegendre) set of points if degree p.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
void Threshold(double eps)
Replace small entries, abs(a_ij) <= eps, with zero.
RT_R2D_TriangleElement(const int p)
Construct the RT_R2D_TriangleElement of order p.
void Set3(const double x1, const double x2, const double x3)
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
double p(const Vector &x, double t)
RT_R2D_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
RT_R2D_SegmentElement(const int p, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R2D_SegmentElement of order p and open BasisType ob_type.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Class for integration point with weight.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
int dof
Number of degrees of freedom.
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
RT_TetrahedronElement(const int p)
Construct the RT_TetrahedronElement of order p.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
RT_TriangleElement(const int p)
Construct the RT_TriangleElement of order p.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Describes the function space on each element.
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
RT_HexahedronElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_HexahedronElement of order p and closed and open BasisType cb_type and ob_type...
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
double u(const Vector &xvec)
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
RT_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type...
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
int order
Order/degree of the shape functions.
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...