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(pp1 + 1), dshape_cy(pp1 + 1);
168 for (
int j = 0; j < pp1; j++)
169 for (
int i = 0; i <= pp1; i++)
174 idx = -1 - idx, s = -1;
180 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
183 for (
int j = 0; j <= pp1; j++)
184 for (
int i = 0; i < pp1; i++)
189 idx = -1 - idx, s = -1;
196 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
203 const int pp1 =
order;
205 #ifdef MFEM_THREAD_SAFE
206 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
207 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
225 for (
int j = 0; j < pp1; j++)
226 for (
int i = 0; i <= pp1; i++)
231 idx = -1 - idx, s = -1;
237 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
239 for (
int j = 0; j <= pp1; j++)
240 for (
int i = 0; i < pp1; i++)
245 idx = -1 - idx, s = -1;
251 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
269 for (
int c = 0; c < 2; c++)
273 for (
int j = 0; j < jm; j++)
274 for (
int i = 0; i < im; i++)
277 if (idx < 0) { idx = -1 - idx; }
278 int ic = (c == 0) ? j : i;
279 const double h = cp[ic+1] - cp[ic];
281 for (
int k = 0; k < nqpt; k++)
284 if (c == 0) { ip2d.
Set2(cp[i], cp[j] + (h*ip1d.
x)); }
285 else { ip2d.
Set2(cp[i] + (h*ip1d.
x), cp[j]); }
287 vc.
Eval(xk, Trans, ip2d);
290 nk + dof2nk[idx]*
dim);
299 const double RT_HexahedronElement::nk[18] =
300 { 0.,0.,-1., 0.,-1.,0., 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,0.,1. };
308 cp(
poly1d.ClosedPoints(p + 1, cb_type))
315 const int dof3 =
dof/3;
317 #ifndef MFEM_THREAD_SAFE
331 for (
int j = 0; j <=
p; j++)
332 for (
int i = 0; i <=
p; i++)
334 dof_map[2*dof3 + i + ((p - j) + 0*(p + 1))*(p + 1)] = o++;
336 for (
int j = 0; j <=
p; j++)
337 for (
int i = 0; i <=
p; i++)
339 dof_map[1*dof3 + i + (0 + j*(p + 2))*(p + 1)] = o++;
341 for (
int j = 0; j <=
p; j++)
342 for (
int i = 0; i <=
p; i++)
344 dof_map[0*dof3 + (p + 1) + (i + j*(p + 1))*(p + 2)] = o++;
346 for (
int j = 0; j <=
p; j++)
347 for (
int i = 0; i <=
p; i++)
349 dof_map[1*dof3 + (p - i) + ((p + 1) + j*(p + 2))*(p + 1)] = o++;
351 for (
int j = 0; j <=
p; j++)
352 for (
int i = 0; i <=
p; i++)
354 dof_map[0*dof3 + 0 + ((p - i) + j*(p + 1))*(p + 2)] = o++;
356 for (
int j = 0; j <=
p; j++)
357 for (
int i = 0; i <=
p; i++)
359 dof_map[2*dof3 + i + (j + (p + 1)*(p + 1))*(p + 1)] = o++;
364 for (
int k = 0; k <=
p; k++)
365 for (
int j = 0; j <=
p; j++)
366 for (
int i = 1; i <=
p; i++)
368 dof_map[0*dof3 + i + (j + k*(p + 1))*(p + 2)] = o++;
371 for (
int k = 0; k <=
p; k++)
372 for (
int j = 1; j <=
p; j++)
373 for (
int i = 0; i <=
p; i++)
375 dof_map[1*dof3 + i + (j + k*(p + 2))*(p + 1)] = o++;
378 for (
int k = 1; k <=
p; k++)
379 for (
int j = 0; j <=
p; j++)
380 for (
int i = 0; i <=
p; i++)
382 dof_map[2*dof3 + i + (j + k*(p + 1))*(p + 1)] = o++;
390 for (
int k = 0; k <=
p; k++)
391 for (
int j = 0; j <=
p; j++)
392 for (
int i = 0; i <= p/2; i++)
394 int idx = 0*dof3 + i + (j + k*(p + 1))*(p + 2);
395 dof_map[idx] = -1 - dof_map[idx];
398 for (
int k = 0; k <=
p; k++)
399 for (
int j = 0; j <= p/2; j++)
400 for (
int i = 0; i <=
p; i++)
402 int idx = 1*dof3 + i + (j + k*(p + 2))*(p + 1);
403 dof_map[idx] = -1 - dof_map[idx];
406 for (
int k = 0; k <= p/2; k++)
407 for (
int j = 0; j <=
p; j++)
408 for (
int i = 0; i <=
p; i++)
410 int idx = 2*dof3 + i + (j + k*(p + 1))*(p + 1);
411 dof_map[idx] = -1 - dof_map[idx];
416 for (
int k = 0; k <=
p; k++)
417 for (
int j = 0; j <=
p; j++)
418 for (
int i = 0; i <= p + 1; i++)
433 for (
int k = 0; k <=
p; k++)
434 for (
int j = 0; j <= p + 1; j++)
435 for (
int i = 0; i <=
p; i++)
450 for (
int k = 0; k <= p + 1; k++)
451 for (
int j = 0; j <=
p; j++)
452 for (
int i = 0; i <=
p; i++)
471 const int pp1 =
order;
473 #ifdef MFEM_THREAD_SAFE
474 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
475 Vector shape_cz(pp1 + 1), shape_oz(pp1);
476 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
501 for (
int k = 0; k < pp1; k++)
502 for (
int j = 0; j < pp1; j++)
503 for (
int i = 0; i <= pp1; i++)
508 idx = -1 - idx, s = -1;
514 shape(idx,0) = s*shape_cx(i)*shape_oy(j)*shape_oz(k);
519 for (
int k = 0; k < pp1; k++)
520 for (
int j = 0; j <= pp1; j++)
521 for (
int i = 0; i < pp1; i++)
526 idx = -1 - idx, s = -1;
533 shape(idx,1) = s*shape_ox(i)*shape_cy(j)*shape_oz(k);
537 for (
int k = 0; k <= pp1; k++)
538 for (
int j = 0; j < pp1; j++)
539 for (
int i = 0; i < pp1; i++)
544 idx = -1 - idx, s = -1;
552 shape(idx,2) = s*shape_ox(i)*shape_oy(j)*shape_cz(k);
559 const int pp1 =
order;
561 #ifdef MFEM_THREAD_SAFE
562 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
563 Vector shape_cz(pp1 + 1), shape_oz(pp1);
564 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
586 for (
int k = 0; k < pp1; k++)
587 for (
int j = 0; j < pp1; j++)
588 for (
int i = 0; i <= pp1; i++)
593 idx = -1 - idx, s = -1;
599 divshape(idx) = s*dshape_cx(i)*shape_oy(j)*shape_oz(k);
602 for (
int k = 0; k < pp1; k++)
603 for (
int j = 0; j <= pp1; j++)
604 for (
int i = 0; i < pp1; i++)
609 idx = -1 - idx, s = -1;
615 divshape(idx) = s*shape_ox(i)*dshape_cy(j)*shape_oz(k);
618 for (
int k = 0; k <= pp1; k++)
619 for (
int j = 0; j < pp1; j++)
620 for (
int i = 0; i < pp1; i++)
625 idx = -1 - idx, s = -1;
631 divshape(idx) = s*shape_ox(i)*shape_oy(j)*dshape_cz(k);
649 for (
int c = 0; c < 3; c++)
654 for (
int k = 0; k < km; k++)
655 for (
int j = 0; j < jm; j++)
656 for (
int i = 0; i < im; i++)
659 if (idx < 0) { idx = -1 - idx; }
661 if (c == 0) { ic1 = j; ic2 = k; }
662 else if (c == 1) { ic1 = i; ic2 = k; }
663 else { ic1 = i; ic2 = j; }
664 const double h1 = cp[ic1+1] - cp[ic1];
665 const double h2 = cp[ic2+1] - cp[ic2];
667 for (
int q = 0; q < nqpt; q++)
670 if (c == 0) { ip3d.
Set3(cp[i], cp[j] + h1*ip2d.
x, cp[k] + h2*ip2d.
y); }
671 else if (c == 1) { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j], cp[k] + h2*ip2d.
y); }
672 else { ip3d.
Set3(cp[i] + h1*ip2d.
x, cp[j] + h2*ip2d.
y, cp[k]); }
674 vc.
Eval(xq, Trans, ip3d);
680 dofs(idx) = val*h1*h2;
686 const double RT_TriangleElement::nk[6] =
687 { 0., -1., 1., 1., -1., 0. };
689 const double RT_TriangleElement::c = 1./3.;
699 #ifndef MFEM_THREAD_SAFE
709 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
714 for (
int i = 0; i <=
p; i++)
719 for (
int i = 0; i <=
p; i++)
724 for (
int i = 0; i <=
p; i++)
731 for (
int j = 0; j <
p; j++)
732 for (
int i = 0; i + j <
p; i++)
734 double w = iop[i] + iop[j] + iop[p-1-i-j];
742 for (
int k = 0; k <
dof; k++)
748 const double *n_k = nk + 2*dof2nk[k];
751 for (
int j = 0; j <=
p; j++)
752 for (
int i = 0; i + j <=
p; i++)
754 double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
755 T(o++, k) = s*n_k[0];
756 T(o++, k) = s*n_k[1];
758 for (
int i = 0; i <=
p; i++)
760 double s = shape_x(i)*shape_y(p-i);
761 T(o++, k) = s*((ip.
x - c)*n_k[0] + (ip.
y - c)*n_k[1]);
774 #ifdef MFEM_THREAD_SAFE
775 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
784 for (
int j = 0; j <=
p; j++)
785 for (
int i = 0; i + j <=
p; i++)
787 double s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
788 u(o,0) =
s; u(o,1) = 0; o++;
789 u(o,0) = 0; u(o,1) =
s; o++;
791 for (
int i = 0; i <=
p; i++)
793 double s = shape_x(i)*shape_y(p-i);
794 u(o,0) = (ip.
x - c)*s;
795 u(o,1) = (ip.
y - c)*s;
807 #ifdef MFEM_THREAD_SAFE
808 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
809 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_l(p + 1);
818 for (
int j = 0; j <=
p; j++)
819 for (
int i = 0; i + j <=
p; i++)
822 divu(o++) = (dshape_x(i)*shape_l(k) -
823 shape_x(i)*dshape_l(k))*shape_y(j);
824 divu(o++) = (dshape_y(j)*shape_l(k) -
825 shape_y(j)*dshape_l(k))*shape_x(i);
827 for (
int i = 0; i <=
p; i++)
830 divu(o++) = ((shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j) +
831 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i));
834 Ti.
Mult(divu, divshape);
838 const double RT_TetrahedronElement::nk[12] =
839 { 1,1,1, -1,0,0, 0,-1,0, 0,0,-1 };
842 const double RT_TetrahedronElement::c = 1./4.;
852 #ifndef MFEM_THREAD_SAFE
864 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
870 for (
int j = 0; j <=
p; j++)
871 for (
int i = 0; i + j <=
p; i++)
873 double w = bop[i] + bop[j] + bop[p-i-j];
877 for (
int j = 0; j <=
p; j++)
878 for (
int i = 0; i + j <=
p; i++)
880 double w = bop[i] + bop[j] + bop[p-i-j];
884 for (
int j = 0; j <=
p; j++)
885 for (
int i = 0; i + j <=
p; i++)
887 double w = bop[i] + bop[j] + bop[p-i-j];
891 for (
int j = 0; j <=
p; j++)
892 for (
int i = 0; i + j <=
p; i++)
894 double w = bop[i] + bop[j] + bop[p-i-j];
900 for (
int k = 0; k <
p; k++)
901 for (
int j = 0; j + k <
p; j++)
902 for (
int i = 0; i + j + k <
p; i++)
904 double w = iop[i] + iop[j] + iop[k] + iop[p-1-i-j-k];
914 for (
int m = 0; m <
dof; m++)
921 const double *nm = nk + 3*dof2nk[m];
924 for (
int k = 0; k <=
p; k++)
925 for (
int j = 0; j + k <=
p; j++)
926 for (
int i = 0; i + j + k <=
p; i++)
928 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
929 T(o++, m) = s * nm[0];
930 T(o++, m) = s * nm[1];
931 T(o++, m) = s * nm[2];
933 for (
int j = 0; j <=
p; j++)
934 for (
int i = 0; i + j <=
p; i++)
936 double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
937 T(o++, m) = s*((ip.
x - c)*nm[0] + (ip.
y - c)*nm[1] +
951 #ifdef MFEM_THREAD_SAFE
952 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
962 for (
int k = 0; k <=
p; k++)
963 for (
int j = 0; j + k <=
p; j++)
964 for (
int i = 0; i + j + k <=
p; i++)
966 double s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
967 u(o,0) =
s; u(o,1) = 0; u(o,2) = 0; o++;
968 u(o,0) = 0; u(o,1) =
s; u(o,2) = 0; o++;
969 u(o,0) = 0; u(o,1) = 0; u(o,2) =
s; o++;
971 for (
int j = 0; j <=
p; j++)
972 for (
int i = 0; i + j <=
p; i++)
974 double s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
975 u(o,0) = (ip.
x - c)*s; u(o,1) = (ip.
y - c)*s; u(o,2) = (ip.
z - c)*s;
987 #ifdef MFEM_THREAD_SAFE
988 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
989 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_z(p + 1), dshape_l(p + 1);
999 for (
int k = 0; k <=
p; k++)
1000 for (
int j = 0; j + k <=
p; j++)
1001 for (
int i = 0; i + j + k <=
p; i++)
1003 int l = p - i - j - k;
1004 divu(o++) = (dshape_x(i)*shape_l(l) -
1005 shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
1006 divu(o++) = (dshape_y(j)*shape_l(l) -
1007 shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
1008 divu(o++) = (dshape_z(k)*shape_l(l) -
1009 shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
1011 for (
int j = 0; j <=
p; j++)
1012 for (
int i = 0; i + j <=
p; i++)
1016 (shape_x(i) + (ip.
x - c)*dshape_x(i))*shape_y(j)*shape_z(k) +
1017 (shape_y(j) + (ip.
y - c)*dshape_y(j))*shape_x(i)*shape_z(k) +
1018 (shape_z(k) + (ip.
z - c)*dshape_z(k))*shape_x(i)*shape_y(j);
1021 Ti.
Mult(divu, divshape);
1024 const double RT_WedgeElement::nk[15] =
1025 { 0,0,-1, 0,0,1, 0,-1,0, 1,1,0, -1,0,0};
1029 (p + 2) * ((p + 1) * (p + 2)) / 2 +
1030 (p + 1) * (p + 1) * (p + 3), p + 1,
1040 MFEM_ASSERT(L2TriangleFE.
GetDof() * H1SegmentFE.
GetDof() +
1042 "Mismatch in number of degrees of freedom "
1043 "when building RT_WedgeElement!");
1045 const int pm1 = p - 1;
1047 #ifndef MFEM_THREAD_SAFE
1065 for (
int j = 0; j <=
p; j++)
1066 for (
int i = 0; i + j <=
p; i++)
1068 l = j + i * (2 * p + 3 - i) / 2;
1069 t_dof[o] = l; s_dof[o] = 0; dof2nk[o] = 0;
1076 for (
int j = 0; j <=
p; j++)
1077 for (
int i = 0; i + j <=
p; i++)
1079 t_dof[o] = l; s_dof[o] = 1; dof2nk[o] = 1; l++;
1085 for (
int j = 0; j <=
p; j++)
1086 for (
int i = 0; i <=
p; i++)
1088 t_dof[o] = i; s_dof[o] = j; dof2nk[o] = 2;
1094 for (
int j = 0; j <=
p; j++)
1095 for (
int i = 0; i <=
p; i++)
1097 t_dof[o] = p + 1 + i; s_dof[o] = j; dof2nk[o] = 3;
1103 for (
int j = 0; j <=
p; j++)
1104 for (
int i = 0; i <=
p; i++)
1106 t_dof[o] = 2 * p + 2 + i; s_dof[o] = j; dof2nk[o] = 4;
1113 for (
int k = 0; k < L2SegmentFE.
GetDof(); k++)
1116 for (
int j = 0; j <= pm1; j++)
1117 for (
int i = 0; i + j <= pm1; i++)
1119 t_dof[o] = 3 * (p + 1) + 2 * l; s_dof[o] = k;
1125 t_dof[o] = 3 * (p + 1) + 2 * l + 1; s_dof[o] = k;
1133 for (
int k = 2; k < H1SegmentFE.
GetDof(); k++)
1135 for (l = 0; l < L2TriangleFE.
GetDof(); l++)
1137 t_dof[o] = l; s_dof[o] = k; dof2nk[o] = 1;
1148 #ifdef MFEM_THREAD_SAFE
1162 for (
int i=0; i<
dof; i++)
1164 if ( dof2nk[i] >= 2 )
1166 shape(i, 0) = trt_shape(t_dof[i], 0) * sl2_shape[s_dof[i]];
1167 shape(i, 1) = trt_shape(t_dof[i], 1) * sl2_shape[s_dof[i]];
1172 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1175 shape(i, 2) = s * tl2_shape[t_dof[i]] * sh1_shape(s_dof[i]);
1183 #ifdef MFEM_THREAD_SAFE
1198 for (
int i=0; i<
dof; i++)
1200 if ( dof2nk[i] >= 2 )
1202 divshape(i) = trt_dshape(t_dof[i]) * sl2_shape(s_dof[i]);
1206 double s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1207 divshape(i) = s * tl2_shape(t_dof[i]) * sh1_dshape(s_dof[i], 0);
1212 const double RT_R1D_SegmentElement::nk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1220 cbasis1d(
poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1221 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1229 #ifndef MFEM_THREAD_SAFE
1241 dof_map[0] = o; dof2nk[o++] = 0;
1245 dof_map[p+1] = o; dof2nk[o++] = 0;
1249 for (
int i = 1; i <=
p; i++)
1252 dof_map[i] = o; dof2nk[o++] = 0;
1255 for (
int i = 0; i <=
p; i++)
1258 dof_map[p+i+2] = o; dof2nk[o++] = 1;
1261 for (
int i = 0; i <=
p; i++)
1264 dof_map[2*p+3+i] = o; dof2nk[o++] = 2;
1273 #ifdef MFEM_THREAD_SAFE
1274 Vector shape_cx(p + 1), shape_ox(p);
1277 cbasis1d.
Eval(ip.
x, shape_cx);
1278 obasis1d.
Eval(ip.
x, shape_ox);
1282 for (
int i = 0; i <=
p; i++)
1284 int idx = dof_map[o++];
1285 shape(idx,0) = shape_cx(i);
1290 for (
int i = 0; i <
p; i++)
1292 int idx = dof_map[o++];
1294 shape(idx,1) = shape_ox(i);
1298 for (
int i = 0; i <
p; i++)
1300 int idx = dof_map[o++];
1303 shape(idx,2) = shape_ox(i);
1312 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1313 "RT_R1D_SegmentElement cannot be embedded in "
1314 "2 or 3 dimensional spaces");
1315 for (
int i=0; i<
dof; i++)
1317 shape(i, 0) *= J(0,0);
1319 shape *= (1.0 / Trans.
Weight());
1327 #ifdef MFEM_THREAD_SAFE
1332 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
1336 for (
int i = 0; i <=
p; i++)
1338 int idx = dof_map[o++];
1339 divshape(idx) = dshape_cx(i);
1342 for (
int i = 0; i <
p; i++)
1344 int idx = dof_map[o++];
1348 for (
int i = 0; i <
p; i++)
1350 int idx = dof_map[o++];
1363 double * nk_ptr =
const_cast<double*
>(nk);
1365 for (
int k = 0; k <
dof; k++)
1371 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1372 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1375 Trans.
Weight() * vk3(1) * n3(1) +
1376 Trans.
Weight() * vk3(2) * n3(2);
1389 double * nk_ptr =
const_cast<double*
>(nk);
1392 for (
int k = 0; k <
dof; k++)
1396 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1397 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1404 vk[1] = n3[1] * Trans.
Weight();
1405 vk[2] = n3[2] * Trans.
Weight();
1408 double w = 1.0/Trans.
Weight();
1409 for (
int d = 0; d < 1; d++)
1415 for (
int j = 0; j < shape.Size(); j++)
1417 double s = shape(j);
1418 if (fabs(s) < 1e-12)
1424 for (
int d = 0; d <
vdim; d++)
1426 I(k,j+d*shape.Size()) = s*vk[d];
1436 double * nk_ptr =
const_cast<double*
>(nk);
1439 for (
int k = 0; k <
dof; k++)
1443 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
1444 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
1456 I(k, j) +=
vshape(j, 0) * vk[0];
1474 double * nk_ptr =
const_cast<double*
>(nk);
1477 for (
int k = 0; k <
dof; k++)
1480 curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
1481 for (
int j = 0; j < curl_k.Size(); j++)
1483 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1488 const double RT_R2D_SegmentElement::nk[2] = { 0.,1.};
1495 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1502 #ifndef MFEM_THREAD_SAFE
1511 for (
int i = 0; i <=
p; i++)
1514 dof_map[i] = o; dof2nk[o++] = 0;
1523 #ifdef MFEM_THREAD_SAFE
1527 obasis1d.
Eval(ip.
x, shape_ox);
1531 for (
int i = 0; i <=
p; i++)
1533 int idx = dof_map[o++];
1534 shape(idx,0) = shape_ox(i);
1544 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
1545 "RT_R2D_SegmentElement cannot be embedded in "
1546 "2 or 3 dimensional spaces");
1547 for (
int i=0; i<
dof; i++)
1549 shape(i, 0) *= J(0,0);
1551 shape *= (1.0 / Trans.
Weight());
1569 double * nk_ptr =
const_cast<double*
>(nk);
1576 for (
int k = 0; k <
dof; k++)
1578 Vector n2(&nk_ptr[dof2nk[k] * 2], 2);
1596 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1602 const double *nk_fe)
1618 MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
1619 "RT_R2D_FiniteElement cannot be embedded in "
1620 "3 dimensional spaces");
1621 for (
int i=0; i<
dof; i++)
1623 double sx = shape(i, 0);
1624 double sy = shape(i, 1);
1625 shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
1626 shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
1628 shape *= (1.0 / Trans.
Weight());
1641 double * nk_ptr =
const_cast<double*
>(
nk);
1648 for (
int k = 0; k <
dof; k++)
1662 for (
int i = 0; i <
dim; i++)
1664 Ikj +=
vshape(j, i) * vk[i];
1667 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1679 #ifdef MFEM_THREAD_SAFE
1683 double * nk_ptr =
const_cast<double*
>(
nk);
1687 const double weight = Trans.
Weight();
1688 for (
int j = 0; j <
dof; j++)
1700 for (
int k = 0; k <
dof; k++)
1703 for (
int d = 0; d <
dim; d++)
1705 R_jk +=
vshape(k,d)*pt_data[d];
1707 R_jk +=
vshape(k,2) * n3(2);
1728 double * nk_ptr =
const_cast<double*
>(
nk);
1730 for (
int k = 0; k <
dof; k++)
1740 Trans.
Weight() * vk3(2) * n3(2);
1753 double * nk_ptr =
const_cast<double*
>(
nk);
1756 for (
int k = 0; k <
dof; k++)
1768 vk[2] = n3[2] * Trans.
Weight();
1771 double w = 1.0/Trans.
Weight();
1772 for (
int d = 0; d < 2; d++)
1778 for (
int j = 0; j < shape.Size(); j++)
1780 double s = shape(j);
1781 if (fabs(s) < 1e-12)
1787 for (
int d = 0; d <
vdim; d++)
1789 I(k,j+d*shape.Size()) = s*vk[d];
1799 double * nk_ptr =
const_cast<double*
>(
nk);
1802 for (
int k = 0; k <
dof; k++)
1819 for (
int i=0; i<2; i++)
1821 I(k, j) +=
vshape(j, i) * vk[i];
1839 double * nk_ptr =
const_cast<double*
>(
nk);
1842 for (
int k = 0; k <
dof; k++)
1845 curl_shape.Mult(nk_ptr +
dof2nk[k] * 3, curl_k);
1846 for (
int j = 0; j < curl_k.Size(); j++)
1848 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1853 const double RT_R2D_TriangleElement::nk_t[12] =
1854 { 0.,-1.,0., 1.,1.,0., -1.,0.,0., 0.,0.,1. };
1863 #ifndef MFEM_THREAD_SAFE
1874 for (
int e=0; e<3; e++)
1877 for (
int i=0; i<=
p; i++)
1884 for (
int j = 0; j <
p; j++)
1885 for (
int i = 0; i + j <
p; i++)
1892 for (
int j = 0; j <=
p; j++)
1893 for (
int i = 0; i + j <=
p; i++)
1898 MFEM_VERIFY(r == RT_FE.
GetDof(),
1899 "RT_R2D_Triangle incorrect number of RT dofs.");
1900 MFEM_VERIFY(l == L2_FE.
GetDof(),
1901 "RT_R2D_Triangle incorrect number of L2 dofs.");
1902 MFEM_VERIFY(o ==
GetDof(),
1903 "RT_R2D_Triangle incorrect number of dofs.");
1908 for (
int i=0; i<
dof; i++)
1927 #ifdef MFEM_THREAD_SAFE
1935 for (
int i=0; i<
dof; i++)
1940 shape(i, 0) = rt_shape(idx, 0);
1941 shape(i, 1) = rt_shape(idx, 1);
1948 shape(i, 2) = l2_shape(-idx-1);
1956 #ifdef MFEM_THREAD_SAFE
1962 for (
int i=0; i<
dof; i++)
1967 div_shape(i) = rt_dshape(idx);
1976 const double RT_R2D_QuadrilateralElement::nk_q[15] =
1977 { 0., -1., 0., 1., 0., 0., 0., 1., 0., -1., 0., 0., 0., 0., 1. };
1983 cbasis1d(
poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
1984 obasis1d(
poly1d.GetBasis(p, VerifyOpen(ob_type)))
1988 const int dofx = (p + 1)*(p + 2);
1989 const int dofy = (p + 1)*(p + 2);
1990 const int dofxy = dofx + dofy;
1992 #ifndef MFEM_THREAD_SAFE
2003 for (
int i = 0; i <=
p; i++)
2005 dof_map[dofx + i + 0*(p + 1)] = o++;
2007 for (
int i = 0; i <=
p; i++)
2009 dof_map[(p + 1) + i*(p + 2)] = o++;
2011 for (
int i = 0; i <=
p; i++)
2013 dof_map[dofx + (p - i) + (p + 1)*(p + 1)] = o++;
2015 for (
int i = 0; i <=
p; i++)
2017 dof_map[0 + (p - i)*(p + 2)] = o++;
2021 for (
int j = 0; j <=
p; j++)
2022 for (
int i = 1; i <=
p; i++)
2026 for (
int j = 1; j <=
p; j++)
2027 for (
int i = 0; i <=
p; i++)
2029 dof_map[dofx + i + j*(p + 1)] = o++;
2031 for (
int j = 0; j <=
p; j++)
2032 for (
int i = 0; i <=
p; i++)
2034 dof_map[dofxy + i + j*(p + 1)] = o++;
2039 for (
int j = 0; j <=
p; j++)
2040 for (
int i = 0; i <= p/2; i++)
2042 int idx = i + j*(p + 2);
2043 dof_map[idx] = -1 - dof_map[idx];
2046 for (
int j = p/2 + 1; j <=
p; j++)
2048 int idx = (p/2 + 1) + j*(p + 2);
2052 for (
int j = 0; j <= p/2; j++)
2053 for (
int i = 0; i <=
p; i++)
2055 int idx = dofx + i + j*(p + 1);
2056 dof_map[idx] = -1 - dof_map[idx];
2059 for (
int i = 0; i <= p/2; i++)
2061 int idx = dofx + i + (p/2 + 1)*(p + 1);
2066 for (
int j = 0; j <=
p; j++)
2067 for (
int i = 0; i <= p + 1; i++)
2081 for (
int j = 0; j <= p + 1; j++)
2082 for (
int i = 0; i <=
p; i++)
2096 for (
int j = 0; j <=
p; j++)
2097 for (
int i = 0; i <=
p; i++)
2108 const int pp1 =
order;
2110 #ifdef MFEM_THREAD_SAFE
2111 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2114 cbasis1d.
Eval(ip.
x, shape_cx);
2115 obasis1d.
Eval(ip.
x, shape_ox);
2116 cbasis1d.
Eval(ip.
y, shape_cy);
2117 obasis1d.
Eval(ip.
y, shape_oy);
2120 for (
int j = 0; j < pp1; j++)
2121 for (
int i = 0; i <= pp1; i++)
2126 idx = -1 - idx, s = -1;
2132 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
2136 for (
int j = 0; j <= pp1; j++)
2137 for (
int i = 0; i < pp1; i++)
2142 idx = -1 - idx, s = -1;
2149 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
2152 for (
int j = 0; j < pp1; j++)
2153 for (
int i = 0; i < pp1; i++)
2158 shape(idx,2) = shape_ox(i)*shape_oy(j);
2165 const int pp1 =
order;
2167 #ifdef MFEM_THREAD_SAFE
2168 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2169 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
2172 cbasis1d.
Eval(ip.
x, shape_cx, dshape_cx);
2173 obasis1d.
Eval(ip.
x, shape_ox);
2174 cbasis1d.
Eval(ip.
y, shape_cy, dshape_cy);
2175 obasis1d.
Eval(ip.
y, shape_oy);
2178 for (
int j = 0; j < pp1; j++)
2179 for (
int i = 0; i <= pp1; i++)
2184 idx = -1 - idx, s = -1;
2190 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
2192 for (
int j = 0; j <= pp1; j++)
2193 for (
int i = 0; i < pp1; i++)
2198 idx = -1 - idx, s = -1;
2204 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
2206 for (
int j = 0; j < pp1; j++)
2207 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...
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.
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...