15 #include "../coefficient.hpp"
22 const double RT_QuadrilateralElement::nk[8] =
23 { 0., -1., 1., 0., 0., 1., -1., 0. };
31 cp(
poly1d.ClosedPoints(p + 1, cb_type))
38 const int dof2 =
dof/2;
40 #ifndef MFEM_THREAD_SAFE
51 for (
int i = 0; i <=
p; i++)
53 dof_map[1*dof2 + i + 0*(p + 1)] = o++;
55 for (
int i = 0; i <=
p; i++)
57 dof_map[0*dof2 + (p + 1) + i*(p + 2)] = o++;
59 for (
int i = 0; i <=
p; i++)
61 dof_map[1*dof2 + (p - i) + (p + 1)*(p + 1)] = o++;
63 for (
int i = 0; i <=
p; i++)
65 dof_map[0*dof2 + 0 + (p - i)*(p + 2)] = o++;
69 for (
int j = 0; j <=
p; j++)
70 for (
int i = 1; i <=
p; i++)
72 dof_map[0*dof2 + i + j*(p + 2)] = o++;
74 for (
int j = 1; j <=
p; j++)
75 for (
int i = 0; i <=
p; i++)
77 dof_map[1*dof2 + i + j*(p + 1)] = o++;
82 for (
int j = 0; j <=
p; j++)
83 for (
int i = 0; i <= p/2; i++)
85 int idx = 0*dof2 + i + j*(p + 2);
86 dof_map[idx] = -1 - dof_map[idx];
89 for (
int j = p/2 + 1; j <=
p; j++)
91 int idx = 0*dof2 + (p/2 + 1) + j*(p + 2);
95 for (
int j = 0; j <= p/2; j++)
96 for (
int i = 0; i <=
p; i++)
98 int idx = 1*dof2 + i + j*(p + 1);
99 dof_map[idx] = -1 - dof_map[idx];
102 for (
int i = 0; i <= p/2; i++)
104 int idx = 1*dof2 + i + (p/2 + 1)*(p + 1);
109 for (
int j = 0; j <=
p; j++)
110 for (
int i = 0; i <= p + 1; i++)
124 for (
int j = 0; j <= p + 1; j++)
125 for (
int i = 0; i <=
p; i++)
144 const int pp1 =
order;
146 #ifdef MFEM_THREAD_SAFE
147 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
148 Vector dshape_cx, dshape_cy;
167 for (
int j = 0; j < pp1; j++)
168 for (
int i = 0; i <= pp1; i++)
173 idx = -1 - idx, s = -1;
179 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
182 for (
int j = 0; j <= pp1; j++)
183 for (
int i = 0; i < pp1; i++)
188 idx = -1 - idx, s = -1;
195 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
202 const int pp1 =
order;
204 #ifdef MFEM_THREAD_SAFE
205 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
206 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
223 for (
int j = 0; j < pp1; j++)
224 for (
int i = 0; i <= pp1; i++)
229 idx = -1 - idx, s = -1;
235 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
237 for (
int j = 0; j <= pp1; j++)
238 for (
int i = 0; i < pp1; i++)
243 idx = -1 - idx, s = -1;
249 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
267 for (
int c = 0; c < 2; c++)
271 for (
int j = 0; j < jm; j++)
272 for (
int i = 0; i < im; i++)
275 if (idx < 0) { idx = -1 - idx; }
276 int ic = (c == 0) ? j : i;
277 const double h = cp[ic+1] - cp[ic];
279 for (
int k = 0; k < nqpt; k++)
282 if (c == 0) { ip2d.
Set2(cp[i], cp[j] + (h*ip1d.
x)); }
283 else { ip2d.
Set2(cp[i] + (h*ip1d.
x), cp[j]); }
285 vc.
Eval(xk, Trans, ip2d);
288 nk + dof2nk[idx]*
dim);
297 const double RT_HexahedronElement::nk[18] =
298 { 0.,0.,-1., 0.,-1.,0., 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,0.,1. };
306 cp(
poly1d.ClosedPoints(p + 1, cb_type))
313 const int dof3 =
dof/3;
315 #ifndef MFEM_THREAD_SAFE
329 for (
int j = 0; j <=
p; j++)
330 for (
int i = 0; i <=
p; i++)
332 dof_map[2*dof3 + i + ((p - j) + 0*(p + 1))*(p + 1)] = o++;
334 for (
int j = 0; j <=
p; j++)
335 for (
int i = 0; i <=
p; i++)
337 dof_map[1*dof3 + i + (0 + j*(p + 2))*(p + 1)] = o++;
339 for (
int j = 0; j <=
p; j++)
340 for (
int i = 0; i <=
p; i++)
342 dof_map[0*dof3 + (p + 1) + (i + j*(p + 1))*(p + 2)] = o++;
344 for (
int j = 0; j <=
p; j++)
345 for (
int i = 0; i <=
p; i++)
347 dof_map[1*dof3 + (p - i) + ((p + 1) + j*(p + 2))*(p + 1)] = o++;
349 for (
int j = 0; j <=
p; j++)
350 for (
int i = 0; i <=
p; i++)
352 dof_map[0*dof3 + 0 + ((p - i) + j*(p + 1))*(p + 2)] = o++;
354 for (
int j = 0; j <=
p; j++)
355 for (
int i = 0; i <=
p; i++)
357 dof_map[2*dof3 + i + (j + (p + 1)*(p + 1))*(p + 1)] = o++;
362 for (
int k = 0; k <=
p; k++)
363 for (
int j = 0; j <=
p; j++)
364 for (
int i = 1; i <=
p; i++)
366 dof_map[0*dof3 + i + (j + k*(p + 1))*(p + 2)] = o++;
369 for (
int k = 0; k <=
p; k++)
370 for (
int j = 1; j <=
p; j++)
371 for (
int i = 0; i <=
p; i++)
373 dof_map[1*dof3 + i + (j + k*(p + 2))*(p + 1)] = o++;
376 for (
int k = 1; k <=
p; k++)
377 for (
int j = 0; j <=
p; j++)
378 for (
int i = 0; i <=
p; i++)
380 dof_map[2*dof3 + i + (j + k*(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 = 0; i <= p/2; i++)
392 int idx = 0*dof3 + i + (j + k*(p + 1))*(p + 2);
393 dof_map[idx] = -1 - dof_map[idx];
396 for (
int k = 0; k <=
p; k++)
397 for (
int j = 0; j <= p/2; j++)
398 for (
int i = 0; i <=
p; i++)
400 int idx = 1*dof3 + i + (j + k*(p + 2))*(p + 1);
401 dof_map[idx] = -1 - dof_map[idx];
404 for (
int k = 0; k <= p/2; k++)
405 for (
int j = 0; j <=
p; j++)
406 for (
int i = 0; i <=
p; i++)
408 int idx = 2*dof3 + i + (j + k*(p + 1))*(p + 1);
409 dof_map[idx] = -1 - dof_map[idx];
414 for (
int k = 0; k <=
p; k++)
415 for (
int j = 0; j <=
p; j++)
416 for (
int i = 0; i <= p + 1; i++)
431 for (
int k = 0; k <=
p; k++)
432 for (
int j = 0; j <= p + 1; j++)
433 for (
int i = 0; i <=
p; i++)
448 for (
int k = 0; k <= p + 1; k++)
449 for (
int j = 0; j <=
p; j++)
450 for (
int i = 0; i <=
p; i++)
469 const int pp1 =
order;
471 #ifdef MFEM_THREAD_SAFE
472 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
473 Vector shape_cz(pp1 + 1), shape_oz(pp1);
474 Vector dshape_cx, dshape_cy, dshape_cz;
498 for (
int k = 0; k < pp1; k++)
499 for (
int j = 0; j < pp1; j++)
500 for (
int i = 0; i <= pp1; i++)
505 idx = -1 - idx, s = -1;
511 shape(idx,0) = s*shape_cx(i)*shape_oy(j)*shape_oz(k);
516 for (
int k = 0; k < pp1; k++)
517 for (
int j = 0; j <= pp1; j++)
518 for (
int i = 0; i < pp1; i++)
523 idx = -1 - idx, s = -1;
530 shape(idx,1) = s*shape_ox(i)*shape_cy(j)*shape_oz(k);
534 for (
int k = 0; k <= pp1; k++)
535 for (
int j = 0; j < pp1; j++)
536 for (
int i = 0; i < pp1; i++)
541 idx = -1 - idx, s = -1;
549 shape(idx,2) = s*shape_ox(i)*shape_oy(j)*shape_cz(k);
556 const int pp1 =
order;
558 #ifdef MFEM_THREAD_SAFE
559 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
560 Vector shape_cz(pp1 + 1), shape_oz(pp1);
561 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
582 for (
int k = 0; k < pp1; k++)
583 for (
int j = 0; j < pp1; j++)
584 for (
int i = 0; i <= pp1; i++)
589 idx = -1 - idx, s = -1;
595 divshape(idx) = s*dshape_cx(i)*shape_oy(j)*shape_oz(k);
598 for (
int k = 0; k < pp1; k++)
599 for (
int j = 0; j <= pp1; j++)
600 for (
int i = 0; i < pp1; i++)
605 idx = -1 - idx, s = -1;
611 divshape(idx) = s*shape_ox(i)*dshape_cy(j)*shape_oz(k);
614 for (
int k = 0; k <= pp1; k++)
615 for (
int j = 0; j < pp1; j++)
616 for (
int i = 0; i < pp1; i++)
621 idx = -1 - idx, s = -1;
627 divshape(idx) = s*shape_ox(i)*shape_oy(j)*dshape_cz(k);
645 for (
int c = 0; c < 3; c++)
650 for (
int k = 0; k < km; k++)
651 for (
int j = 0; j < jm; j++)
652 for (
int i = 0; i < im; i++)
655 if (idx < 0) { idx = -1 - idx; }
657 if (c == 0) { ic1 = j; ic2 = k; }
658 else if (c == 1) { ic1 = i; ic2 = k; }
659 else { ic1 = i; ic2 = j; }
660 const double h1 = cp[ic1+1] - cp[ic1];
661 const double h2 = cp[ic2+1] - cp[ic2];
663 for (
int q = 0; q < nqpt; q++)
666 if (c == 0) { ip3d.
Set3(cp[i], cp[j] + h1*ip2d.
x, cp[k] + h2*ip2d.
y); }
667 else if (c == 1) { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j], cp[k] + h2*ip2d.
y); }
668 else { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j] + h2*ip2d.
y, cp[k]); }
670 vc.
Eval(xq, Trans, ip3d);
676 dofs(idx) = val*h1*h2;
682 const double RT_TriangleElement::nk[6] =
683 { 0., -1., 1., 1., -1., 0. };
685 const double RT_TriangleElement::c = 1./3.;
695 #ifndef MFEM_THREAD_SAFE
705 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
710 for (
int i = 0; i <=
p; i++)
715 for (
int i = 0; i <=
p; i++)
720 for (
int i = 0; i <=
p; i++)
727 for (
int j = 0; j <
p; j++)
728 for (
int i = 0; i + j <
p; i++)
730 double w = iop[i] + iop[j] + iop[p-1-i-j];
738 for (
int k = 0; k <
dof; k++)
744 const double *n_k = nk + 2*dof2nk[k];
747 for (
int j = 0; j <=
p; j++)
748 for (
int i = 0; i + j <=
p; i++)
750 double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
751 T(o++, k) = s*n_k[0];
752 T(o++, k) = s*n_k[1];
754 for (
int i = 0; i <=
p; i++)
756 double s = shape_x(i)*shape_y(p-i);
757 T(o++, k) = s*((ip.
x - c)*n_k[0] + (ip.
y - c)*n_k[1]);
770 #ifdef MFEM_THREAD_SAFE
771 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
780 for (
int j = 0; j <=
p; j++)
781 for (
int i = 0; i + j <=
p; i++)
783 double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
784 u(o,0) =
s; u(o,1) = 0; o++;
785 u(o,0) = 0; u(o,1) =
s; o++;
787 for (
int i = 0; i <=
p; i++)
789 double s = shape_x(i)*shape_y(p-i);
790 u(o,0) = (ip.
x - c)*s;
791 u(o,1) = (ip.
y - c)*s;
803 #ifdef MFEM_THREAD_SAFE
804 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
805 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_l(p + 1);
814 for (
int j = 0; j <=
p; j++)
815 for (
int i = 0; i + j <=
p; i++)
818 divu(o++) = (dshape_x(i)*shape_l(k) -
819 shape_x(i)*dshape_l(k))*shape_y(j);
820 divu(o++) = (dshape_y(j)*shape_l(k) -
821 shape_y(j)*dshape_l(k))*shape_x(i);
823 for (
int i = 0; i <=
p; i++)
826 divu(o++) = ((shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j) +
827 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i));
830 Ti.
Mult(divu, divshape);
834 const double RT_TetrahedronElement::nk[12] =
835 { 1,1,1, -1,0,0, 0,-1,0, 0,0,-1 };
838 const double RT_TetrahedronElement::c = 1./4.;
848 #ifndef MFEM_THREAD_SAFE
860 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
866 for (
int j = 0; j <=
p; j++)
867 for (
int i = 0; i + j <=
p; i++)
869 double w = bop[i] + bop[j] + bop[p-i-j];
873 for (
int j = 0; j <=
p; j++)
874 for (
int i = 0; i + j <=
p; i++)
876 double w = bop[i] + bop[j] + bop[p-i-j];
880 for (
int j = 0; j <=
p; j++)
881 for (
int i = 0; i + j <=
p; i++)
883 double w = bop[i] + bop[j] + bop[p-i-j];
887 for (
int j = 0; j <=
p; j++)
888 for (
int i = 0; i + j <=
p; i++)
890 double w = bop[i] + bop[j] + bop[p-i-j];
896 for (
int k = 0; k <
p; k++)
897 for (
int j = 0; j + k <
p; j++)
898 for (
int i = 0; i + j + k <
p; i++)
900 double w = iop[i] + iop[j] + iop[k] + iop[p-1-i-j-k];
910 for (
int m = 0; m <
dof; m++)
917 const double *nm = nk + 3*dof2nk[m];
920 for (
int k = 0; k <=
p; k++)
921 for (
int j = 0; j + k <=
p; j++)
922 for (
int i = 0; i + j + k <=
p; i++)
924 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
925 T(o++, m) = s * nm[0];
926 T(o++, m) = s * nm[1];
927 T(o++, m) = s * nm[2];
929 for (
int j = 0; j <=
p; j++)
930 for (
int i = 0; i + j <=
p; i++)
932 double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
933 T(o++, m) = s*((ip.
x - c)*nm[0] + (ip.
y - c)*nm[1] +
947 #ifdef MFEM_THREAD_SAFE
948 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
958 for (
int k = 0; k <=
p; k++)
959 for (
int j = 0; j + k <=
p; j++)
960 for (
int i = 0; i + j + k <=
p; i++)
962 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
963 u(o,0) =
s; u(o,1) = 0; u(o,2) = 0; o++;
964 u(o,0) = 0; u(o,1) =
s; u(o,2) = 0; o++;
965 u(o,0) = 0; u(o,1) = 0; u(o,2) =
s; o++;
967 for (
int j = 0; j <=
p; j++)
968 for (
int i = 0; i + j <=
p; i++)
970 double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
971 u(o,0) = (ip.
x - c)*s; u(o,1) = (ip.
y - c)*s; u(o,2) = (ip.
z - c)*s;
983 #ifdef MFEM_THREAD_SAFE
984 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
985 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_z(p + 1), dshape_l(p + 1);
995 for (
int k = 0; k <=
p; k++)
996 for (
int j = 0; j + k <=
p; j++)
997 for (
int i = 0; i + j + k <=
p; i++)
999 int l = p - i - j - k;
1000 divu(o++) = (dshape_x(i)*shape_l(l) -
1001 shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
1002 divu(o++) = (dshape_y(j)*shape_l(l) -
1003 shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
1004 divu(o++) = (dshape_z(k)*shape_l(l) -
1005 shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
1007 for (
int j = 0; j <=
p; j++)
1008 for (
int i = 0; i + j <=
p; i++)
1012 (shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j)*shape_z(k) +
1013 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i)*shape_z(k) +
1014 (shape_z(k) + (ip.
z - c)*dshape_z(k))*shape_x(i)*shape_y(j);
1017 Ti.
Mult(divu, divshape);
1020 const double RT_WedgeElement::nk[15] =
1021 { 0,0,-1, 0,0,1, 0,-1,0, 1,1,0, -1,0,0};
1025 (p + 2) * ((p + 1) * (p + 2)) / 2 +
1026 (p + 1) * (p + 1) * (p + 3), p + 1,
1036 MFEM_ASSERT(L2TriangleFE.
GetDof() * H1SegmentFE.
GetDof() +
1038 "Mismatch in number of degrees of freedom "
1039 "when building RT_WedgeElement!");
1041 const int pm1 = p - 1;
1043 #ifndef MFEM_THREAD_SAFE
1061 for (
int j = 0; j <=
p; j++)
1062 for (
int i = 0; i + j <=
p; i++)
1064 l = j + i * (2 * p + 3 - i) / 2;
1065 t_dof[o] = l; s_dof[o] = 0; dof2nk[o] = 0;
1072 for (
int j = 0; j <=
p; j++)
1073 for (
int i = 0; i + j <=
p; i++)
1075 t_dof[o] = l; s_dof[o] = 1; dof2nk[o] = 1; l++;
1081 for (
int j = 0; j <=
p; j++)
1082 for (
int i = 0; i <=
p; i++)
1084 t_dof[o] = i; s_dof[o] = j; dof2nk[o] = 2;
1090 for (
int j = 0; j <=
p; j++)
1091 for (
int i = 0; i <=
p; i++)
1093 t_dof[o] = p + 1 + i; s_dof[o] = j; dof2nk[o] = 3;
1099 for (
int j = 0; j <=
p; j++)
1100 for (
int i = 0; i <=
p; i++)
1102 t_dof[o] = 2 * p + 2 + i; s_dof[o] = j; dof2nk[o] = 4;
1109 for (
int k = 0; k < L2SegmentFE.
GetDof(); k++)
1112 for (
int j = 0; j <= pm1; j++)
1113 for (
int i = 0; i + j <= pm1; i++)
1115 t_dof[o] = 3 * (p + 1) + 2 * l; s_dof[o] = k;
1121 t_dof[o] = 3 * (p + 1) + 2 * l + 1; s_dof[o] = k;
1129 for (
int k = 2; k < H1SegmentFE.
GetDof(); k++)
1131 for (l = 0; l < L2TriangleFE.
GetDof(); l++)
1133 t_dof[o] = l; s_dof[o] = k; dof2nk[o] = 1;
1144 #ifdef MFEM_THREAD_SAFE
1158 for (
int i=0; i<
dof; i++)
1160 if ( dof2nk[i] >= 2 )
1162 shape(i, 0) = trt_shape(t_dof[i], 0) * sl2_shape[s_dof[i]];
1163 shape(i, 1) = trt_shape(t_dof[i], 1) * sl2_shape[s_dof[i]];
1168 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1171 shape(i, 2) = s * tl2_shape[t_dof[i]] * sh1_shape(s_dof[i]);
1179 #ifdef MFEM_THREAD_SAFE
1180 Vector trt_dshape(RTSegmentFE.GetDof());
1194 for (
int i=0; i<
dof; i++)
1196 if ( dof2nk[i] >= 2 )
1198 divshape(i) = trt_dshape(t_dof[i]) * sl2_shape(s_dof[i]);
1202 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1203 divshape(i) = s * tl2_shape(t_dof[i]) * sh1_dshape(s_dof[i], 0);
1208 const double RT_R1D_SegmentElement::nk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1216 cbasis1d(
poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1217 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1225 #ifndef MFEM_THREAD_SAFE
1237 dof_map[0] = o; dof2nk[o++] = 0;
1241 dof_map[p+1] = o; dof2nk[o++] = 0;
1245 for (
int i = 1; i <=
p; i++)
1248 dof_map[i] = o; dof2nk[o++] = 0;
1251 for (
int i = 0; i <=
p; i++)
1254 dof_map[p+i+2] = o; dof2nk[o++] = 1;
1257 for (
int i = 0; i <=
p; i++)
1260 dof_map[2*p+3+i] = o; dof2nk[o++] = 2;
1269 #ifdef MFEM_THREAD_SAFE
1270 Vector shape_cx(p + 1), shape_ox(p);
1273 cbasis1d.
Eval(ip.
x, shape_cx);
1274 obasis1d.
Eval(ip.
x, shape_ox);
1278 for (
int i = 0; i <=
p; i++)
1280 int idx = dof_map[o++];
1281 shape(idx,0) = shape_cx(i);
1286 for (
int i = 0; i <
p; i++)
1288 int idx = dof_map[o++];
1290 shape(idx,1) = shape_ox(i);
1294 for (
int i = 0; i <
p; i++)
1296 int idx = dof_map[o++];
1299 shape(idx,2) = shape_ox(i);
1308 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1309 "RT_R1D_SegmentElement cannot be embedded in "
1310 "2 or 3 dimensional spaces");
1311 for (
int i=0; i<
dof; i++)
1313 shape(i, 0) *= J(0,0);
1315 shape *= (1.0 / Trans.
Weight());
1323 #ifdef MFEM_THREAD_SAFE
1328 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
1332 for (
int i = 0; i <=
p; i++)
1334 int idx = dof_map[o++];
1335 divshape(idx) = dshape_cx(i);
1338 for (
int i = 0; i <
p; i++)
1340 int idx = dof_map[o++];
1344 for (
int i = 0; i <
p; i++)
1346 int idx = dof_map[o++];
1359 double * nk_ptr =
const_cast<double*
>(nk);
1361 for (
int k = 0; k <
dof; k++)
1367 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1368 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1371 Trans.
Weight() * vk3(1) * n3(1) +
1372 Trans.
Weight() * vk3(2) * n3(2);
1385 double * nk_ptr =
const_cast<double*
>(nk);
1388 for (
int k = 0; k <
dof; k++)
1392 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1393 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1400 vk[1] = n3[1] * Trans.
Weight();
1401 vk[2] = n3[2] * Trans.
Weight();
1404 double w = 1.0/Trans.
Weight();
1405 for (
int d = 0; d < 1; d++)
1411 for (
int j = 0; j < shape.Size(); j++)
1413 double s = shape(j);
1414 if (fabs(s) < 1e-12)
1420 for (
int d = 0; d <
vdim; d++)
1422 I(k,j+d*shape.Size()) = s*vk[d];
1432 double * nk_ptr =
const_cast<double*
>(nk);
1435 for (
int k = 0; k <
dof; k++)
1439 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1440 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1452 I(k, j) +=
vshape(j, 0) * vk[0];
1470 double * nk_ptr =
const_cast<double*
>(nk);
1473 for (
int k = 0; k <
dof; k++)
1476 curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
1477 for (
int j = 0; j < curl_k.Size(); j++)
1479 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1484 const double RT_R2D_SegmentElement::nk[2] = { 0.,1.};
1491 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1498 #ifndef MFEM_THREAD_SAFE
1507 for (
int i = 0; i <=
p; i++)
1510 dof_map[i] = o; dof2nk[o++] = 0;
1519 #ifdef MFEM_THREAD_SAFE
1523 obasis1d.
Eval(ip.
x, shape_ox);
1527 for (
int i = 0; i <=
p; i++)
1529 int idx = dof_map[o++];
1530 shape(idx,0) = shape_ox(i);
1540 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1541 "RT_R2D_SegmentElement cannot be embedded in "
1542 "2 or 3 dimensional spaces");
1543 for (
int i=0; i<
dof; i++)
1545 shape(i, 0) *= J(0,0);
1547 shape *= (1.0 / Trans.
Weight());
1565 double * nk_ptr =
const_cast<double*
>(nk);
1572 for (
int k = 0; k <
dof; k++)
1574 Vector n2(&nk_ptr[dof2nk[k] * 2], 2);
1592 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1598 const double *nk_fe)
1614 MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
1615 "RT_R2D_FiniteElement cannot be embedded in "
1616 "3 dimensional spaces");
1617 for (
int i=0; i<
dof; i++)
1619 double sx = shape(i, 0);
1620 double sy = shape(i, 1);
1621 shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
1622 shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
1624 shape *= (1.0 / Trans.
Weight());
1637 double * nk_ptr =
const_cast<double*
>(
nk);
1644 for (
int k = 0; k <
dof; k++)
1658 for (
int i = 0; i <
dim; i++)
1660 Ikj +=
vshape(j, i) * vk[i];
1663 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1675 #ifdef MFEM_THREAD_SAFE
1679 double * nk_ptr =
const_cast<double*
>(
nk);
1683 const double weight = Trans.
Weight();
1684 for (
int j = 0; j <
dof; j++)
1696 for (
int k = 0; k <
dof; k++)
1699 for (
int d = 0; d <
dim; d++)
1701 R_jk +=
vshape(k,d)*pt_data[d];
1703 R_jk +=
vshape(k,2) * n3(2);
1724 double * nk_ptr =
const_cast<double*
>(
nk);
1726 for (
int k = 0; k <
dof; k++)
1736 Trans.
Weight() * vk3(2) * n3(2);
1749 double * nk_ptr =
const_cast<double*
>(
nk);
1752 for (
int k = 0; k <
dof; k++)
1764 vk[2] = n3[2] * Trans.
Weight();
1767 double w = 1.0/Trans.
Weight();
1768 for (
int d = 0; d < 2; d++)
1774 for (
int j = 0; j < shape.Size(); j++)
1776 double s = shape(j);
1777 if (fabs(s) < 1e-12)
1783 for (
int d = 0; d <
vdim; d++)
1785 I(k,j+d*shape.Size()) = s*vk[d];
1795 double * nk_ptr =
const_cast<double*
>(
nk);
1798 for (
int k = 0; k <
dof; k++)
1815 for (
int i=0; i<2; i++)
1817 I(k, j) +=
vshape(j, i) * vk[i];
1835 double * nk_ptr =
const_cast<double*
>(
nk);
1838 for (
int k = 0; k <
dof; k++)
1841 curl_shape.Mult(nk_ptr +
dof2nk[k] * 3, curl_k);
1842 for (
int j = 0; j < curl_k.Size(); j++)
1844 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1849 const double RT_R2D_TriangleElement::nk_t[12] =
1850 { 0.,-1.,0., 1.,1.,0., -1.,0.,0., 0.,0.,1. };
1859 #ifndef MFEM_THREAD_SAFE
1870 for (
int e=0; e<3; e++)
1873 for (
int i=0; i<=
p; i++)
1880 for (
int j = 0; j <
p; j++)
1881 for (
int i = 0; i + j <
p; i++)
1888 for (
int j = 0; j <=
p; j++)
1889 for (
int i = 0; i + j <=
p; i++)
1894 MFEM_VERIFY(r == RT_FE.
GetDof(),
1895 "RT_R2D_Triangle incorrect number of RT dofs.");
1896 MFEM_VERIFY(l == L2_FE.
GetDof(),
1897 "RT_R2D_Triangle incorrect number of L2 dofs.");
1898 MFEM_VERIFY(o ==
GetDof(),
1899 "RT_R2D_Triangle incorrect number of dofs.");
1904 for (
int i=0; i<
dof; i++)
1923 #ifdef MFEM_THREAD_SAFE
1931 for (
int i=0; i<
dof; i++)
1936 shape(i, 0) = rt_shape(idx, 0);
1937 shape(i, 1) = rt_shape(idx, 1);
1944 shape(i, 2) = l2_shape(-idx-1);
1952 #ifdef MFEM_THREAD_SAFE
1958 for (
int i=0; i<
dof; i++)
1963 div_shape(i) = rt_dshape(idx);
1972 const double RT_R2D_QuadrilateralElement::nk_q[15] =
1973 { 0., -1., 0., 1., 0., 0., 0., 1., 0., -1., 0., 0., 0., 0., 1. };
1979 cbasis1d(
poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1980 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1984 const int dofx = (p + 1)*(p + 2);
1985 const int dofy = (p + 1)*(p + 2);
1986 const int dofxy = dofx + dofy;
1988 #ifndef MFEM_THREAD_SAFE
1999 for (
int i = 0; i <=
p; i++)
2001 dof_map[dofx + i + 0*(p + 1)] = o++;
2003 for (
int i = 0; i <=
p; i++)
2005 dof_map[(p + 1) + i*(p + 2)] = o++;
2007 for (
int i = 0; i <=
p; i++)
2009 dof_map[dofx + (p - i) + (p + 1)*(p + 1)] = o++;
2011 for (
int i = 0; i <=
p; i++)
2013 dof_map[0 + (p - i)*(p + 2)] = o++;
2017 for (
int j = 0; j <=
p; j++)
2018 for (
int i = 1; i <=
p; i++)
2022 for (
int j = 1; j <=
p; j++)
2023 for (
int i = 0; i <=
p; i++)
2025 dof_map[dofx + i + j*(p + 1)] = o++;
2027 for (
int j = 0; j <=
p; j++)
2028 for (
int i = 0; i <=
p; i++)
2030 dof_map[dofxy + i + j*(p + 1)] = o++;
2035 for (
int j = 0; j <=
p; j++)
2036 for (
int i = 0; i <= p/2; i++)
2038 int idx = i + j*(p + 2);
2039 dof_map[idx] = -1 - dof_map[idx];
2042 for (
int j = p/2 + 1; j <=
p; j++)
2044 int idx = (p/2 + 1) + j*(p + 2);
2048 for (
int j = 0; j <= p/2; j++)
2049 for (
int i = 0; i <=
p; i++)
2051 int idx = dofx + i + j*(p + 1);
2052 dof_map[idx] = -1 - dof_map[idx];
2055 for (
int i = 0; i <= p/2; i++)
2057 int idx = dofx + i + (p/2 + 1)*(p + 1);
2062 for (
int j = 0; j <=
p; j++)
2063 for (
int i = 0; i <= p + 1; i++)
2077 for (
int j = 0; j <= p + 1; j++)
2078 for (
int i = 0; i <=
p; i++)
2092 for (
int j = 0; j <=
p; j++)
2093 for (
int i = 0; i <=
p; i++)
2104 const int pp1 =
order;
2106 #ifdef MFEM_THREAD_SAFE
2107 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2110 cbasis1d.
Eval(ip.
x, shape_cx);
2111 obasis1d.
Eval(ip.
x, shape_ox);
2112 cbasis1d.
Eval(ip.
y, shape_cy);
2113 obasis1d.
Eval(ip.
y, shape_oy);
2116 for (
int j = 0; j < pp1; j++)
2117 for (
int i = 0; i <= pp1; i++)
2122 idx = -1 - idx, s = -1;
2128 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
2132 for (
int j = 0; j <= pp1; j++)
2133 for (
int i = 0; i < pp1; i++)
2138 idx = -1 - idx, s = -1;
2145 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
2148 for (
int j = 0; j < pp1; j++)
2149 for (
int i = 0; i < pp1; i++)
2154 shape(idx,2) = shape_ox(i)*shape_oy(j);
2161 const int pp1 =
order;
2163 #ifdef MFEM_THREAD_SAFE
2164 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2165 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
2168 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
2169 obasis1d.
Eval(ip.
x, shape_ox);
2170 cbasis1d.
Eval(ip.
y, shape_cy, dshape_cy);
2171 obasis1d.
Eval(ip.
y, shape_oy);
2174 for (
int j = 0; j < pp1; j++)
2175 for (
int i = 0; i <= pp1; i++)
2180 idx = -1 - idx, s = -1;
2186 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
2188 for (
int j = 0; j <= pp1; j++)
2189 for (
int i = 0; i < pp1; i++)
2194 idx = -1 - idx, s = -1;
2200 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
2202 for (
int j = 0; j < pp1; j++)
2203 for (
int i = 0; i < pp1; i++)
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
void Set(const double *p, const int dim)
Class for an integration rule - an Array of IntegrationPoint.
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...
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
Base class for vector Coefficients that optionally depend on time and space.
void SetRow(int r, const double *row)
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 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 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.
void SetSize(int s)
Resize the vector to size s.
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 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...
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
int dim
Dimension of reference space.
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Data type dense matrix using column-major storage.
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.
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 GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
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 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 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...
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
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.
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
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 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...
Poly_1D::Basis & obasis1d
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...
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
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 ...
Poly_1D::Basis & cbasis1d
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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 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...
RT_R2D_FiniteElement(int p, Geometry::Type G, int Do, const double *nk_fe)
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...
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
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 Set2(const double x1, const double x2)
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 GetVDim()
Returns dimension of the vector.
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 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...
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.
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
void Set3(const double x1, const double x2, const double x3)
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 ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
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...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
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.
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
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 Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Class for integration point with weight.
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
int dof
Number of degrees of freedom.
RT_TetrahedronElement(const int p)
Construct the RT_TetrahedronElement of order p.
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...
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
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_TriangleElement(const int p)
Construct the RT_TriangleElement of order p.
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 Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
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.
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 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...
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
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...