15 #include "../coefficient.hpp"
34 #ifndef MFEM_THREAD_SAFE
42 MFEM_ABORT(
"method is not implemented for this class");
48 MFEM_ABORT(
"method is not implemented for this class");
54 MFEM_ABORT(
"method is not implemented for this class");
61 div_shape *= (1.0 / Trans.
Weight());
67 MFEM_ABORT(
"method is not implemented for this class");
77 #ifdef MFEM_THREAD_SAFE
82 curl_shape *= (1.0 / Trans.
Weight());
88 curl_shape *= (1.0 / Trans.
Weight());
91 MFEM_ABORT(
"Invalid dimension, Dim = " <<
dim);
97 MFEM_ABORT(
"method is not overloaded");
103 MFEM_ABORT(
"method is not overloaded");
109 MFEM_ABORT(
"method is not overloaded");
115 MFEM_ABORT(
"method is not overloaded");
122 MFEM_ABORT(
"method is not overloaded");
128 MFEM_ABORT(
"method is not overloaded");
134 MFEM_ABORT(
"method is not overloaded");
140 mfem_error (
"FiniteElement::ProjectFromNodes() (vector) is not overloaded!");
146 MFEM_ABORT(
"method is not overloaded");
151 MFEM_ABORT(
"method is not implemented for this element");
157 MFEM_ABORT(
"method is not implemented for this element");
164 MFEM_ABORT(
"method is not implemented for this element");
171 MFEM_ABORT(
"method is not implemented for this element");
178 MFEM_ABORT(
"method is not implemented for this element");
195 #ifdef MFEM_THREAD_SAFE
215 int size = (
dim*(
dim+1))/2;
221 for (
int nd = 0; nd <
dof; nd++)
223 Laplacian[nd] = hess(nd,0) + hess(nd,4) + hess(nd,5);
228 for (
int nd = 0; nd <
dof; nd++)
230 Laplacian[nd] = hess(nd,0) + hess(nd,2);
235 for (
int nd = 0; nd <
dof; nd++)
237 Laplacian[nd] = hess(nd,0);
248 int size = (
dim*(
dim+1))/2;
259 scale[1] = 2*Gij(0,1);
260 scale[2] = 2*Gij(0,2);
262 scale[3] = 2*Gij(1,2);
270 scale[1] = 2*Gij(0,1);
278 for (
int nd = 0; nd <
dof; nd++)
281 for (
int ii = 0; ii < size; ii++)
283 Laplacian[nd] += hess(nd,ii)*scale[ii];
324 int size = (dim*(dim+1))/2;
342 for (
int i = 0; i <
dim; i++)
344 for (
int j = 0; j <
dim; j++)
346 for (
int k = 0; k <
dim; k++)
348 for (
int l = 0; l <
dim; l++)
350 lhm(map[i*dim+j],map[k*dim+l]) += invJ(i,k)*invJ(j,l);
358 for (
int i = 0; i < dim*
dim; i++) { mult[map[i]]++; }
363 Mult( hess, lhm, Hessian);
369 MFEM_ABORT(
"method is not implemented for this element");
390 #ifdef MFEM_THREAD_SAFE
397 for (
int i = 0; i < fine_fe.
dof; i++)
402 for (
int j = 0; j <
dof; j++)
404 if (fabs(I(i,j) =
c_shape(j)) < 1.0e-12)
430 Vector fine_shape(fs), coarse_shape(cs);
431 DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs);
449 fine_mass_inv.
Mult(fine_coarse_mass, I);
470 Vector fine_shape(fs), coarse_shape(cs);
471 DenseMatrix coarse_mass(cs), coarse_fine_mass(cs, fs);
488 coarse_mass_inv.
Mult(coarse_fine_mass, R);
494 R *= 1.0 / Trans.
Weight();
505 if (d2q.
IntRule == &ir && d2q.
mode == mode) {
return d2q; }
519 #ifdef MFEM_THREAD_SAFE
523 for (
int i = 0; i < nqpt; i++)
527 for (
int j = 0; j <
dof; j++)
529 d2q->
B[i+nqpt*j] = d2q->
Bt[j+dof*i] =
c_shape(j);
532 for (
int d = 0; d <
dim; d++)
534 for (
int j = 0; j <
dof; j++)
536 d2q->
G[i+nqpt*(d+dim*j)] = d2q->
Gt[j+dof*(i+nqpt*d)] =
vshape(j,d);
554 if (d2q.
IntRule == &ir && d2q.
mode == mode) {
return d2q; }
559 const int ndof =
order + 1;
560 const int nqpt = (int)floor(pow(ir.
GetNPoints(), 1.0/
dim) + 0.5);
570 Vector val(ndof), grad(ndof);
571 for (
int i = 0; i < nqpt; i++)
576 for (
int j = 0; j < ndof; j++)
578 d2q->
B[i+nqpt*j] = d2q->
Bt[j+ndof*i] = val(j);
579 d2q->
G[i+nqpt*j] = d2q->
Gt[j+ndof*i] = grad(j);
594 for (
int i = 0; i <
dof; i++)
604 for (
int j = 0; j < fe.
GetDof(); j++)
606 curl(i,j) = w * curl_shape(j,0);
633 #ifdef MFEM_THREAD_SAFE
639 for (
int j = 0; j <
dof; j++)
659 for (
int i = 0; i <
dof; i++)
665 dofs(i) = coeff.
Eval (Trans, ip);
668 dofs(i) *= Trans.
Weight();
679 for (
int i = 0; i <
dof; i++)
683 vc.
Eval (x, Trans, ip);
688 for (
int j = 0; j < x.Size(); j++)
690 dofs(dof*j+i) = x(j);
702 for (
int k = 0; k <
dof; k++)
707 for (
int r = 0; r < MQ.Height(); r++)
709 for (
int d = 0; d < MQ.Width(); d++)
711 dofs(k+dof*(d+MQ.Width()*r)) = MQ(r,d);
727 for (
int k = 0; k <
dof; k++)
730 for (
int j = 0; j < shape.Size(); j++)
732 I(k,j) = (fabs(shape(j)) < 1e-12) ? 0.0 : shape(j);
738 for (
int k = 0; k <
dof; k++)
746 for (
int j = 0; j < shape.Size(); j++)
748 I(k,j) = (fabs(shape(j)) < 1e-12) ? 0.0 : shape(j);
759 for (
int k = 0; k <
dof; k++)
767 for (
int j = 0; j < vshape.
Height(); j++)
768 for (
int d = 0; d < vshape.
Width(); d++)
770 I(k+d*dof,j) =
vshape(j,d);
786 for (
int k = 0; k <
dof; k++)
792 Mult(dshape, Jinv, grad_k);
797 for (
int j = 0; j < grad_k.Height(); j++)
798 for (
int d = 0; d <
dim; d++)
800 grad(k+d*dof,j) = grad_k(j,d);
813 for (
int k = 0; k <
dof; k++)
821 for (
int j = 0; j < div_shape.Size(); j++)
823 div(k,j) = (fabs(div_shape(j)) < 1e-12) ? 0.0 : div_shape(j)/detJ;
828 for (
int j = 0; j < div_shape.Size(); j++)
830 div(k,j) = (fabs(div_shape(j)) < 1e-12) ? 0.0 : div_shape(j);
838 int Do,
int O,
int M,
int F)
852 void VectorFiniteElement::CalcShape (
855 mfem_error (
"Error: Cannot use scalar CalcShape(...) function with\n"
856 " VectorFiniteElements!");
859 void VectorFiniteElement::CalcDShape (
860 const IntegrationPoint &ip, DenseMatrix &dshape )
const
862 mfem_error (
"Error: Cannot use scalar CalcDShape(...) function with\n"
863 " VectorFiniteElements!");
895 MFEM_ABORT(
"Invalid dimension, Dim = " <<
dim);
899 MFEM_ABORT(
"Invalid MapType = " <<
map_type);
907 #ifdef MFEM_THREAD_SAFE
912 shape *= (1.0 / Trans.
Weight());
919 #ifdef MFEM_THREAD_SAFE
932 MFEM_ASSERT(vc.
GetVDim() == sdim,
"");
934 const bool square_J = (
dim == sdim);
936 for (
int k = 0; k <
dof; k++)
942 if (!square_J) { dofs(k) /= Trans.
Weight(); }
951 const bool square_J = (
dim == sdim);
953 for (
int k = 0; k <
dof; k++)
959 if (!square_J) { dofs(k) /= Trans.
Weight(); }
970 MFEM_ASSERT(mc.
GetWidth() == sdim,
"");
971 const bool square_J = (
dim == sdim);
973 Vector nk_phys(sdim), dofs_k(MQ.Height());
974 MFEM_ASSERT(dofs.
Size() ==
dof*MQ.Height(),
"");
976 for (
int k = 0; k <
dof; k++)
982 if (!square_J) { nk_phys /= T.
Weight(); }
983 MQ.Mult(nk_phys, dofs_k);
984 for (
int r = 0; r < MQ.Height(); r++)
986 dofs(k+dof*r) = dofs_k(r);
1002 for (
int k = 0; k <
dof; k++)
1013 double w = 1.0/Trans.
Weight();
1014 for (
int d = 0; d <
dim; d++)
1020 for (
int j = 0; j < shape.Size(); j++)
1022 double s = shape(j);
1023 if (fabs(s) < 1e-12)
1029 for (
int d = 0; d < sdim; d++)
1031 I(k,j+d*shape.Size()) = s*vk[d];
1042 const bool square_J = (
dim == sdim);
1045 for (
int k = 0; k <
dof; k++)
1057 if (!square_J) { vshapenk /= Trans.
Weight(); }
1058 for (
int j=0; j<vshapenk.Size(); j++)
1060 I(k,j) = vshapenk(j);
1072 mfem_error(
"VectorFiniteElement::ProjectGrad_RT works only in 2D!");
1080 for (
int k = 0; k <
dof; k++)
1083 tk[0] = nk[d2n[k]*
dim+1];
1084 tk[1] = -nk[d2n[k]*
dim];
1085 dshape.Mult(tk, grad_k);
1086 for (
int j = 0; j < grad_k.Size(); j++)
1088 grad(k,j) = (fabs(grad_k(j)) < 1e-12) ? 0.0 : grad_k(j);
1097 #ifdef MFEM_THREAD_SAFE
1110 for (
int k = 0; k <
dof; k++)
1117 JtJ *= 1.0 / Trans.
Weight();
1124 for (
int j = 0; j < curl_k.Size(); j++)
1126 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1139 for (
int k = 0; k <
dof; k++)
1142 curl_shape.Mult(nk + d2n[k]*
dim, curl_k);
1143 for (
int j = 0; j < curl_k.Size(); j++)
1145 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1157 for (
int k = 0; k <
dof; k++)
1171 for (
int k = 0; k <
dof; k++)
1187 MFEM_ASSERT(mc.
GetWidth() == sdim,
"");
1189 Vector tk_phys(sdim), dofs_k(MQ.Height());
1190 MFEM_ASSERT(dofs.
Size() ==
dof*MQ.Height(),
"");
1192 for (
int k = 0; k <
dof; k++)
1198 MQ.Mult(tk_phys, dofs_k);
1199 for (
int r = 0; r < MQ.Height(); r++)
1201 dofs(k+dof*r) = dofs_k(r);
1217 for (
int k = 0; k <
dof; k++)
1228 double w = 1.0/Trans.
Weight();
1229 for (
int d = 0; d < sdim; d++)
1235 for (
int j = 0; j < shape.Size(); j++)
1237 double s = shape(j);
1238 if (fabs(s) < 1e-12)
1244 for (
int d = 0; d < sdim; d++)
1246 I(k, j + d*shape.Size()) = s*vk[d];
1259 for (
int k = 0; k <
dof; k++)
1271 for (
int j=0; j<vshapetk.Size(); j++)
1273 I(k, j) = vshapetk(j);
1289 for (
int k = 0; k <
dof; k++)
1292 dshape.Mult(tk + d2t[k]*
dim, grad_k);
1293 for (
int j = 0; j < grad_k.Size(); j++)
1295 grad(k,j) = (fabs(grad_k(j)) < 1e-12) ? 0.0 : grad_k(j);
1310 DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs);
1311 const int ir_order =
1327 for (
int k=0; k<fs; ++k)
1329 for (
int j=0; j<cs; ++j)
1332 for (
int d1=0; d1<
dim; ++d1)
1334 for (
int d2=0; d2<
dim; ++d2)
1336 Mkj += w*fine_shape(k,d1)*adjJ(d2,d1)*coarse_shape(j,d2);
1339 fine_coarse_mass(k,j) += (fabs(Mkj) < 1e-12) ? 0.0 : Mkj;
1344 fine_mass_inv.
Mult(fine_coarse_mass, I);
1358 #ifdef MFEM_THREAD_SAFE
1368 for (
int k = 0; k <
dof; k++)
1379 for (
int i = 0; i <
dim; i++)
1381 Ikj +=
vshape(j, i) * vk[i];
1383 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1398 DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs);
1399 const int ir_order =
1414 for (
int k=0; k<fs; ++k)
1416 for (
int j=0; j<cs; ++j)
1419 for (
int d1=0; d1<
dim; ++d1)
1421 for (
int d2=0; d2<
dim; ++d2)
1423 Mkj += ip.
weight*fine_shape(k,d1)*J(d1,d2)*coarse_shape(j,d2);
1426 fine_coarse_mass(k,j) += (fabs(Mkj) < 1e-12) ? 0.0 : Mkj;
1431 fine_mass_inv.
Mult(fine_coarse_mass, I);
1443 #ifdef MFEM_THREAD_SAFE
1453 for (
int k = 0; k <
dof; k++)
1464 for (
int i = 0; i <
dim; i++)
1466 Ikj +=
vshape(j, i) * vk[i];
1468 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1481 #ifdef MFEM_THREAD_SAFE
1487 const double weight = Trans.
Weight();
1488 for (
int j = 0; j <
dof; j++)
1497 for (
int k = 0; k <
dof; k++)
1500 for (
int d = 0; d <
dim; d++)
1502 R_jk +=
vshape(k,d)*pt_data[d];
1524 #ifdef MFEM_THREAD_SAFE
1530 for (
int j = 0; j <
dof; j++)
1537 Jinv.
Mult(tk+
dim*d2t[j], pt_data);
1538 for (
int k = 0; k <
dof; k++)
1541 for (
int d = 0; d <
dim; d++)
1543 R_jk +=
vshape(k,d)*pt_data[d];
1559 : etype(etype), auxiliary_basis(NULL), scale_integrated(false)
1568 for (
int i = 0; i <=
p; i++)
1582 for (
int i = 0; i <=
p; i++)
1584 for (
int j = 0; j < i; j++)
1586 double xij = x(i) - x(j);
1591 for (
int i = 0; i <=
p; i++)
1598 for (
int i = 0; i <
p; i++)
1602 mfem_error(
"Poly_1D::Basis::Basis : nodes are not increasing!");
1612 auxiliary_basis =
new Basis(
1634 int i, k,
p = x.Size() - 1;
1644 for (k = 0; k <
p; k++)
1646 if (y >= (x(k) + x(k+1))/2)
1652 for (i = k+1; i <=
p; i++)
1659 l = lk * (y - x(k));
1661 for (i = 0; i < k; i++)
1663 u(i) = l * w(i) / (y - x(i));
1666 for (i++; i <=
p; i++)
1668 u(i) = l * w(i) / (y - x(i));
1676 auxiliary_basis->Eval(y, u_aux, d_aux);
1677 EvalIntegrated(d_aux, u);
1696 int i, k,
p = x.Size() - 1;
1697 double l, lp, lk, sk, si;
1707 for (k = 0; k <
p; k++)
1709 if (y >= (x(k) + x(k+1))/2)
1715 for (i = k+1; i <=
p; i++)
1722 l = lk * (y - x(k));
1725 for (i = 0; i < k; i++)
1727 si = 1.0/(y - x(i));
1729 u(i) = l * si * w(i);
1732 for (i++; i <=
p; i++)
1734 si = 1.0/(y - x(i));
1736 u(i) = l * si * w(i);
1740 for (i = 0; i < k; i++)
1742 d(i) = (lp * w(i) -
u(i))/(y - x(i));
1745 for (i++; i <=
p; i++)
1747 d(i) = (lp * w(i) -
u(i))/(y - x(i));
1755 auxiliary_basis->Eval(y, u_aux, d_aux, d2_aux);
1756 EvalIntegrated(d_aux,u);
1757 EvalIntegrated(d2_aux,d);
1767 "Basis::Eval with second order derivatives not implemented for"
1768 " etype = " << etype);
1781 int i, k,
p = x.Size() - 1;
1782 double l, lp, lp2, lk, sk, si, sk2;
1793 for (k = 0; k <
p; k++)
1795 if (y >= (x(k) + x(k+1))/2)
1801 for (i = k+1; i <=
p; i++)
1808 l = lk * (y - x(k));
1812 for (i = 0; i < k; i++)
1814 si = 1.0/(y - x(i));
1817 u(i) = l * si * w(i);
1820 for (i++; i <=
p; i++)
1822 si = 1.0/(y - x(i));
1825 u(i) = l * si * w(i);
1828 lp2 = lp * sk + l * sk2 + sk * lk;
1830 for (i = 0; i < k; i++)
1832 d(i) = (lp * w(i) -
u(i))/(y - x(i));
1833 d2(i) = (lp2 * w(i) - 2 * d(i))/(y - x(i));
1836 d2(k) = sk2 *
u(k) + sk * d(k);
1837 for (i++; i <=
p; i++)
1839 d(i) = (lp * w(i) -
u(i))/(y - x(i));
1840 d2(i) = (lp2 * w(i) - 2 * d(i))/(y - x(i));
1848 MFEM_ABORT(
"Integrated basis must be evaluated with EvalIntegrated");
1857 "EvalIntegrated is only valid for Integrated basis type");
1858 int p = d_aux_.
Size() - 1;
1862 for (
int j=1; j<
p; ++j)
1864 u[j] = u[j-1] - d_aux_[j];
1869 if (scale_integrated)
1871 Vector &aux_nodes = auxiliary_basis->x;
1872 for (
int j=0; j<aux_nodes.
Size()-1; ++j)
1874 u[j] *= aux_nodes[j+1] - aux_nodes[j];
1881 scale_integrated = scale_integrated_;
1886 delete auxiliary_basis;
1894 for (
int i = 0; i <=
p; i++)
1896 binom(i,0) = binom(i,i) = 1;
1897 for (
int j = 1; j < i; j++)
1899 binom(i,j) = binom(i-1,j) + binom(i-1,j-1);
1908 for (
int i = 0; i <=
p; i++)
1911 double s = sin(M_PI_2*(i + 0.5)/(p + 1));
1916 void Poly_1D::CalcMono(
const int p,
const double x,
double *
u)
1920 for (
int n = 1; n <=
p; n++)
1926 void Poly_1D::CalcMono(
const int p,
const double x,
double *u,
double *d)
1931 for (
int n = 1; n <=
p; n++)
1951 for (i = 1; i <
p; i++)
1958 for (i--; i > 0; i--)
1968 double *u,
double *d)
1979 const double xpy = x + y, ptx = p*x;
1982 for (i = 1; i <
p; i++)
1984 d[i] = b[i]*z*(i*xpy - ptx);
1991 for (i--; i > 0; i--)
2013 const double xpy = x + y, ptx = p*x;
2016 for (i = 1; i <
p; i++)
2018 d[i] = b[i]*z*(i*xpy - ptx);
2023 for (i--; i > 0; i--)
2038 if (p == 0) {
return; }
2039 u[1] = z = 2.*x - 1.;
2040 for (
int n = 1; n <
p; n++)
2042 u[n+1] = ((2*n + 1)*z*u[n] - n*u[n-1])/(n + 1);
2055 if (p == 0) {
return; }
2056 u[1] = z = 2.*x - 1.;
2058 for (
int n = 1; n <
p; n++)
2060 u[n+1] = ((2*n + 1)*z*u[n] - n*u[n-1])/(n + 1);
2061 d[n+1] = (4*n + 2)*u[n] + d[n-1];
2065 void Poly_1D::CalcChebyshev(
const int p,
const double x,
double *u)
2072 if (p == 0) {
return; }
2073 u[1] = z = 2.*x - 1.;
2074 for (
int n = 1; n <
p; n++)
2076 u[n+1] = 2*z*u[n] - u[n-1];
2080 void Poly_1D::CalcChebyshev(
const int p,
const double x,
double *u,
double *d)
2093 if (p == 0) {
return; }
2094 u[1] = z = 2.*x - 1.;
2096 for (
int n = 1; n <
p; n++)
2098 u[n+1] = 2*z*u[n] - u[n-1];
2099 d[n+1] = (n + 1)*(z*d[n]/n + 2*u[n]);
2103 void Poly_1D::CalcChebyshev(
const int p,
const double x,
double *u,
double *d,
2119 if (p == 0) {
return; }
2120 u[1] = z = 2.*x - 1.;
2123 for (
int n = 1; n <
p; n++)
2125 u[n+1] = 2*z*u[n] - u[n-1];
2126 d[n+1] = (n + 1)*(z*d[n]/n + 2*u[n]);
2127 dd[n+1] = (n + 1)*(2.*(n + 1)*d[n] + z*dd[n])/n;
2138 if (points_container.find(btype) == points_container.end())
2143 if (pts.
Size() <=
p)
2149 pts[
p] =
new double[p + 1];
2159 if ( bases_container.find(btype) == bases_container.end() )
2165 if (bases.
Size() <=
p)
2169 if (bases[p] == NULL)
2182 for (PointsMap::iterator it = points_container.begin();
2183 it != points_container.end() ; ++it)
2186 for (
int i = 0 ; i < pts.
Size() ; ++i )
2193 for (BasisMap::iterator it = bases_container.begin();
2194 it != bases_container.end() ; ++it)
2197 for (
int i = 0 ; i < bases.
Size() ; ++i )
2220 for (
int i = 1; i <
p; i++)
2228 const int p1 = p + 1;
2239 for (
int i = 1; i <
p; i++)
2243 for (
int i = 1; i <
p; i++)
2247 for (
int i = 1; i <
p; i++)
2251 for (
int i = 1; i <
p; i++)
2257 for (
int j = 1; j <
p; j++)
2259 for (
int i = 1; i <
p; i++)
2268 const int p1 = p + 1;
2272 dof_map[0 + (0 + 0*p1)*p1] = 0;
2273 dof_map[p + (0 + 0*p1)*p1] = 1;
2274 dof_map[p + (p + 0*p1)*p1] = 2;
2275 dof_map[0 + (p + 0*p1)*p1] = 3;
2276 dof_map[0 + (0 + p*p1)*p1] = 4;
2277 dof_map[p + (0 + p*p1)*p1] = 5;
2278 dof_map[p + (p + p*p1)*p1] = 6;
2279 dof_map[0 + (p + p*p1)*p1] = 7;
2284 for (
int i = 1; i <
p; i++)
2286 dof_map[i + (0 + 0*p1)*p1] = o++;
2288 for (
int i = 1; i <
p; i++)
2290 dof_map[p + (i + 0*p1)*p1] = o++;
2292 for (
int i = 1; i <
p; i++)
2294 dof_map[i + (p + 0*p1)*p1] = o++;
2296 for (
int i = 1; i <
p; i++)
2298 dof_map[0 + (i + 0*p1)*p1] = o++;
2300 for (
int i = 1; i <
p; i++)
2302 dof_map[i + (0 + p*p1)*p1] = o++;
2304 for (
int i = 1; i <
p; i++)
2306 dof_map[p + (i + p*p1)*p1] = o++;
2308 for (
int i = 1; i <
p; i++)
2310 dof_map[i + (p + p*p1)*p1] = o++;
2312 for (
int i = 1; i <
p; i++)
2314 dof_map[0 + (i + p*p1)*p1] = o++;
2316 for (
int i = 1; i <
p; i++)
2318 dof_map[0 + (0 + i*p1)*p1] = o++;
2320 for (
int i = 1; i <
p; i++)
2322 dof_map[p + (0 + i*p1)*p1] = o++;
2324 for (
int i = 1; i <
p; i++)
2326 dof_map[p + (p + i*p1)*p1] = o++;
2328 for (
int i = 1; i <
p; i++)
2330 dof_map[0 + (p + i*p1)*p1] = o++;
2334 for (
int j = 1; j <
p; j++)
2336 for (
int i = 1; i <
p; i++)
2338 dof_map[i + ((p-j) + 0*p1)*p1] = o++;
2341 for (
int j = 1; j <
p; j++)
2343 for (
int i = 1; i <
p; i++)
2345 dof_map[i + (0 + j*p1)*p1] = o++;
2348 for (
int j = 1; j <
p; j++)
2350 for (
int i = 1; i <
p; i++)
2352 dof_map[p + (i + j*p1)*p1] = o++;
2355 for (
int j = 1; j <
p; j++)
2357 for (
int i = 1; i <
p; i++)
2359 dof_map[(p-i) + (p + j*p1)*p1] = o++;
2362 for (
int j = 1; j <
p; j++)
2364 for (
int i = 1; i <
p; i++)
2366 dof_map[0 + ((p-i) + j*p1)*p1] = o++;
2369 for (
int j = 1; j <
p; j++)
2371 for (
int i = 1; i <
p; i++)
2373 dof_map[i + (j + p*p1)*p1] = o++;
2378 for (
int k = 1; k <
p; k++)
2380 for (
int j = 1; j <
p; j++)
2382 for (
int i = 1; i <
p; i++)
2384 dof_map[i + (j + k*p1)*p1] = o++;
2391 MFEM_ABORT(
"invalid dimension: " << dims);
2402 MFEM_ABORT(
"invalid DofMapType: " << dmtype);
2441 cbasis1d(
poly1d.GetBasis(p, VerifyClosed(cbtype))),
2442 obasis1d(
poly1d.GetBasis(p - 1, VerifyOpen(obtype)))
2444 MFEM_VERIFY(dims > 1,
"Constructor for VectorTensorFiniteElement with both "
2445 "open and closed bases is not valid for 1D elements.");
2457 cbasis1d(
poly1d.GetBasis(p, VerifyOpen(obtype))),
2458 obasis1d(
poly1d.GetBasis(p, VerifyOpen(obtype)))
2460 MFEM_VERIFY(dims == 1,
"Constructor for VectorTensorFiniteElement without "
2461 "closed basis is only valid for 1D elements.");
2485 const bool closed)
const
2490 i < (closed ?
dof2quad_array.Size() : dof2quad_array_open.Size());
2494 if (d2q.
IntRule == &ir && d2q.
mode == mode) {
return d2q; }
2499 const int nqpt = (int)floor(pow(ir.
GetNPoints(), 1.0/
dim) + 0.5);
2509 Vector val(ndof), grad(ndof);
2510 for (
int i = 0; i < nqpt; i++)
2524 for (
int j = 0; j < ndof; j++)
2526 d2q->
B[i+nqpt*j] = d2q->
Bt[j+ndof*i] = val(j);
2527 d2q->
G[i+nqpt*j] = d2q->
Gt[j+ndof*i] = grad(j);
2537 dof2quad_array_open.Append(d2q);
2545 for (
int i = 0; i < dof2quad_array_open.Size(); i++)
2547 delete dof2quad_array_open[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)
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
int Size() const
Return the logical size of the array.
void Get(double *p, const int dim) const
int GetDim() const
Returns the reference space dimension for the finite element.
void trans(const Vector &u, Vector &x)
Array< int > lex_ordering
Class for an integration rule - an Array of IntegrationPoint.
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
virtual void CalcPhysLinLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
void ProjectMatrixCoefficient_RT(const double *nk, const Array< int > &d2n, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an RT space.
EvalType
One-dimensional basis evaluation type.
void ProjectMatrixCoefficient_ND(const double *tk, const Array< int > &d2t, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an ND 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...
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
NodalTensorFiniteElement(const int dims, const int p, const int btype, const DofMapType dmtype)
void LocalInterpolation_RT(const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Basis(const int p, const double *nodes, EvalType etype=Barycentric)
Create a nodal or positive (Bernstein) basis of degree p.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Base class for vector Coefficients that optionally depend on time and space.
void SetRow(int r, const double *row)
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
virtual void CalcPhysHessian(ElementTransformation &Trans, DenseMatrix &Hessian) const
Evaluate the Hessian of all shape functions of a scalar finite element in reference space at the give...
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in physical space at the po...
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 ...
static void GivePolyPoints(const int np, double *pts, const int type)
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
void SetSize(int s)
Resize the vector to size s.
virtual void CalcPhysLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the gi...
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
FiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct FiniteElement with given.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
int dim
Dimension of reference space.
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
TensorBasisElement(const int dims, const int p, const int btype, const DofMapType dmtype)
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
Array< double > Gt
Transpose of G.
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
~VectorTensorFiniteElement()
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
VectorFiniteElement(int D, Geometry::Type G, int Do, int O, int M, int F=FunctionSpace::Pk)
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
virtual void SetMapType(const int map_type_)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
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.
static void CalcLegendre(const int p, const double x, double *u)
void NodalLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get the matrix I that defines nodal interpolation between this element and the refined element fine_f...
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Use change of basis, O(p^2) Evals.
const DofToQuad & GetTensorDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode, const bool closed) const
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
int cdim
Dimension of curl for vector-valued basis functions.
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 ProjectFromNodes(Vector &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector of values at the finite element nodes and a transformation, compute its projection (ap...
int GetHeight() const
Get the height of the matrix.
static void CalcBinomTerms(const int p, const double x, const double y, double *u)
Compute the p terms in the expansion of the binomial (x + y)^p and store them in the already allocate...
const double * GetPoints(const int p, const int btype)
Get the coordinates of the points of the given BasisType, btype.
Poly_1D::Basis & obasis1d
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.
Class for standard nodal finite elements.
void ScalarLocalRestriction(ElementTransformation &Trans, DenseMatrix &R, const ScalarFiniteElement &coarse_fe) const
Get restriction matrix R defined through local L2-projection in the space defined by the coarse_fe...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
VectorTensorFiniteElement(const int dims, const int d, const int p, const int cbtype, const int obtype, const int M, const DofMapType dmtype)
Class for finite elements with basis functions that return scalar values.
void SetSize(int m, int n)
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...
Poly_1D::Basis & cbasis1d
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Fast evaluation of Bernstein polynomials.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
static void ChebyshevPoints(const int p, double *x)
Compute the points for the Chebyshev polynomials of order p and place them in the already allocated x...
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
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...
void Invert()
Replaces the current matrix with its inverse.
class FiniteElement * FE
The FiniteElement that created and owns this object.
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
int orders[Geometry::MaxDim]
Anisotropic orders.
virtual ~FiniteElement()
Deconstruct the FiniteElement.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Use barycentric Lagrangian interpolation, O(p) Evals.
int GetVDim()
Returns dimension of the vector.
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
void GetColumn(int c, Vector &col) const
void ScalarLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get matrix I "Interpolation" defined through local L2-projection in the space defined by the fine_fe...
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
void Threshold(double eps)
Replace small entries, abs(a_ij) <= eps, with zero.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
void Set3(const double x1, const double x2, const double x3)
Array< double > Bt
Transpose of B.
double * Data() const
Returns the matrix data array.
void LocalInterpolation_ND(const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
double p(const Vector &x, double t)
Integrated indicator functions (cf. Gerritsma)
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &Hessian) const
Evaluate the Hessians of all shape functions of a scalar finite element in reference space at the giv...
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Base class for Matrix Coefficients that optionally depend on time and space.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
static const int * Binom(const int p)
Get a pointer to an array containing the binomial coefficients "p choose k" for k=0,...,p for the given 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 Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Integrated GLL indicator functions.
Array< double > B
Basis functions evaluated at quadrature points.
virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const
Get the dofs associated with the given face. *dofs is set to an internal array of the local dofc on t...
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
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...
void LocalRestriction_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
Class for integration point with weight.
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Full multidimensional representation which does not use tensor product structure. The ordering of the...
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
const Poly_1D::Basis & GetBasis1D() const
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
void Project_RT(const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the RT basis functions.
int dof
Number of degrees of freedom.
Array< DofToQuad * > dof2quad_array
Container for all DofToQuad objects created by the FiniteElement.
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto 'this' FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Implements CalcDivShape methods.
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
void ProjectCurl_2D(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
void LocalL2Projection_ND(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
void Mult(const double *x, double *y) const
Matrix vector multiplication.
static void CalcBernstein(const int p, const double x, double *u)
Compute the values of the Bernstein basis functions of order p at coordinate x and store the results ...
int GetWidth() const
Get the width of the matrix.
void LocalRestriction_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
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.
void pts(int iphi, int t, double x[])
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
double u(const Vector &xvec)
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
void Project_ND(const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the ND basis functions.
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
double FNorm2() const
Compute the square of the Frobenius norm of the matrix.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
static void CalcDBinomTerms(const int p, const double x, const double y, double *d)
Compute the derivatives (w.r.t. x) of the terms in the expansion of the binomial (x + y)^p assuming t...
const IntegrationRule * IntRule
IntegrationRule that defines the quadrature points at which the basis functions of the FE are evaluat...
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
const DofToQuad & GetTensorDofToQuad(const class TensorBasisElement &tb, const IntegrationRule &ir, DofToQuad::Mode mode) const
int order
Order/degree of the shape functions.
void LocalL2Projection_RT(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
No derivatives implemented.
Implements CalcCurlShape methods.
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const