26 MFEM_ABORT(
"BilinearFormIntegrator::AssemblePA(fes)\n"
27 " is not implemented for this class.");
32 mfem_error (
"BilinearFormIntegrator::AssembleNURBSPA(fes)\n"
33 " is not implemented for this class.");
39 MFEM_ABORT(
"BilinearFormIntegrator::AssemblePA(fes, fes)\n"
40 " is not implemented for this class.");
45 MFEM_ABORT(
"BilinearFormIntegrator::AssemblePABoundary(fes)\n"
46 " is not implemented for this class.");
51 MFEM_ABORT(
"BilinearFormIntegrator::AssemblePAInteriorFaces(...)\n"
52 " is not implemented for this class.");
57 MFEM_ABORT(
"BilinearFormIntegrator::AssemblePABoundaryFaces(...)\n"
58 " is not implemented for this class.");
63 MFEM_ABORT(
"BilinearFormIntegrator::AssembleDiagonalPA(...)\n"
64 " is not implemented for this class.");
71 MFEM_ABORT(
"BilinearFormIntegrator::AssembleEA(...)\n"
72 " is not implemented for this class.");
81 MFEM_ABORT(
"BilinearFormIntegrator::AssembleEAInteriorFaces(...)\n"
82 " is not implemented for this class.");
90 MFEM_ABORT(
"BilinearFormIntegrator::AssembleEABoundaryFaces(...)\n"
91 " is not implemented for this class.");
96 MFEM_ABORT(
"BilinearFormIntegrator::AssembleDiagonalPA_ADAt(...)\n"
97 " is not implemented for this class.");
102 MFEM_ABORT(
"BilinearFormIntegrator:AddMultPA:(...)\n"
103 " is not implemented for this class.");
108 MFEM_ABORT(
"BilinearFormIntegrator::AddMultNURBSPA(...)\n"
109 " is not implemented for this class.");
114 MFEM_ABORT(
"BilinearFormIntegrator::AddMultTransposePA(...)\n"
115 " is not implemented for this class.");
120 MFEM_ABORT(
"BilinearFormIntegrator::AssembleMF(...)\n"
121 " is not implemented for this class.");
126 MFEM_ABORT(
"BilinearFormIntegrator::AddMultMF(...)\n"
127 " is not implemented for this class.");
132 MFEM_ABORT(
"BilinearFormIntegrator::AddMultTransposeMF(...)\n"
133 " is not implemented for this class.");
138 MFEM_ABORT(
"BilinearFormIntegrator::AssembleDiagonalMF(...)\n"
139 " is not implemented for this class.");
146 MFEM_ABORT(
"BilinearFormIntegrator::AssembleElementMatrix(...)\n"
147 " is not implemented for this class.");
154 MFEM_ABORT(
"BilinearFormIntegrator::AssembleElementMatrix2(...)\n"
155 " is not implemented for this class.");
161 mfem_error (
"BilinearFormIntegrator::AssemblePatchMatrix(...)\n"
162 " is not implemented for this class.");
169 MFEM_ABORT(
"BilinearFormIntegrator::AssembleFaceMatrix(...)\n"
170 " is not implemented for this class.");
178 MFEM_ABORT(
"AssembleFaceMatrix (mixed form) is not implemented for this"
179 " Integrator class.");
188 MFEM_ABORT(
"AssembleTraceFaceMatrix (DPG form) is not implemented for this"
189 " Integrator class.");
195 MFEM_ABORT(
"Not implemented.");
206 elmat.
Mult(elfun, elvect);
217 elmat.
Mult(elfun, elvect);
281 for (
int i = 0; i < integrators.Size(); i++)
283 integrators[i]->SetIntRule(ir);
290 MFEM_ASSERT(integrators.Size() > 0,
"empty SumIntegrator.");
292 integrators[0]->AssembleElementMatrix(el, Trans, elmat);
293 for (
int i = 1; i < integrators.Size(); i++)
295 integrators[i]->AssembleElementMatrix(el, Trans, elem_mat);
304 MFEM_ASSERT(integrators.Size() > 0,
"empty SumIntegrator.");
306 integrators[0]->AssembleElementMatrix2(el1, el2, Trans, elmat);
307 for (
int i = 1; i < integrators.Size(); i++)
309 integrators[i]->AssembleElementMatrix2(el1, el2, Trans, elem_mat);
318 MFEM_ASSERT(integrators.Size() > 0,
"empty SumIntegrator.");
320 integrators[0]->AssembleFaceMatrix(el1, el2, Trans, elmat);
321 for (
int i = 1; i < integrators.Size(); i++)
323 integrators[i]->AssembleFaceMatrix(el1, el2, Trans, elem_mat);
333 MFEM_ASSERT(integrators.Size() > 0,
"empty SumIntegrator.");
335 integrators[0]->AssembleFaceMatrix(tr_fe, te_fe1, te_fe2, Trans, elmat);
336 for (
int i = 1; i < integrators.Size(); i++)
338 integrators[i]->AssembleFaceMatrix(tr_fe, te_fe1, te_fe2, Trans, elem_mat);
345 for (
int i = 0; i < integrators.Size(); i++)
347 integrators[i]->AssemblePA(fes);
353 for (
int i = 0; i < integrators.Size(); i++)
355 integrators[i]->AssembleDiagonalPA(diag);
361 for (
int i = 0; i < integrators.Size(); i++)
363 integrators[i]->AssemblePAInteriorFaces(fes);
369 for (
int i = 0; i < integrators.Size(); i++)
371 integrators[i]->AssemblePABoundaryFaces(fes);
377 for (
int i = 0; i < integrators.Size(); i++)
379 integrators[i]->AddMultPA(x, y);
385 for (
int i = 0; i < integrators.Size(); i++)
387 integrators[i]->AddMultTransposePA(x, y);
393 for (
int i = 0; i < integrators.Size(); i++)
395 integrators[i]->AssembleMF(fes);
401 for (
int i = 0; i < integrators.Size(); i++)
403 integrators[i]->AddMultTransposeMF(x, y);
409 for (
int i = 0; i < integrators.Size(); i++)
411 integrators[i]->AddMultMF(x, y);
417 for (
int i = 0; i < integrators.Size(); i++)
419 integrators[i]->AssembleDiagonalMF(diag);
426 for (
int i = 0; i < integrators.Size(); i++)
428 integrators[i]->AssembleEA(fes, emat,
add);
437 for (
int i = 0; i < integrators.Size(); i++)
439 integrators[i]->AssembleEAInteriorFaces(fes,ea_data_int,ea_data_ext,
add);
447 for (
int i = 0; i < integrators.Size(); i++)
449 integrators[i]->AssembleEABoundaryFaces(fes, ea_data_bdr,
add);
457 for (
int i = 0; i < integrators.Size(); i++)
459 delete integrators[i];
471 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
474#ifdef MFEM_THREAD_SAFE
475 Vector test_shape(test_nd);
489 elmat.
SetSize(test_nd, trial_nd);
511 w *=
Q->
Eval(Trans, ip);
515#ifndef MFEM_THREAD_SAFE
531 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
539 "Dimension mismatch in height of matrix coefficient.");
541 "Dimension mismatch in width of matrix coefficient.");
545 MFEM_VERIFY(trial_vdim == test_vdim,
546 "Diagonal matrix coefficient requires matching "
547 "test and trial vector dimensions.");
549 "Dimension mismatch in diagonal matrix coefficient.");
553 MFEM_VERIFY(
VQ->
GetVDim() == 3,
"Vector coefficient must have "
554 "dimension equal to three.");
557#ifdef MFEM_THREAD_SAFE
568 test_shape.SetSize(test_nd, test_vdim);
569 shape_tmp.
SetSize(test_nd, trial_vdim);
573 trial_shape.
Reset(test_shape.Data(), trial_nd, trial_vdim);
577 trial_shape.
SetSize(trial_nd, trial_vdim);
580 elmat.
SetSize(test_nd, trial_nd);
607 Mult(test_shape, M, shape_tmp);
621 for (
int j=0; j<test_nd; j++)
627 if (test_vdim == 3 && trial_vdim == 3)
629 shape_tmp(j,0) = test_shape(j,1) * V(2) -
630 test_shape(j,2) * V(1);
631 shape_tmp(j,1) = test_shape(j,2) * V(0) -
632 test_shape(j,0) * V(2);
633 shape_tmp(j,2) = test_shape(j,0) * V(1) -
634 test_shape(j,1) * V(0);
636 else if (test_vdim == 3 && trial_vdim == 2)
638 shape_tmp(j,0) = test_shape(j,1) * V(2) -
639 test_shape(j,2) * V(1);
640 shape_tmp(j,1) = test_shape(j,2) * V(0) -
641 test_shape(j,0) * V(2);
643 else if (test_vdim == 3 && trial_vdim == 1)
645 shape_tmp(j,0) = test_shape(j,1) * V(2) -
646 test_shape(j,2) * V(1);
648 else if (test_vdim == 2 && trial_vdim == 3)
650 shape_tmp(j,0) = test_shape(j,1) * V(2);
651 shape_tmp(j,1) = -test_shape(j,0) * V(2);
652 shape_tmp(j,2) = test_shape(j,0) * V(1) -
653 test_shape(j,1) * V(0);
655 else if (test_vdim == 2 && trial_vdim == 2)
657 shape_tmp(j,0) = test_shape(j,1) * V(2);
658 shape_tmp(j,1) = -test_shape(j,0) * V(2);
660 else if (test_vdim == 1 && trial_vdim == 3)
662 shape_tmp(j,0) = 0.0;
663 shape_tmp(j,1) = -test_shape(j,0) * V(2);
664 shape_tmp(j,2) = test_shape(j,0) * V(1);
666 else if (test_vdim == 1 && trial_vdim == 1)
668 shape_tmp(j,0) = 0.0;
677 w *=
Q -> Eval (Trans, ip);
689#ifndef MFEM_THREAD_SAFE
704 MFEM_VERIFY(
VQ,
"MixedScalarVectorIntegrator: "
705 "VectorCoefficient must be set");
711 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
712 int sca_nd = sca_fe->
GetDof();
713 int vec_nd = vec_fe->
GetDof();
717 MFEM_VERIFY(
VQ->
GetVDim() == vdim,
"MixedScalarVectorIntegrator: "
718 "Dimensions of VectorCoefficient and Vector-valued basis "
719 "functions must match");
721#ifdef MFEM_THREAD_SAFE
725 Vector vshape_tmp(vec_nd);
736 elmat.
SetSize(test_nd, trial_nd);
766 vshape.
Mult(V,vshape_tmp);
777 int trial_dof = trial_fe.
GetDof();
778 int test_dof = test_fe.
GetDof();
782 dshape.
SetSize(trial_dof, dim);
783 gshape.
SetSize(trial_dof, dim);
786 elmat.
SetSize(dim * test_dof, trial_dof);
792 elmat_comp.
SetSize(test_dof, trial_dof);
804 Mult(dshape, Jadj, gshape);
809 c *=
Q->
Eval(Trans, ip);
813 for (
int d = 0; d < dim; ++d)
816 MultVWt(shape, d_col, elmat_comp);
817 for (
int jj = 0; jj < trial_dof; ++jj)
819 for (
int ii = 0; ii < test_dof; ++ii)
821 elmat(d * test_dof + ii, jj) += elmat_comp(ii, jj);
845 bool square = (dim == spaceDim);
851 "Unexpected dimension for VectorCoefficient");
856 "Unexpected width for MatrixCoefficient");
858 "Unexpected height for MatrixCoefficient");
861#ifdef MFEM_THREAD_SAFE
862 DenseMatrix dshape(nd, dim), dshapedxt(nd, spaceDim);
868 dshapedxt.
SetSize(nd, spaceDim);
869 dshapedxt_m.
SetSize(nd,
MQ ? spaceDim : 0);
870 M.SetSize(
MQ ? spaceDim : 0);
880 bool deleteRule =
false;
883 const int patch = NURBSFE->
GetPatch();
884 const int* ijk = NURBSFE->
GetIJK();
898 w = ip.
weight / (square ? w : w*w*w);
906 Mult(dshapedxt, M, dshapedxt_m);
919 w *=
Q->
Eval(Trans, ip);
935 int tr_nd = trial_fe.
GetDof();
936 int te_nd = test_fe.
GetDof();
939 bool square = (dim == spaceDim);
945 "Unexpected dimension for VectorCoefficient");
950 "Unexpected width for MatrixCoefficient");
952 "Unexpected height for MatrixCoefficient");
955#ifdef MFEM_THREAD_SAFE
956 DenseMatrix dshape(tr_nd, dim), dshapedxt(tr_nd, spaceDim);
957 DenseMatrix te_dshape(te_nd, dim), te_dshapedxt(te_nd, spaceDim);
964 dshapedxt.
SetSize(tr_nd, spaceDim);
966 te_dshapedxt.
SetSize(te_nd, spaceDim);
967 invdfdx.
SetSize(dim, spaceDim);
968 dshapedxt_m.
SetSize(te_nd,
MQ ? spaceDim : 0);
969 M.SetSize(
MQ ? spaceDim : 0);
986 w = ip.
weight / (square ? w : w*w*w);
987 Mult(dshape, invdfdx, dshapedxt);
988 Mult(te_dshape, invdfdx, te_dshapedxt);
994 Mult(te_dshapedxt, M, dshapedxt_m);
1007 w *=
Q->
Eval(Trans, ip);
1027 "Unexpected dimension for VectorCoefficient");
1032 "Unexpected width for MatrixCoefficient");
1034 "Unexpected height for MatrixCoefficient");
1037#ifdef MFEM_THREAD_SAFE
1038 DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim), M(
MQ ? spaceDim : 0);
1042 invdfdx.
SetSize(dim, spaceDim);
1070 w *=
Q->
Eval(Tr, ip);
1080 M.
Mult(vecdxt, pointflux);
1085 for (
int j=0; j<spaceDim; ++j)
1087 pointflux[j] = D[j] * vecdxt[j];
1092 invdfdx.
Mult(pointflux, vec);
1102 int nd, spaceDim, fnd;
1111 "Unexpected dimension for VectorCoefficient");
1116 "Unexpected width for MatrixCoefficient");
1118 "Unexpected height for MatrixCoefficient");
1121#ifdef MFEM_THREAD_SAFE
1122 DenseMatrix dshape(nd,dim), invdfdx(dim, spaceDim);
1127 invdfdx.
SetSize(dim, spaceDim);
1140 flux.
SetSize( fnd * spaceDim );
1142 for (
int i = 0; i < fnd; i++)
1158 vecdxt *=
Q->
Eval(Trans,ip);
1160 for (
int j = 0; j < spaceDim; j++)
1162 flux(fnd*j+i) = vecdxt(j);
1170 M.
Mult(vecdxt, pointflux);
1175 for (
int j=0; j<spaceDim; ++j)
1177 pointflux[j] = D[j] * vecdxt[j];
1180 for (
int j = 0; j < spaceDim; j++)
1182 flux(fnd*j+i) = pointflux(j);
1188 for (
int j = 0; j < spaceDim; j++)
1190 flux(fnd*j+i) = vecdxt(j);
1200 int nd = fluxelem.
GetDof();
1204#ifdef MFEM_THREAD_SAFE
1213 if (d_energy) { vec.
SetSize(spaceDim); }
1216 int order = 2 * fluxelem.
GetOrder();
1220 if (d_energy) { *d_energy = 0.0; }
1228 for (
int k = 0; k < spaceDim; k++)
1230 for (
int j = 0; j < nd; j++)
1232 pointflux(k) += flux(k*nd+j)*shape(j);
1248 energy += w * (D * pointflux);
1252 real_t e = (pointflux * pointflux);
1253 if (
Q) { e *=
Q->
Eval(Trans, ip); }
1261 for (
int k = 0; k < dim; k++)
1263 (*d_energy)[k] += w * vec[k] * vec[k];
1302#ifdef MFEM_THREAD_SAFE
1321 w *=
Q -> Eval(Trans, ip);
1332 int tr_nd = trial_fe.
GetDof();
1333 int te_nd = test_fe.
GetDof();
1336#ifdef MFEM_THREAD_SAFE
1344 &
GetRule(trial_fe, test_fe, Trans);
1358 w *=
Q -> Eval(Trans, ip);
1385 MFEM_ASSERT(Trans.
Elem2No < 0,
1386 "support for interior faces is not implemented");
1391#ifdef MFEM_THREAD_SAFE
1420 w *=
Q -> Eval(Trans, ip);
1433#ifdef MFEM_THREAD_SAFE
1435 Vector shape, vec2, BdFidxT;
1453 Q->
Eval(Q_ir, Trans, *ir);
1467 adjJ.
Mult(vec1, vec2);
1468 dshape.
Mult(vec2, BdFidxT);
1506 Mult(dshape, adjJ, grad);
1511 for (
int k = 0; k < nd; k++)
1514 for (
int l = 0; l < nd; l++)
1517 for (
int s = 0;
s <
dim;
s++)
1519 a += Q_nodal(
s,k)*grad(l,
s);
1521 elmat(k,l) += wsk*
a;
1552 vdim = (vdim == -1) ? spaceDim : vdim;
1594 VQ->
Eval(vec, Trans, ip);
1595 for (
int k = 0; k < vdim; k++)
1597 elmat.
AddMatrix(norm*vec(k), partelmat, nd*k, nd*k);
1602 MQ->
Eval(mcoeff, Trans, ip);
1603 for (
int i = 0; i < vdim; i++)
1604 for (
int j = 0; j < vdim; j++)
1606 elmat.
AddMatrix(norm*mcoeff(i,j), partelmat, nd*i, nd*j);
1613 norm *=
Q->
Eval(Trans, ip);
1616 for (
int k = 0; k < vdim; k++)
1628 int tr_nd = trial_fe.
GetDof();
1629 int te_nd = test_fe.
GetDof();
1636 elmat.
SetSize(te_nd*vdim, tr_nd*vdim);
1639 partelmat.
SetSize(te_nd, tr_nd);
1653 Trans.
OrderW() + Q_order);
1675 MultVWt(te_shape, shape, partelmat);
1679 VQ->
Eval(vec, Trans, ip);
1680 for (
int k = 0; k < vdim; k++)
1682 elmat.
AddMatrix(norm*vec(k), partelmat, te_nd*k, tr_nd*k);
1687 MQ->
Eval(mcoeff, Trans, ip);
1688 for (
int i = 0; i < vdim; i++)
1689 for (
int j = 0; j < vdim; j++)
1691 elmat.
AddMatrix(norm*mcoeff(i,j), partelmat, te_nd*i, tr_nd*j);
1698 norm *=
Q->
Eval(Trans, ip);
1701 for (
int k = 0; k < vdim; k++)
1703 elmat.
AddMatrix(partelmat, te_nd*k, tr_nd*k);
1713 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
1715#ifdef MFEM_THREAD_SAFE
1716 Vector divshape(trial_nd), shape(test_nd);
1722 elmat.
SetSize(test_nd, trial_nd);
1742 w *=
Q->
Eval(Trans, ip);
1753 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
1759 "Trial space must be H(Curl) and test space must be H_1");
1761#ifdef MFEM_THREAD_SAFE
1773 elmat.
SetSize(test_nd, trial_nd);
1815 Mult(dshape, invdfdx, dshapedxt);
1823 w *=
Q->
Eval(Trans, ip);
1835 int trial_nd = trial_fe.
GetDof(), test_nd = test_fe.
GetDof(), i;
1837 int dimc = (
dim == 3) ? 3 : 1;
1841 "At least one of the finite elements must be in H(Curl)");
1843 int curl_nd, vec_nd;
1855#ifdef MFEM_THREAD_SAFE
1866 elmat.
SetSize(test_nd, trial_nd);
1913 w *=
Q->
Eval(Trans, ip);
1919 AddMultABt(vshapeTest, curlshapeTrial_dFT, elmat);
1923 AddMultABt(curlshapeTrial_dFT, vshapeTest, elmat);
1935 int trial_nd = trial_fe.
GetDof();
1936 int test_nd = test_fe.
GetDof();
1942 elmat.
SetSize (test_nd,trial_nd);
1944 dshapedxt.
SetSize(trial_nd, spaceDim);
1982 Mult (dshape, invdfdx, dshapedxt);
1986 for (l = 0; l < trial_nd; l++)
1988 dshapedxi(l) = dshapedxt(l,xi);
2005#ifdef MFEM_THREAD_SAFE
2046 Mult(curlshape_dFt, M, curlshape);
2057 w *=
Q->
Eval(Trans, ip);
2072 int tr_nd = trial_fe.
GetDof();
2073 int te_nd = test_fe.
GetDof();
2078#ifdef MFEM_THREAD_SAFE
2123 Mult(te_curlshape_dFt, M, te_curlshape);
2124 AddMultABt(te_curlshape, curlshape_dFt, elmat);
2130 AddMultADBt(te_curlshape_dFt,D,curlshape_dFt,elmat);
2136 w *=
Q->
Eval(Trans, ip);
2139 AddMultABt(te_curlshape_dFt, curlshape_dFt, elmat);
2144void CurlCurlIntegrator
2149#ifdef MFEM_THREAD_SAFE
2153 MFEM_VERIFY(ir == NULL,
"Integration rule (ir) must be NULL")
2158 projcurl.
Mult(
u, flux);
2167 int nd = fluxelem.
GetDof();
2170#ifdef MFEM_THREAD_SAFE
2177 int order = 2 * fluxelem.
GetOrder();
2181 if (d_energy) { *d_energy = 0.0; }
2200 real_t e = w * (pointflux * pointflux);
2209#if ANISO_EXPERIMENTAL
2229#if ANISO_EXPERIMENTAL
2233 int n = (int) round(pow(ir.
GetNPoints(), 1.0/3.0));
2237 for (
int k = 0; k < n; k++)
2238 for (
int l = 0; l < n; l++)
2239 for (
int m = 0; m < n; m++)
2241 Vector &vec = pfluxes[(k*n + l)*n + m];
2244 tmp = vec; tmp -= pfluxes[(k*n + l)*n + (m-1)];
2245 (*d_energy)[0] += (tmp * tmp);
2249 tmp = vec; tmp -= pfluxes[(k*n + (l-1))*n + m];
2250 (*d_energy)[1] += (tmp * tmp);
2254 tmp = vec; tmp -= pfluxes[((k-1)*n + l)*n + m];
2255 (*d_energy)[2] += (tmp * tmp);
2273 int cld = (
dim*(
dim-1))/2;
2275#ifdef MFEM_THREAD_SAFE
2304 Mult(dshape_hat, Jadj, dshape);
2309 w *=
Q->
Eval(Trans, ip);
2322#ifdef MFEM_THREAD_SAFE
2347 MultAtB(elfun_mat, dshape_hat, grad_hat);
2353 Mult(grad_hat, Jadj, grad);
2357 real_t curl = grad(0,1) - grad(1,0);
2362 real_t curl_x = grad(2,1) - grad(1,2);
2363 real_t curl_y = grad(0,2) - grad(2,0);
2364 real_t curl_z = grad(1,0) - grad(0,1);
2365 w *= curl_x * curl_x + curl_y * curl_y + curl_z * curl_z;
2370 w *=
Q->
Eval(Tr, ip);
2378 return 0.5 * energy;
2386 int trial_dof = trial_fe.
GetDof();
2387 int test_dof = test_fe.
GetDof();
2388 int dimc = (
dim == 3) ? 3 : 1;
2392 "Trial finite element must be either 2D/3D H(Curl) or 2D H1");
2395 "Test finite element must be in H1/L2");
2408 elmat_comp.
SetSize(test_dof, trial_dof);
2441 c *=
Q->
Eval(Trans, ip);
2445 for (
int d = 0; d <
dimc; ++d)
2448 for (
int jj = 0; jj < trial_dof; ++jj)
2450 for (
int ii = 0; ii < test_dof; ++ii)
2452 elmat(d * test_dof + ii, jj) += shape(ii) * curldata[jj];
2471#ifdef MFEM_THREAD_SAFE
2476 trial_vshape.
SetSize(dof, vdim);
2506 Mult(trial_vshape,K,tmp);
2519 w *=
Q -> Eval (Trans, ip);
2535 int vdim = std::max(spaceDim, trial_fe.
GetRangeDim());
2536 int trial_dof = trial_fe.
GetDof();
2537 int test_dof = test_fe.
GetDof();
2540#ifdef MFEM_THREAD_SAFE
2546 trial_vshape.
SetSize(trial_dof, spaceDim);
2552 elmat.
SetSize(vdim*test_dof, trial_dof);
2576 for (
int d = 0; d < vdim; d++)
2578 for (
int j = 0; j < test_dof; j++)
2580 for (
int k = 0; k < trial_dof; k++)
2582 elmat(d * test_dof + j, k) +=
2583 shape(j) * D(d) * trial_vshape(k, d);
2592 for (
int d = 0; d < vdim; d++)
2594 for (
int j = 0; j < test_dof; j++)
2596 for (
int k = 0; k < trial_dof; k++)
2599 for (
int vd = 0; vd < spaceDim; vd++)
2601 Kv += K(d, vd) * trial_vshape(k, vd);
2603 elmat(d * test_dof + j, k) += shape(j) * Kv;
2612 w *=
Q->
Eval(Trans, ip);
2614 for (
int d = 0; d < vdim; d++)
2616 for (
int j = 0; j < test_dof; j++)
2618 for (
int k = 0; k < trial_dof; k++)
2620 elmat(d * test_dof + j, k) +=
2621 w * shape(j) * trial_vshape(k, d);
2633 int trial_vdim = std::max(spaceDim, trial_fe.
GetRangeDim());
2634 int test_vdim = std::max(spaceDim, test_fe.
GetRangeDim());
2635 int trial_dof = trial_fe.
GetDof();
2636 int test_dof = test_fe.
GetDof();
2639#ifdef MFEM_THREAD_SAFE
2645 trial_vshape.
SetSize(trial_dof,trial_vdim);
2646 test_vshape.
SetSize(test_dof,test_vdim);
2652 elmat.
SetSize (test_dof, trial_dof);
2676 Mult(test_vshape,K,tmp);
2689 w *=
Q -> Eval (Trans, ip);
2697 MFEM_ABORT(
"VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
2698 " is not implemented for given trial and test bases.");
2709 int trial_dof = trial_fe.
GetDof();
2710 int test_dof = test_fe.
GetDof();
2713 dshape.
SetSize (trial_dof, dim);
2714 gshape.
SetSize (trial_dof, dim);
2716 divshape.
SetSize (dim*trial_dof);
2719 elmat.
SetSize (test_dof, dim*trial_dof);
2726 for (
int i = 0; i < ir -> GetNPoints(); i++)
2736 Mult (dshape, Jadj, gshape);
2743 c *=
Q -> Eval (Trans, ip);
2770#ifdef MFEM_THREAD_SAFE
2786 for (
int i = 0; i < ir -> GetNPoints(); i++)
2797 c *=
Q -> Eval (Trans, ip);
2811 int tr_nd = trial_fe.
GetDof();
2812 int te_nd = test_fe.
GetDof();
2815#ifdef MFEM_THREAD_SAFE
2817 Vector te_divshape(te_nd);
2827 int order = 2 * max(test_fe.
GetOrder(),
2834 for (
int i = 0; i < ir -> GetNPoints(); i++)
2846 c *=
Q -> Eval (Trans, ip);
2859 const int dof = el.
GetDof();
2864 vdim = (vdim <= 0) ?
sdim : vdim;
2865 const bool square = (
dim ==
sdim);
2890 for (
int i = 0; i < ir -> GetNPoints(); i++)
2898 w = ip.
weight / (square ? w : w*w*w);
2905 VQ->
Eval(vcoeff, Trans, ip);
2906 for (
int k = 0; k < vdim; ++k)
2914 MQ->
Eval(mcoeff, Trans, ip);
2915 for (
int ii = 0; ii < vdim; ++ii)
2917 for (
int jj = 0; jj < vdim; ++jj)
2919 Mult_a_AAt(w*mcoeff(ii,jj), dshapedxt, pelmat);
2920 elmat.
AddMatrix(pelmat, dof*ii, dof*jj);
2926 if (
Q) { w *=
Q->
Eval(Trans, ip); }
2928 for (
int k = 0; k < vdim; ++k)
2940 const int dof = el.
GetDof();
2945 vdim = (vdim <= 0) ?
sdim : vdim;
2946 const bool square = (
dim ==
sdim);
2983 w = ip.
weight / (square ? w : w*w*w);
2989 VQ->
Eval(vcoeff, Tr, ip);
2990 for (
int k = 0; k < vdim; ++k)
2992 pelmat *= w*vcoeff(k);
2995 pelmat.
AddMult(vec_in, vec_out);
3000 MQ->
Eval(mcoeff, Tr, ip);
3001 for (
int ii = 0; ii < vdim; ++ii)
3004 for (
int jj = 0; jj < vdim; ++jj)
3006 pelmat *= w*mcoeff(ii,jj);
3008 pelmat.
Mult(vec_in, vec_out);
3014 if (
Q) { w *=
Q->
Eval(Tr, ip); }
3016 for (
int k = 0; k < vdim; ++k)
3020 pelmat.
AddMult(vec_in, vec_out);
3042#ifdef MFEM_THREAD_SAFE
3063 for (
int i = 0; i < ir -> GetNPoints(); i++)
3093 for (
int d = 0; d <
dim; d++)
3095 for (
int k = 0; k < dof; k++)
3096 for (
int l = 0; l < dof; l++)
3098 elmat (dof*d+k, dof*d+l) += (M * w) * pelmat(k, l);
3101 for (
int ii = 0; ii <
dim; ii++)
3102 for (
int jj = 0; jj <
dim; jj++)
3104 for (
int kk = 0; kk < dof; kk++)
3105 for (
int ll = 0; ll < dof; ll++)
3107 elmat(dof*ii+kk, dof*jj+ll) +=
3108 (M * w) * gshape(kk, jj) * gshape(ll, ii);
3120 const int dof = el.
GetDof();
3122 const int tdim =
dim*(
dim+1)/2;
3125 MFEM_ASSERT(
dim == 2 ||
dim == 3,
3126 "dimension is not supported: dim = " <<
dim);
3131#ifdef MFEM_THREAD_SAFE
3137 real_t gh_data[9], grad_data[9];
3149 for (
int i = 0; i < fnd; i++)
3153 MultAtB(loc_data_mat, dshape, gh);
3174 L *= (grad(0,0) + grad(1,1));
3176 flux(i+fnd*0) = M2*grad(0,0) + L;
3177 flux(i+fnd*1) = M2*grad(1,1) + L;
3178 flux(i+fnd*2) = M*(grad(0,1) + grad(1,0));
3182 L *= (grad(0,0) + grad(1,1) + grad(2,2));
3184 flux(i+fnd*0) = M2*grad(0,0) + L;
3185 flux(i+fnd*1) = M2*grad(1,1) + L;
3186 flux(i+fnd*2) = M2*grad(2,2) + L;
3187 flux(i+fnd*3) = M*(grad(0,1) + grad(1,0));
3188 flux(i+fnd*4) = M*(grad(0,2) + grad(2,0));
3189 flux(i+fnd*5) = M*(grad(1,2) + grad(2,1));
3198 const int dof = fluxelem.
GetDof();
3200 const int tdim =
dim*(
dim+1)/2;
3205 MFEM_ASSERT(d_energy == NULL,
"anisotropic estimates are not supported");
3206 MFEM_ASSERT(flux.
Size() == dof*tdim,
"invalid 'flux' vector");
3208#ifndef MFEM_THREAD_SAFE
3213 real_t pointstress_data[6];
3214 Vector pointstress(pointstress_data, tdim);
3225 int order = 2 * Trans.
OrderGrad(&fluxelem);
3264 const real_t *
s = pointstress_data;
3268 const real_t tr_e = (
s[0] +
s[1])/(2*(M + L));
3270 pt_e = (0.25/M)*(
s[0]*(
s[0] - L) +
s[1]*(
s[1] - L) + 2*
s[2]*
s[2]);
3275 const real_t tr_e = (
s[0] +
s[1] +
s[2])/(2*M + 3*L);
3277 pt_e = (0.25/M)*(
s[0]*(
s[0] - L) +
s[1]*(
s[1] - L) +
s[2]*(
s[2] - L) +
3278 2*(
s[3]*
s[3] +
s[4]*
s[4] +
s[5]*
s[5]));
3351 nor(0) = 2*eip1.
x - 1.0;
3360 b =
beta * fabs(un);
3368 if (un >= 0.0 && ndof2)
3383 for (
int i = 0; i < ndof1; i++)
3384 for (
int j = 0; j < ndof1; j++)
3386 elmat(i, j) += w * shape1(i) * shape1(j);
3395 for (
int i = 0; i < ndof2; i++)
3396 for (
int j = 0; j < ndof1; j++)
3398 elmat(ndof1+i, j) -= w * shape2(i) * shape1(j);
3404 for (
int i = 0; i < ndof2; i++)
3405 for (
int j = 0; j < ndof2; j++)
3407 elmat(ndof1+i, ndof1+j) += w * shape2(i) * shape2(j);
3410 for (
int i = 0; i < ndof1; i++)
3411 for (
int j = 0; j < ndof2; j++)
3413 elmat(i, ndof1+j) -= w * shape1(i) * shape2(j);
3432 int ndof1, ndof2, ndofs;
3433 bool kappa_is_nonzero = (
kappa != 0.);
3463 ndofs = ndof1 + ndof2;
3466 if (kappa_is_nonzero)
3475 const int order = (ndof2) ? max(el1.
GetOrder(),
3496 nor(0) = 2*eip1.
x - 1.0;
3526 if (kappa_is_nonzero)
3542 for (
int i = 0; i < ndof1; i++)
3543 for (
int j = 0; j < ndof1; j++)
3569 if (kappa_is_nonzero)
3576 for (
int i = 0; i < ndof1; i++)
3577 for (
int j = 0; j < ndof2; j++)
3582 for (
int i = 0; i < ndof2; i++)
3583 for (
int j = 0; j < ndof1; j++)
3588 for (
int i = 0; i < ndof2; i++)
3589 for (
int j = 0; j < ndof2; j++)
3595 if (kappa_is_nonzero)
3599 for (
int i = 0; i < ndof1; i++)
3602 for (
int j = 0; j <= i; j++)
3609 for (
int i = 0; i < ndof2; i++)
3611 const int i2 = ndof1 + i;
3613 for (
int j = 0; j < ndof1; j++)
3617 for (
int j = 0; j <= i; j++)
3627 if (kappa_is_nonzero)
3629 for (
int i = 0; i < ndofs; i++)
3631 for (
int j = 0; j < i; j++)
3633 real_t aij = elmat(i,j), aji = elmat(j,i), mij =
jmat(i,j);
3634 elmat(i,j) =
sigma*aji - aij + mij;
3635 elmat(j,i) =
sigma*aij - aji + mij;
3637 elmat(i,i) = (
sigma - 1.)*elmat(i,i) +
jmat(i,i);
3642 for (
int i = 0; i < ndofs; i++)
3644 for (
int j = 0; j < i; j++)
3646 real_t aij = elmat(i,j), aji = elmat(j,i);
3647 elmat(i,j) =
sigma*aji - aij;
3648 elmat(j,i) =
sigma*aij - aji;
3650 elmat(i,i) *= (
sigma - 1.);
3665 const int dim,
const int row_ndofs,
const int col_ndofs,
3666 const int row_offset,
const int col_offset,
3672 for (
int jm = 0, j = col_offset; jm <
dim; ++jm)
3674 for (
int jdof = 0; jdof < col_ndofs; ++jdof, ++j)
3676 const real_t t2 = col_dshape_dnM(jdof);
3677 for (
int im = 0, i = row_offset; im <
dim; ++im)
3679 const real_t t1 = col_dshape(jdof, jm) * col_nL(im);
3680 const real_t t3 = col_dshape(jdof, im) * col_nM(jm);
3681 const real_t tt = t1 + ((im == jm) ? t2 : 0.0) + t3;
3682 for (
int idof = 0; idof < row_ndofs; ++idof, ++i)
3684 elmat(i, j) += row_shape(idof) * tt;
3690 if (jmatcoef == 0.0) {
return; }
3692 for (
int d = 0; d <
dim; ++d)
3694 const int jo = col_offset + d*col_ndofs;
3695 const int io = row_offset + d*row_ndofs;
3696 for (
int jdof = 0, j = jo; jdof < col_ndofs; ++jdof, ++j)
3698 const real_t sj = jmatcoef * col_shape(jdof);
3699 for (
int i = max(io,j), idof = i - io; idof < row_ndofs; ++idof, ++i)
3701 jmat(i, j) += row_shape(idof) * sj;
3711#ifdef MFEM_THREAD_SAFE
3725 const int ndofs1 = el1.
GetDof();
3727 const int nvdofs =
dim*(ndofs1 + ndofs2);
3737 const bool kappa_is_nonzero = (
kappa != 0.0);
3738 if (kappa_is_nonzero)
3771 for (
int pind = 0; pind < ir->
GetNPoints(); ++pind)
3791 nor(0) = 2*eip1.
x - 1.0;
3812 wLM = (wL2 + 2.0*wM2);
3827 wLM += (wL1 + 2.0*wM1);
3835 dim, ndofs1, ndofs1, 0, 0, jmatcoef,
nL1,
nM1,
3838 if (ndofs2 == 0) {
continue; }
3845 dim, ndofs1, ndofs2, 0,
dim*ndofs1, jmatcoef,
nL2,
nM2,
3849 dim, ndofs2, ndofs1,
dim*ndofs1, 0, jmatcoef,
nL1,
nM1,
3858 if (kappa_is_nonzero)
3860 for (
int i = 0; i < nvdofs; ++i)
3862 for (
int j = 0; j < i; ++j)
3864 real_t aij = elmat(i,j), aji = elmat(j,i), mij =
jmat(i,j);
3865 elmat(i,j) =
alpha*aji - aij + mij;
3866 elmat(j,i) =
alpha*aij - aji + mij;
3868 elmat(i,i) = (
alpha - 1.)*elmat(i,i) +
jmat(i,i);
3873 for (
int i = 0; i < nvdofs; ++i)
3875 for (
int j = 0; j < i; ++j)
3877 real_t aij = elmat(i,j), aji = elmat(j,i);
3878 elmat(i,j) =
alpha*aji - aij;
3879 elmat(j,i) =
alpha*aij - aji;
3881 elmat(i,i) *= (
alpha - 1.);
3892 int i, j, face_ndof, ndof1, ndof2;
3897 face_ndof = trial_face_fe.
GetDof();
3898 ndof1 = test_fe1.
GetDof();
3900 face_shape.
SetSize(face_ndof);
3905 ndof2 = test_fe2.
GetDof();
3913 elmat.
SetSize(ndof1 + ndof2, face_ndof);
3948 trial_face_fe.
CalcShape(ip, face_shape);
3962 for (i = 0; i < ndof1; i++)
3963 for (j = 0; j < face_ndof; j++)
3965 elmat(i, j) += shape1(i) * face_shape(j);
3970 for (i = 0; i < ndof2; i++)
3971 for (j = 0; j < face_ndof; j++)
3973 elmat(ndof1+i, j) -= shape2(i) * face_shape(j);
3984 int i, j, face_ndof, ndof1, ndof2,
dim;
3989 face_ndof = trial_face_fe.
GetDof();
3990 ndof1 = test_fe1.
GetDof();
3993 face_shape.
SetSize(face_ndof);
4000 ndof2 = test_fe2.
GetDof();
4009 elmat.
SetSize(ndof1 + ndof2, face_ndof);
4032 trial_face_fe.
CalcShape(ip, face_shape);
4038 shape1.
Mult(normal, shape1_n);
4046 shape2.
Mult(normal, shape2_n);
4049 for (i = 0; i < ndof1; i++)
4050 for (j = 0; j < face_ndof; j++)
4052 elmat(i, j) += shape1_n(i) * face_shape(j);
4057 for (i = 0; i < ndof2; i++)
4058 for (j = 0; j < face_ndof; j++)
4060 elmat(ndof1+i, j) -= shape2_n(i) * face_shape(j);
4073 "TraceIntegrator::AssembleTraceFaceMatrix: Test space should be H1");
4075 "TraceIntegrator::AssembleTraceFaceMatrix: Trial space should be RT trace");
4077 int i, j, face_ndof, ndof;
4080 face_ndof = trial_face_fe.
GetDof();
4083 face_shape.
SetSize(face_ndof);
4086 elmat.
SetSize(ndof, face_ndof);
4100 MFEM_VERIFY(elem == Trans.
Elem2->
ElementNo,
"Elem != Trans.Elem2->ElementNo");
4104 if (iel != elem) { scale = -1.; }
4119 for (i = 0; i < ndof; i++)
4121 for (j = 0; j < face_ndof; j++)
4123 elmat(i, j) += shape(i) * face_shape(j);
4135 int i, j, face_ndof, ndof,
dim;
4139 "NormalTraceIntegrator::AssembleTraceFaceMatrix: Test space should be RT");
4141 "NormalTraceIntegrator::AssembleTraceFaceMatrix: Trial space should be H1 (trace)");
4143 face_ndof = trial_face_fe.
GetDof();
4147 face_shape.
SetSize(face_ndof);
4152 elmat.
SetSize(ndof, face_ndof);
4166 MFEM_VERIFY(elem == Trans.
Elem2->
ElementNo,
"Elem != Trans.Elem2->ElementNo");
4170 if (iel != elem) { scale = -1.; }
4180 shape.
Mult(normal, shape_n);
4181 face_shape *= ip.
weight*scale;
4183 for (i = 0; i < ndof; i++)
4185 for (j = 0; j < face_ndof; j++)
4187 elmat(i, j) += shape_n(i) * face_shape(j);
4201 "TangentTraceIntegrator::AssembleTraceFaceMatrix: Test space should be ND");
4203 int face_ndof, ndof,
dim;
4209 "Trial space should be ND face trace and test space should be a ND vector field in 3D ";
4211 trial_face_fe.
GetDim() == 2 && test_fe.
GetDim() == 3, msg);
4216 "Trial space should be H1 edge trace and test space should be a ND vector field in 2D";
4218 trial_face_fe.
GetDim() == 1 && test_fe.
GetDim() == 2, msg);
4220 face_ndof = trial_face_fe.
GetDof();
4223 int dimc = (
dim == 3) ? 3 : 1;
4231 elmat.
SetSize(ndof, face_ndof);
4245 MFEM_VERIFY(elem == Trans.
Elem2->
ElementNo,
"Elem != Trans.Elem2->ElementNo");
4249 if (iel != elem) { scale = -1.; }
4270 cross_product(normal, shape, shape_n);
4286 for (
int i = 0; i < ran_nodes.
Size(); i++)
4292 for (
int j = 0; j < shape.
Size(); j++)
4294 for (
int d = 0; d < spaceDim; d++)
4296 elmat(i, j+d*shape.
Size()) = shape(j)*n(d);
4317 virtual void Eval(Vector &V, ElementTransformation &T,
4318 const IntegrationPoint &ip)
4334 internal::ShapeCoefficient dom_shape_coeff(*
Q, dom_fe);
4340 ran_fe.
Project(dom_shape_coeff, Trans, elmat_as_vec);
4369 VShapeCoefficient dom_shape_coeff(*
Q, dom_fe, Trans.
GetSpaceDim());
4395 vc(width), shape(height) { }
4407 VecShapeCoefficient dom_shape_coeff(*
VQ, dom_fe);
4436 using VectorCoefficient::Eval;
4443 for (
int k = 0; k < vdim; k++)
4445 V(k) = vc(0) * vshape(k,1) - vc(1) * vshape(k,0);
4450 VCrossVShapeCoefficient dom_shape_coeff(*
VQ, dom_fe);
4456 ran_fe.
Project(dom_shape_coeff, Trans, elmat_as_vec);
4476 vshape(height, width), vc(width)
4478 MFEM_ASSERT(width == 3,
"");
4487 for (
int k = 0; k < height; k++)
4489 M(k,0) = vc(1) * vshape(k,2) - vc(2) * vshape(k,1);
4490 M(k,1) = vc(2) * vshape(k,0) - vc(0) * vshape(k,2);
4491 M(k,2) = vc(0) * vshape(k,1) - vc(1) * vshape(k,0);
4496 VCrossVShapeCoefficient dom_shape_coeff(*
VQ, dom_fe);
4530 virtual void Eval(Vector &V, ElementTransformation &T,
4531 const IntegrationPoint &ip)
4549 internal::VDotVShapeCoefficient dom_shape_coeff(*
VQ, dom_fe);
4555 ran_fe.
Project(dom_shape_coeff, Trans, elmat_as_vec);
int Size() const
Return the logical size of the array.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
static const IntegrationRule & GetRule(const FiniteElement &el, ElementTransformation &Trans)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual real_t ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
DiagonalMatrixCoefficient * DQ
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat) override
const IntegrationRule & GetRule(int order, FaceElementTransformations &T)
static void AssembleBlock(const int dim, const int row_ndofs, const int col_ndofs, const int row_offset, const int col_offset, const real_t jmatcoef, const Vector &col_nL, const Vector &col_nM, const Vector &row_shape, const Vector &col_shape, const Vector &col_dshape_dnM, const DenseMatrix &col_dshape, DenseMatrix &elmat, DenseMatrix &jmat)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
const FaceGeometricFactors * geom
Not owned.
static const IntegrationRule & GetRule(Geometry::Type geom, int order, FaceElementTransformations &T)
Data type dense matrix using column-major storage.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
void Transpose()
(*this) = (*this)^t
void GetColumnReference(int c, Vector &col)
real_t InnerProduct(const real_t *x, const real_t *y) const
Compute y^t A x.
real_t * Data() const
Returns the matrix data array.
real_t * GetData() const
Returns the matrix data array.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
virtual void AddMult(const Vector &x, Vector &y, const real_t a=1.0) const
y += a * A.x
void Reset(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
void Invert()
Replaces the current matrix with its inverse.
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i.
void GradToCurl(DenseMatrix &curl)
void GetColumn(int c, Vector &col) const
void GradToVectorCurl2D(DenseMatrix &curl)
void GradToDiv(Vector &div)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
virtual real_t ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
ElasticityComponentIntegrator(ElasticityIntegrator &parent_, int i_, int j_)
Given an ElasticityIntegrator, create an integrator that represents the th component block.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Tr, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual real_t ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Abstract class for all finite elements.
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 ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
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...
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...
int GetDim() const
Returns the reference space dimension for the finite element.
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
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...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference 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...
int Space() const
Returns the type of FunctionSpace on the element.
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 CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
int GetCurlDim() const
Returns the dimension of the curl for vector-valued finite elements.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement....
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 GetDof() const
Returns the number of degrees of freedom in the finite element.
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 ...
@ Pk
Polynomials of order k.
@ rQk
Refined tensor products of polynomials of order k.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Base class for Matrix Coefficients that optionally depend on time and space.
int GetVDim() const
For backward compatibility get the width of the matrix.
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 GetWidth() const
Get the width of the matrix.
int GetHeight() const
Get the height of the matrix.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape_)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual int GetVDim(const FiniteElement &vector_fe)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape_)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual int GetTestVDim(const FiniteElement &test_fe)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual int GetTrialVDim(const FiniteElement &trial_fe)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
DiagonalMatrixCoefficient * DQ
An arbitrary order and dimension NURBS element.
Array< const KnotVector * > & KnotVectors() const
const int * GetIJK() const
IntegrationRule & GetElementRule(const int elem, const int patch, const int *ijk, Array< const KnotVector * > const &kv, bool &deleteRule) const
Returns a rule for the element.
Class for standard nodal finite elements.
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleTraceFaceMatrix(int ielem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void AssembleElementMatrix2(const FiniteElement &nd_fe, const FiniteElement &l2_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
virtual void AddMultMF(const Vector &x, Vector &y) const
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
virtual void AssembleMF(const FiniteElementSpace &fes)
Method defining matrix-free assembly.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleEAInteriorFaces(const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add)
virtual void AddMultTransposeMF(const Vector &x, Vector &y) const
virtual void AssembleDiagonalMF(Vector &diag)
Assemble diagonal and add it to Vector diag.
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
virtual void AssembleEA(const FiniteElementSpace &fes, Vector &emat, const bool add)
Method defining element assembly.
virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes)
virtual void AssembleEABoundaryFaces(const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add)
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes)
void AssembleTraceFaceMatrix(int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
void AssembleTraceFaceMatrix(int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
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 AssembleElementMatrix2(const FiniteElement &nd_fe, const FiniteElement &rt_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Assemble an element matrix.
virtual real_t GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute element energy: .
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the b...
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
DiagonalMatrixCoefficient * DQ
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &rt_fe, const FiniteElement &l2_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
void Neg()
(*this) = -(*this)
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
void NewDataAndSize(real_t *d, int s)
Set the Vector data and size, deleting the old data, if owned.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
void CalcOrtho(const DenseMatrix &J, Vector &n)
void AddMult_a_ABt(real_t a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
real_t u(const Vector &xvec)
void mfem_error(const char *msg)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void add(const Vector &v1, const Vector &v2, Vector &v)
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
void AddMult_a_VWt(const real_t a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
void AddMult_a_VVt(const real_t a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
void AddMult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
void AddMultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
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 Mult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)