12 #include "../general/forall.hpp"
15 #include "libceed/mass.hpp"
27 static void PAHcurlSetup2D(
const int Q1D,
34 const int NQ = Q1D*Q1D;
43 for (
int q = 0; q < NQ; ++q)
45 const double J11 = J(q,0,0,e);
46 const double J21 = J(q,1,0,e);
47 const double J12 = J(q,0,1,e);
48 const double J22 = J(q,1,1,e);
49 const double c_detJ = W[q] * coeff(q, e) / ((J11*J22)-(J21*J12));
50 y(q,0,e) = c_detJ * (J12*J12 + J22*J22);
51 y(q,1,e) = -c_detJ * (J12*J11 + J22*J21);
52 y(q,2,e) = c_detJ * (J11*J11 + J21*J21);
58 static void PAHcurlSetup3D(
const int Q1D,
60 const Array<double> &w,
65 const int NQ = Q1D*Q1D*Q1D;
67 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
68 auto coeff =
Reshape(_coeff.Read(), NQ, NE);
69 auto y =
Reshape(op.Write(), NQ, 6, NE);
73 for (
int q = 0; q < NQ; ++q)
75 const double J11 = J(q,0,0,e);
76 const double J21 = J(q,1,0,e);
77 const double J31 = J(q,2,0,e);
78 const double J12 = J(q,0,1,e);
79 const double J22 = J(q,1,1,e);
80 const double J32 = J(q,2,1,e);
81 const double J13 = J(q,0,2,e);
82 const double J23 = J(q,1,2,e);
83 const double J33 = J(q,2,2,e);
84 const double detJ = J11 * (J22 * J33 - J32 * J23) -
85 J21 * (J12 * J33 - J32 * J13) +
86 J31 * (J12 * J23 - J22 * J13);
87 const double c_detJ = W[q] * coeff(q, e) / detJ;
89 const double A11 = (J22 * J33) - (J23 * J32);
90 const double A12 = (J32 * J13) - (J12 * J33);
91 const double A13 = (J12 * J23) - (J22 * J13);
92 const double A21 = (J31 * J23) - (J21 * J33);
93 const double A22 = (J11 * J33) - (J13 * J31);
94 const double A23 = (J21 * J13) - (J11 * J23);
95 const double A31 = (J21 * J32) - (J31 * J22);
96 const double A32 = (J31 * J12) - (J11 * J32);
97 const double A33 = (J11 * J22) - (J12 * J21);
99 y(q,0,e) = c_detJ * (A11*A11 + A12*A12 + A13*A13);
100 y(q,1,e) = c_detJ * (A11*A21 + A12*A22 + A13*A23);
101 y(q,2,e) = c_detJ * (A11*A31 + A12*A32 + A13*A33);
102 y(q,3,e) = c_detJ * (A21*A21 + A22*A22 + A23*A23);
103 y(q,4,e) = c_detJ * (A21*A31 + A22*A32 + A23*A33);
104 y(q,5,e) = c_detJ * (A31*A31 + A32*A32 + A33*A33);
117 MFEM_VERIFY(el != NULL,
"Only VectorTensorFiniteElement is supported!");
120 = IntRule ? IntRule : &MassIntegrator::GetRule(*el, *el,
122 const int dims = el->
GetDim();
123 MFEM_VERIFY(dims == 2 || dims == 3,
"");
125 const int symmDims = (dims * (dims + 1)) / 2;
128 MFEM_VERIFY(
dim == 2 ||
dim == 3,
"");
134 dofs1D = mapsC->
ndof;
135 quad1D = mapsC->nqpt;
137 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt,
"");
145 for (
int e=0; e<ne; ++e)
148 for (
int p=0; p<nq; ++p)
150 coeff[p + (e * nq)] = Q->Eval(*tr, ir->
IntPoint(p));
167 MFEM_ABORT(
"Unknown kernel.");
171 static void PAHcurlMassApply2D(
const int D1D,
182 constexpr
static int VDIM = 2;
196 for (
int qy = 0; qy < Q1D; ++qy)
198 for (
int qx = 0; qx < Q1D; ++qx)
200 for (
int c = 0; c < VDIM; ++c)
202 mass[qy][qx][c] = 0.0;
209 for (
int c = 0; c < VDIM; ++c)
211 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
212 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
214 for (
int dy = 0; dy < D1Dy; ++dy)
217 for (
int qx = 0; qx < Q1D; ++qx)
222 for (
int dx = 0; dx < D1Dx; ++dx)
224 const double t = x(dx + (dy * D1Dx) + osc, e);
225 for (
int qx = 0; qx < Q1D; ++qx)
227 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
231 for (
int qy = 0; qy < Q1D; ++qy)
233 const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
234 for (
int qx = 0; qx < Q1D; ++qx)
236 mass[qy][qx][c] += massX[qx] * wy;
245 for (
int qy = 0; qy < Q1D; ++qy)
247 for (
int qx = 0; qx < Q1D; ++qx)
249 const double O11 = op(qx,qy,0,e);
250 const double O12 = op(qx,qy,1,e);
251 const double O22 = op(qx,qy,2,e);
252 const double massX = mass[qy][qx][0];
253 const double massY = mass[qy][qx][1];
254 mass[qy][qx][0] = (O11*massX)+(O12*massY);
255 mass[qy][qx][1] = (O12*massX)+(O22*massY);
259 for (
int qy = 0; qy < Q1D; ++qy)
263 for (
int c = 0; c < VDIM; ++c)
265 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
266 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
269 for (
int dx = 0; dx < D1Dx; ++dx)
273 for (
int qx = 0; qx < Q1D; ++qx)
275 for (
int dx = 0; dx < D1Dx; ++dx)
277 massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
281 for (
int dy = 0; dy < D1Dy; ++dy)
283 const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
285 for (
int dx = 0; dx < D1Dx; ++dx)
287 y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
297 static void PAHcurlMassAssembleDiagonal2D(
const int D1D,
300 const Array<double> &_Bo,
301 const Array<double> &_Bc,
305 constexpr
static int VDIM = 2;
307 auto Bo =
Reshape(_Bo.Read(), Q1D, D1D-1);
308 auto Bc =
Reshape(_Bc.Read(), Q1D, D1D);
309 auto op =
Reshape(_op.Read(), Q1D, Q1D, 3, NE);
310 auto diag =
Reshape(_diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
316 for (
int c = 0; c < VDIM; ++c)
318 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
319 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
323 for (
int dy = 0; dy < D1Dy; ++dy)
325 for (
int qx = 0; qx < Q1D; ++qx)
328 for (
int qy = 0; qy < Q1D; ++qy)
330 const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
332 mass[qx] += wy * wy * ((c == 0) ? op(qx,qy,0,e) : op(qx,qy,2,e));
336 for (
int dx = 0; dx < D1Dx; ++dx)
338 for (
int qx = 0; qx < Q1D; ++qx)
340 const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
341 diag(dx + (dy * D1Dx) + osc, e) += mass[qx] * wx * wx;
351 template<
int MAX_D1D = HCURL_MAX_D1D,
int MAX_Q1D = HCURL_MAX_Q1D>
352 static void PAHcurlMassAssembleDiagonal3D(
const int D1D,
355 const Array<double> &_Bo,
356 const Array<double> &_Bc,
360 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D >
MAX_D1D");
361 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D >
MAX_Q1D");
362 constexpr static int VDIM = 3;
364 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
365 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
366 auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
367 auto diag = Reshape(_diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
373 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
375 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
376 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
377 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
379 const int opc = (c == 0) ? 0 : ((c == 1) ? 3 : 5);
381 double mass[MAX_Q1D];
383 for (int dz = 0; dz < D1Dz; ++dz)
385 for (int dy = 0; dy < D1Dy; ++dy)
387 for (int qx = 0; qx < Q1D; ++qx)
390 for (int qy = 0; qy < Q1D; ++qy)
392 const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
394 for (int qz = 0; qz < Q1D; ++qz)
396 const double wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
398 mass[qx] += wy * wy * wz * wz * op(qx,qy,qz,opc,e);
403 for (int dx = 0; dx < D1Dx; ++dx)
405 for (int qx = 0; qx < Q1D; ++qx)
407 const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
408 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += mass[qx] * wx * wx;
414 osc += D1Dx * D1Dy * D1Dz;
416 }); // end of element loop
419 void VectorFEMassIntegrator::AssembleDiagonalPA(Vector& diag)
422 PAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne,
423 mapsO->B, mapsC->B, pa_data, diag);
425 PAHcurlMassAssembleDiagonal2D(dofs1D, quad1D, ne,
426 mapsO->B, mapsC->B, pa_data, diag);
429 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
430 static void PAHcurlMassApply3D(const int D1D,
433 const Array<double> &_Bo,
434 const Array<double> &_Bc,
435 const Array<double> &_Bot,
436 const Array<double> &_Bct,
441 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D >
MAX_D1D");
442 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D >
MAX_Q1D");
443 constexpr static int VDIM = 3;
445 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
446 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
447 auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
448 auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
449 auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
450 auto x = Reshape(_x.Read(), 3*(D1D-1)*D1D*D1D, NE);
451 auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
455 double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
457 for (int qz = 0; qz < Q1D; ++qz)
459 for (int qy = 0; qy < Q1D; ++qy)
461 for (int qx = 0; qx < Q1D; ++qx)
463 for (int c = 0; c < VDIM; ++c)
465 mass[qz][qy][qx][c] = 0.0;
473 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
475 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
476 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
477 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
479 for (int dz = 0; dz < D1Dz; ++dz)
481 double massXY[MAX_Q1D][MAX_Q1D];
482 for (int qy = 0; qy < Q1D; ++qy)
484 for (int qx = 0; qx < Q1D; ++qx)
486 massXY[qy][qx] = 0.0;
490 for (int dy = 0; dy < D1Dy; ++dy)
492 double massX[MAX_Q1D];
493 for (int qx = 0; qx < Q1D; ++qx)
498 for (int dx = 0; dx < D1Dx; ++dx)
500 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
501 for (int qx = 0; qx < Q1D; ++qx)
503 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
507 for (int qy = 0; qy < Q1D; ++qy)
509 const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
510 for (int qx = 0; qx < Q1D; ++qx)
512 const double wx = massX[qx];
513 massXY[qy][qx] += wx * wy;
518 for (int qz = 0; qz < Q1D; ++qz)
520 const double wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
521 for (int qy = 0; qy < Q1D; ++qy)
523 for (int qx = 0; qx < Q1D; ++qx)
525 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
531 osc += D1Dx * D1Dy * D1Dz;
532 } // loop (c) over components
535 for (int qz = 0; qz < Q1D; ++qz)
537 for (int qy = 0; qy < Q1D; ++qy)
539 for (int qx = 0; qx < Q1D; ++qx)
541 const double O11 = op(qx,qy,qz,0,e);
542 const double O12 = op(qx,qy,qz,1,e);
543 const double O13 = op(qx,qy,qz,2,e);
544 const double O22 = op(qx,qy,qz,3,e);
545 const double O23 = op(qx,qy,qz,4,e);
546 const double O33 = op(qx,qy,qz,5,e);
547 const double massX = mass[qz][qy][qx][0];
548 const double massY = mass[qz][qy][qx][1];
549 const double massZ = mass[qz][qy][qx][2];
550 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
551 mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ);
552 mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ);
557 for (int qz = 0; qz < Q1D; ++qz)
559 double massXY[MAX_D1D][MAX_D1D];
563 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
565 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
566 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
567 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
569 for (int dy = 0; dy < D1Dy; ++dy)
571 for (int dx = 0; dx < D1Dx; ++dx)
576 for (int qy = 0; qy < Q1D; ++qy)
578 double massX[MAX_D1D];
579 for (int dx = 0; dx < D1Dx; ++dx)
583 for (int qx = 0; qx < Q1D; ++qx)
585 for (int dx = 0; dx < D1Dx; ++dx)
587 massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
590 for (int dy = 0; dy < D1Dy; ++dy)
592 const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
593 for (int dx = 0; dx < D1Dx; ++dx)
595 massXY[dy][dx] += massX[dx] * wy;
600 for (int dz = 0; dz < D1Dz; ++dz)
602 const double wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
603 for (int dy = 0; dy < D1Dy; ++dy)
605 for (int dx = 0; dx < D1Dx; ++dx)
607 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
612 osc += D1Dx * D1Dy * D1Dz;
615 }); // end of element loop
618 void VectorFEMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
622 PAHcurlMassApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
623 mapsC->Bt, pa_data, x, y);
627 PAHcurlMassApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
628 mapsC->Bt, pa_data, x, y);
632 // PA H(curl) curl-curl assemble 2D kernel
633 static void PACurlCurlSetup2D(const int Q1D,
635 const Array<double> &w,
640 const int NQ = Q1D*Q1D;
642 auto J = Reshape(j.Read(), NQ, 2, 2, NE);
643 auto coeff = Reshape(_coeff.Read(), NQ, NE);
644 auto y = Reshape(op.Write(), NQ, NE);
647 for (int q = 0; q < NQ; ++q)
649 const double J11 = J(q,0,0,e);
650 const double J21 = J(q,1,0,e);
651 const double J12 = J(q,0,1,e);
652 const double J22 = J(q,1,1,e);
653 const double detJ = (J11*J22)-(J21*J12);
654 y(q,e) = W[q] * coeff(q,e) / detJ;
659 // PA H(curl) curl-curl assemble 3D kernel
660 static void PACurlCurlSetup3D(const int Q1D,
662 const Array<double> &w,
667 const int NQ = Q1D*Q1D*Q1D;
669 auto J = Reshape(j.Read(), NQ, 3, 3, NE);
670 auto coeff = Reshape(_coeff.Read(), NQ, NE);
671 auto y = Reshape(op.Write(), NQ, 6, NE);
674 for (int q = 0; q < NQ; ++q)
676 const double J11 = J(q,0,0,e);
677 const double J21 = J(q,1,0,e);
678 const double J31 = J(q,2,0,e);
679 const double J12 = J(q,0,1,e);
680 const double J22 = J(q,1,1,e);
681 const double J32 = J(q,2,1,e);
682 const double J13 = J(q,0,2,e);
683 const double J23 = J(q,1,2,e);
684 const double J33 = J(q,2,2,e);
685 const double detJ = J11 * (J22 * J33 - J32 * J23) -
686 /* */ J21 * (J12 * J33 - J32 * J13) +
687 /* */ J31 * (J12 * J23 - J22 * J13);
689 // set y to the 6 entries of J^T J / det^2
690 const double c_detJ = W[q] * coeff(q,e) / detJ;
692 y(q,0,e) = c_detJ * (J11*J11 + J21*J21 + J31*J31); // 1,1
693 y(q,1,e) = c_detJ * (J11*J12 + J21*J22 + J31*J32); // 1,2
694 y(q,2,e) = c_detJ * (J11*J13 + J21*J23 + J31*J33); // 1,3
695 y(q,3,e) = c_detJ * (J12*J12 + J22*J22 + J32*J32); // 2,2
696 y(q,4,e) = c_detJ * (J12*J13 + J22*J23 + J32*J33); // 2,3
697 y(q,5,e) = c_detJ * (J13*J13 + J23*J23 + J33*J33); // 3,3
702 void CurlCurlIntegrator::AssemblePA(const FiniteElementSpace &fes)
704 // Assumes tensor-product elements
705 Mesh *mesh = fes.GetMesh();
706 const FiniteElement *fel = fes.GetFE(0);
708 const VectorTensorFiniteElement *el =
709 dynamic_cast<const VectorTensorFiniteElement*>(fel);
712 const IntegrationRule *ir
713 = IntRule ? IntRule : &MassIntegrator::GetRule(*el, *el,
714 *mesh->GetElementTransformation(0));
715 const int dims = el->GetDim();
716 MFEM_VERIFY(dims == 2 || dims == 3, "");
718 const int nq = ir->GetNPoints();
719 dim = mesh->Dimension();
720 MFEM_VERIFY(dim == 2 || dim == 3, "");
723 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
724 mapsC = &el->GetDofToQuad(*ir, DofToQuad::TENSOR);
725 mapsO = &el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
726 dofs1D = mapsC->ndof;
727 quad1D = mapsC->nqpt;
729 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
731 const int ndata = (dim == 2) ? 1 : 6;
732 pa_data.SetSize(ndata * nq * ne, Device::GetMemoryType());
734 Vector coeff(ne * nq);
738 for (int e=0; e<ne; ++e)
740 ElementTransformation *tr = mesh->GetElementTransformation(e);
741 for (int p=0; p<nq; ++p)
743 coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
748 if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
750 // pa_data_2.SetSize(6 * nq * ne, Device::GetMemoryType());
752 PACurlCurlSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
755 else if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
757 PACurlCurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
762 MFEM_ABORT("Unknown kernel.
");
766 static void PACurlCurlApply2D(const int D1D,
769 const Array<double> &_Bo,
770 const Array<double> &_Bot,
771 const Array<double> &_Gc,
772 const Array<double> &_Gct,
777 constexpr static int VDIM = 2;
779 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
780 auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
781 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
782 auto Gct = Reshape(_Gct.Read(), D1D, Q1D);
783 auto op = Reshape(_op.Read(), Q1D, Q1D, NE);
784 auto x = Reshape(_x.Read(), 2*(D1D-1)*D1D, NE);
785 auto y = Reshape(_y.ReadWrite(), 2*(D1D-1)*D1D, NE);
789 double curl[MAX_Q1D][MAX_Q1D];
791 // curl[qy][qx] will be computed as du_y/dx - du_x/dy
793 for (int qy = 0; qy < Q1D; ++qy)
795 for (int qx = 0; qx < Q1D; ++qx)
803 for (int c = 0; c < VDIM; ++c) // loop over x, y components
805 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
806 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
808 for (int dy = 0; dy < D1Dy; ++dy)
810 double gradX[MAX_Q1D];
811 for (int qx = 0; qx < Q1D; ++qx)
816 for (int dx = 0; dx < D1Dx; ++dx)
818 const double t = x(dx + (dy * D1Dx) + osc, e);
819 for (int qx = 0; qx < Q1D; ++qx)
821 gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
825 for (int qy = 0; qy < Q1D; ++qy)
827 const double wy = (c == 0) ? -Gc(qy,dy) : Bo(qy,dy);
828 for (int qx = 0; qx < Q1D; ++qx)
830 curl[qy][qx] += gradX[qx] * wy;
836 } // loop (c) over components
839 for (int qy = 0; qy < Q1D; ++qy)
841 for (int qx = 0; qx < Q1D; ++qx)
843 curl[qy][qx] *= op(qx,qy,e);
847 for (int qy = 0; qy < Q1D; ++qy)
851 for (int c = 0; c < VDIM; ++c) // loop over x, y components
853 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
854 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
856 double gradX[MAX_D1D];
857 for (int dx = 0; dx < D1Dx; ++dx)
861 for (int qx = 0; qx < Q1D; ++qx)
863 for (int dx = 0; dx < D1Dx; ++dx)
865 gradX[dx] += curl[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
868 for (int dy = 0; dy < D1Dy; ++dy)
870 const double wy = (c == 0) ? -Gct(dy,qy) : Bot(dy,qy);
872 for (int dx = 0; dx < D1Dx; ++dx)
874 y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
881 }); // end of element loop
884 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
885 static void PACurlCurlApply3D(const int D1D,
888 const Array<double> &_Bo,
889 const Array<double> &_Bc,
890 const Array<double> &_Bot,
891 const Array<double> &_Bct,
892 const Array<double> &_Gc,
893 const Array<double> &_Gct,
898 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D >
MAX_D1D");
899 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D >
MAX_Q1D");
900 // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
901 // (\nabla\times u) \cdot (\nabla\times v) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{v}
902 // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
903 // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
904 // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
906 constexpr static int VDIM = 3;
908 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
909 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
910 auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
911 auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
912 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
913 auto Gct = Reshape(_Gct.Read(), D1D, Q1D);
914 auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
915 auto x = Reshape(_x.Read(), 3*(D1D-1)*D1D*D1D, NE);
916 auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
920 double curl[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
921 // curl[qz][qy][qx] will be computed as the vector curl at each quadrature point.
923 for (int qz = 0; qz < Q1D; ++qz)
925 for (int qy = 0; qy < Q1D; ++qy)
927 for (int qx = 0; qx < Q1D; ++qx)
929 for (int c = 0; c < VDIM; ++c)
931 curl[qz][qy][qx][c] = 0.0;
937 // We treat x, y, z components separately for optimization specific to each.
943 const int D1Dz = D1D;
944 const int D1Dy = D1D;
945 const int D1Dx = D1D - 1;
947 for (int dz = 0; dz < D1Dz; ++dz)
949 double gradXY[MAX_Q1D][MAX_Q1D][2];
950 for (int qy = 0; qy < Q1D; ++qy)
952 for (int qx = 0; qx < Q1D; ++qx)
954 for (int d = 0; d < 2; ++d)
956 gradXY[qy][qx][d] = 0.0;
961 for (int dy = 0; dy < D1Dy; ++dy)
963 double massX[MAX_Q1D];
964 for (int qx = 0; qx < Q1D; ++qx)
969 for (int dx = 0; dx < D1Dx; ++dx)
971 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
972 for (int qx = 0; qx < Q1D; ++qx)
974 massX[qx] += t * Bo(qx,dx);
978 for (int qy = 0; qy < Q1D; ++qy)
980 const double wy = Bc(qy,dy);
981 const double wDy = Gc(qy,dy);
982 for (int qx = 0; qx < Q1D; ++qx)
984 const double wx = massX[qx];
985 gradXY[qy][qx][0] += wx * wDy;
986 gradXY[qy][qx][1] += wx * wy;
991 for (int qz = 0; qz < Q1D; ++qz)
993 const double wz = Bc(qz,dz);
994 const double wDz = Gc(qz,dz);
995 for (int qy = 0; qy < Q1D; ++qy)
997 for (int qx = 0; qx < Q1D; ++qx)
999 // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1000 curl[qz][qy][qx][1] += gradXY[qy][qx][1] * wDz; // (u_0)_{x_2}
1001 curl[qz][qy][qx][2] -= gradXY[qy][qx][0] * wz; // -(u_0)_{x_1}
1007 osc += D1Dx * D1Dy * D1Dz;
1012 const int D1Dz = D1D;
1013 const int D1Dy = D1D - 1;
1014 const int D1Dx = D1D;
1016 for (int dz = 0; dz < D1Dz; ++dz)
1018 double gradXY[MAX_Q1D][MAX_Q1D][2];
1019 for (int qy = 0; qy < Q1D; ++qy)
1021 for (int qx = 0; qx < Q1D; ++qx)
1023 for (int d = 0; d < 2; ++d)
1025 gradXY[qy][qx][d] = 0.0;
1030 for (int dx = 0; dx < D1Dx; ++dx)
1032 double massY[MAX_Q1D];
1033 for (int qy = 0; qy < Q1D; ++qy)
1038 for (int dy = 0; dy < D1Dy; ++dy)
1040 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
1041 for (int qy = 0; qy < Q1D; ++qy)
1043 massY[qy] += t * Bo(qy,dy);
1047 for (int qx = 0; qx < Q1D; ++qx)
1049 const double wx = Bc(qx,dx);
1050 const double wDx = Gc(qx,dx);
1051 for (int qy = 0; qy < Q1D; ++qy)
1053 const double wy = massY[qy];
1054 gradXY[qy][qx][0] += wDx * wy;
1055 gradXY[qy][qx][1] += wx * wy;
1060 for (int qz = 0; qz < Q1D; ++qz)
1062 const double wz = Bc(qz,dz);
1063 const double wDz = Gc(qz,dz);
1064 for (int qy = 0; qy < Q1D; ++qy)
1066 for (int qx = 0; qx < Q1D; ++qx)
1068 // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1069 curl[qz][qy][qx][0] -= gradXY[qy][qx][1] * wDz; // -(u_1)_{x_2}
1070 curl[qz][qy][qx][2] += gradXY[qy][qx][0] * wz; // (u_1)_{x_0}
1076 osc += D1Dx * D1Dy * D1Dz;
1081 const int D1Dz = D1D - 1;
1082 const int D1Dy = D1D;
1083 const int D1Dx = D1D;
1085 for (int dx = 0; dx < D1Dx; ++dx)
1087 double gradYZ[MAX_Q1D][MAX_Q1D][2];
1088 for (int qz = 0; qz < Q1D; ++qz)
1090 for (int qy = 0; qy < Q1D; ++qy)
1092 for (int d = 0; d < 2; ++d)
1094 gradYZ[qz][qy][d] = 0.0;
1099 for (int dy = 0; dy < D1Dy; ++dy)
1101 double massZ[MAX_Q1D];
1102 for (int qz = 0; qz < Q1D; ++qz)
1107 for (int dz = 0; dz < D1Dz; ++dz)
1109 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
1110 for (int qz = 0; qz < Q1D; ++qz)
1112 massZ[qz] += t * Bo(qz,dz);
1116 for (int qy = 0; qy < Q1D; ++qy)
1118 const double wy = Bc(qy,dy);
1119 const double wDy = Gc(qy,dy);
1120 for (int qz = 0; qz < Q1D; ++qz)
1122 const double wz = massZ[qz];
1123 gradYZ[qz][qy][0] += wz * wy;
1124 gradYZ[qz][qy][1] += wz * wDy;
1129 for (int qx = 0; qx < Q1D; ++qx)
1131 const double wx = Bc(qx,dx);
1132 const double wDx = Gc(qx,dx);
1134 for (int qy = 0; qy < Q1D; ++qy)
1136 for (int qz = 0; qz < Q1D; ++qz)
1138 // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1139 curl[qz][qy][qx][0] += gradYZ[qz][qy][1] * wx; // (u_2)_{x_1}
1140 curl[qz][qy][qx][1] -= gradYZ[qz][qy][0] * wDx; // -(u_2)_{x_0}
1147 // Apply D operator.
1148 for (int qz = 0; qz < Q1D; ++qz)
1150 for (int qy = 0; qy < Q1D; ++qy)
1152 for (int qx = 0; qx < Q1D; ++qx)
1154 const double O11 = op(qx,qy,qz,0,e);
1155 const double O12 = op(qx,qy,qz,1,e);
1156 const double O13 = op(qx,qy,qz,2,e);
1157 const double O22 = op(qx,qy,qz,3,e);
1158 const double O23 = op(qx,qy,qz,4,e);
1159 const double O33 = op(qx,qy,qz,5,e);
1161 const double c1 = (O11 * curl[qz][qy][qx][0]) + (O12 * curl[qz][qy][qx][1]) +
1162 (O13 * curl[qz][qy][qx][2]);
1163 const double c2 = (O12 * curl[qz][qy][qx][0]) + (O22 * curl[qz][qy][qx][1]) +
1164 (O23 * curl[qz][qy][qx][2]);
1165 const double c3 = (O13 * curl[qz][qy][qx][0]) + (O23 * curl[qz][qy][qx][1]) +
1166 (O33 * curl[qz][qy][qx][2]);
1168 curl[qz][qy][qx][0] = c1;
1169 curl[qz][qy][qx][1] = c2;
1170 curl[qz][qy][qx][2] = c3;
1178 const int D1Dz = D1D;
1179 const int D1Dy = D1D;
1180 const int D1Dx = D1D - 1;
1182 for (int qz = 0; qz < Q1D; ++qz)
1184 double gradXY12[MAX_D1D][MAX_D1D];
1185 double gradXY21[MAX_D1D][MAX_D1D];
1187 for (int dy = 0; dy < D1Dy; ++dy)
1189 for (int dx = 0; dx < D1Dx; ++dx)
1191 gradXY12[dy][dx] = 0.0;
1192 gradXY21[dy][dx] = 0.0;
1195 for (int qy = 0; qy < Q1D; ++qy)
1197 double massX[MAX_D1D][2];
1198 for (int dx = 0; dx < D1Dx; ++dx)
1200 for (int n = 0; n < 2; ++n)
1205 for (int qx = 0; qx < Q1D; ++qx)
1207 for (int dx = 0; dx < D1Dx; ++dx)
1209 const double wx = Bot(dx,qx);
1211 massX[dx][0] += wx * curl[qz][qy][qx][1];
1212 massX[dx][1] += wx * curl[qz][qy][qx][2];
1215 for (int dy = 0; dy < D1Dy; ++dy)
1217 const double wy = Bct(dy,qy);
1218 const double wDy = Gct(dy,qy);
1220 for (int dx = 0; dx < D1Dx; ++dx)
1222 gradXY21[dy][dx] += massX[dx][0] * wy;
1223 gradXY12[dy][dx] += massX[dx][1] * wDy;
1228 for (int dz = 0; dz < D1Dz; ++dz)
1230 const double wz = Bct(dz,qz);
1231 const double wDz = Gct(dz,qz);
1232 for (int dy = 0; dy < D1Dy; ++dy)
1234 for (int dx = 0; dx < D1Dx; ++dx)
1236 // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1237 // (u_0)_{x_2} * (op * curl)_1 - (u_0)_{x_1} * (op * curl)_2
1238 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1239 e) += (gradXY21[dy][dx] * wDz) - (gradXY12[dy][dx] * wz);
1245 osc += D1Dx * D1Dy * D1Dz;
1250 const int D1Dz = D1D;
1251 const int D1Dy = D1D - 1;
1252 const int D1Dx = D1D;
1254 for (int qz = 0; qz < Q1D; ++qz)
1256 double gradXY02[MAX_D1D][MAX_D1D];
1257 double gradXY20[MAX_D1D][MAX_D1D];
1259 for (int dy = 0; dy < D1Dy; ++dy)
1261 for (int dx = 0; dx < D1Dx; ++dx)
1263 gradXY02[dy][dx] = 0.0;
1264 gradXY20[dy][dx] = 0.0;
1267 for (int qx = 0; qx < Q1D; ++qx)
1269 double massY[MAX_D1D][2];
1270 for (int dy = 0; dy < D1Dy; ++dy)
1275 for (int qy = 0; qy < Q1D; ++qy)
1277 for (int dy = 0; dy < D1Dy; ++dy)
1279 const double wy = Bot(dy,qy);
1281 massY[dy][0] += wy * curl[qz][qy][qx][2];
1282 massY[dy][1] += wy * curl[qz][qy][qx][0];
1285 for (int dx = 0; dx < D1Dx; ++dx)
1287 const double wx = Bct(dx,qx);
1288 const double wDx = Gct(dx,qx);
1290 for (int dy = 0; dy < D1Dy; ++dy)
1292 gradXY02[dy][dx] += massY[dy][0] * wDx;
1293 gradXY20[dy][dx] += massY[dy][1] * wx;
1298 for (int dz = 0; dz < D1Dz; ++dz)
1300 const double wz = Bct(dz,qz);
1301 const double wDz = Gct(dz,qz);
1302 for (int dy = 0; dy < D1Dy; ++dy)
1304 for (int dx = 0; dx < D1Dx; ++dx)
1306 // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1307 // -(u_1)_{x_2} * (op * curl)_0 + (u_1)_{x_0} * (op * curl)_2
1308 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1309 e) += (-gradXY20[dy][dx] * wDz) + (gradXY02[dy][dx] * wz);
1315 osc += D1Dx * D1Dy * D1Dz;
1320 const int D1Dz = D1D - 1;
1321 const int D1Dy = D1D;
1322 const int D1Dx = D1D;
1324 for (int qx = 0; qx < Q1D; ++qx)
1326 double gradYZ01[MAX_D1D][MAX_D1D];
1327 double gradYZ10[MAX_D1D][MAX_D1D];
1329 for (int dy = 0; dy < D1Dy; ++dy)
1331 for (int dz = 0; dz < D1Dz; ++dz)
1333 gradYZ01[dz][dy] = 0.0;
1334 gradYZ10[dz][dy] = 0.0;
1337 for (int qy = 0; qy < Q1D; ++qy)
1339 double massZ[MAX_D1D][2];
1340 for (int dz = 0; dz < D1Dz; ++dz)
1342 for (int n = 0; n < 2; ++n)
1347 for (int qz = 0; qz < Q1D; ++qz)
1349 for (int dz = 0; dz < D1Dz; ++dz)
1351 const double wz = Bot(dz,qz);
1353 massZ[dz][0] += wz * curl[qz][qy][qx][0];
1354 massZ[dz][1] += wz * curl[qz][qy][qx][1];
1357 for (int dy = 0; dy < D1Dy; ++dy)
1359 const double wy = Bct(dy,qy);
1360 const double wDy = Gct(dy,qy);
1362 for (int dz = 0; dz < D1Dz; ++dz)
1364 gradYZ01[dz][dy] += wy * massZ[dz][1];
1365 gradYZ10[dz][dy] += wDy * massZ[dz][0];
1370 for (int dx = 0; dx < D1Dx; ++dx)
1372 const double wx = Bct(dx,qx);
1373 const double wDx = Gct(dx,qx);
1375 for (int dy = 0; dy < D1Dy; ++dy)
1377 for (int dz = 0; dz < D1Dz; ++dz)
1379 // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1380 // (u_2)_{x_1} * (op * curl)_0 - (u_2)_{x_0} * (op * curl)_1
1381 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1382 e) += (gradYZ10[dz][dy] * wx) - (gradYZ01[dz][dy] * wDx);
1389 }); // end of element loop
1392 void CurlCurlIntegrator::AddMultPA(const Vector &x, Vector &y) const
1396 PACurlCurlApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1397 mapsC->Bt, mapsC->G, mapsC->Gt, pa_data, x, y);
1401 PACurlCurlApply2D(dofs1D, quad1D, ne, mapsO->B, mapsO->Bt,
1402 mapsC->G, mapsC->Gt, pa_data, x, y);
1406 MFEM_ABORT("Unsupported dimension!
");
1410 static void PACurlCurlAssembleDiagonal2D(const int D1D,
1413 const Array<double> &_Bo,
1414 const Array<double> &_Gc,
1418 constexpr static int VDIM = 2;
1420 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
1421 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1422 auto op = Reshape(_op.Read(), Q1D, Q1D, NE);
1423 auto diag = Reshape(_diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
1429 for (int c = 0; c < VDIM; ++c) // loop over x, y components
1431 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1432 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1436 for (int dy = 0; dy < D1Dy; ++dy)
1438 for (int qx = 0; qx < Q1D; ++qx)
1441 for (int qy = 0; qy < Q1D; ++qy)
1443 const double wy = (c == 1) ? Bo(qy,dy) : -Gc(qy,dy);
1444 t[qx] += wy * wy * op(qx,qy,e);
1448 for (int dx = 0; dx < D1Dx; ++dx)
1450 for (int qx = 0; qx < Q1D; ++qx)
1452 const double wx = ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
1453 diag(dx + (dy * D1Dx) + osc, e) += t[qx] * wx * wx;
1460 }); // end of element loop
1463 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
1464 static void PACurlCurlAssembleDiagonal3D(const int D1D,
1467 const Array<double> &_Bo,
1468 const Array<double> &_Bc,
1469 const Array<double> &_Go,
1470 const Array<double> &_Gc,
1474 constexpr static int VDIM = 3;
1475 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D >
MAX_D1D");
1476 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D >
MAX_Q1D");
1478 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
1479 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1480 auto Go = Reshape(_Go.Read(), Q1D, D1D-1);
1481 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1482 auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
1483 auto diag = Reshape(_diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
1487 // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
1488 // (\nabla\times u) \cdot (\nabla\times u) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{u}
1489 // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1490 // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1491 // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1493 // For each c, we will keep 6 arrays for derivatives multiplied by the 6 entries of the symmetric 3x3 matrix (dF^T dF).
1497 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
1499 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
1500 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1501 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1503 double zt[MAX_Q1D][MAX_Q1D][MAX_D1D][6][3];
1506 for (int qx = 0; qx < Q1D; ++qx)
1508 for (int qy = 0; qy < Q1D; ++qy)
1510 for (int dz = 0; dz < D1Dz; ++dz)
1512 for (int i=0; i<6; ++i)
1514 for (int d=0; d<3; ++d)
1516 zt[qx][qy][dz][i][d] = 0.0;
1520 for (int qz = 0; qz < Q1D; ++qz)
1522 const double wz = ((c == 2) ? Bo(qz,dz) : Bc(qz,dz));
1523 const double wDz = ((c == 2) ? Go(qz,dz) : Gc(qz,dz));
1525 for (int i=0; i<6; ++i)
1527 zt[qx][qy][dz][i][0] += wz * wz * op(qx,qy,qz,i,e);
1528 zt[qx][qy][dz][i][1] += wDz * wz * op(qx,qy,qz,i,e);
1529 zt[qx][qy][dz][i][2] += wDz * wDz * op(qx,qy,qz,i,e);
1534 } // end of z contraction
1536 double yt[MAX_Q1D][MAX_D1D][MAX_D1D][6][3][3];
1539 for (int qx = 0; qx < Q1D; ++qx)
1541 for (int dz = 0; dz < D1Dz; ++dz)
1543 for (int dy = 0; dy < D1Dy; ++dy)
1545 for (int i=0; i<6; ++i)
1547 for (int d=0; d<3; ++d)
1548 for (int j=0; j<3; ++j)
1550 yt[qx][dy][dz][i][d][j] = 0.0;
1554 for (int qy = 0; qy < Q1D; ++qy)
1556 const double wy = ((c == 1) ? Bo(qy,dy) : Bc(qy,dy));
1557 const double wDy = ((c == 1) ? Go(qy,dy) : Gc(qy,dy));
1559 for (int i=0; i<6; ++i)
1561 for (int d=0; d<3; ++d)
1563 yt[qx][dy][dz][i][d][0] += wy * wy * zt[qx][qy][dz][i][d];
1564 yt[qx][dy][dz][i][d][1] += wDy * wy * zt[qx][qy][dz][i][d];
1565 yt[qx][dy][dz][i][d][2] += wDy * wDy * zt[qx][qy][dz][i][d];
1571 } // end of y contraction
1574 for (int dz = 0; dz < D1Dz; ++dz)
1576 for (int dy = 0; dy < D1Dy; ++dy)
1578 for (int dx = 0; dx < D1Dx; ++dx)
1580 for (int qx = 0; qx < Q1D; ++qx)
1582 const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
1583 const double wDx = ((c == 0) ? Go(qx,dx) : Gc(qx,dx));
1585 // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
1586 // (\nabla\times u) \cdot (\nabla\times u) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{u}
1587 // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1588 // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1589 // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1592 const double O11 = op(q,0,e);
1593 const double O12 = op(q,1,e);
1594 const double O13 = op(q,2,e);
1595 const double O22 = op(q,3,e);
1596 const double O23 = op(q,4,e);
1597 const double O33 = op(q,5,e);
1602 // (u_0)_{x_2} (O22 (u_0)_{x_2} - O23 (u_0)_{x_1}) - (u_0)_{x_1} (O32 (u_0)_{x_2} - O33 (u_0)_{x_1})
1604 // (u_0)_{x_2} O22 (u_0)_{x_2}
1605 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1606 e) += yt[qx][dy][dz][3][2][0] * wx * wx;
1608 // -(u_0)_{x_2} O23 (u_0)_{x_1} - (u_0)_{x_1} O32 (u_0)_{x_2}
1609 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1610 e) += -2.0 * yt[qx][dy][dz][4][1][1] * wx * wx;
1612 // (u_0)_{x_1} O33 (u_0)_{x_1}
1613 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1614 e) += yt[qx][dy][dz][5][0][2] * wx * wx;
1618 // (u_1)_{x_2} (O11 (u_1)_{x_2} - O13 (u_1)_{x_0}) + (u_1)_{x_0} (-O31 (u_1)_{x_2} + O33 (u_1)_{x_0})
1620 // (u_1)_{x_2} O11 (u_1)_{x_2}
1621 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1622 e) += yt[qx][dy][dz][0][2][0] * wx * wx;
1624 // -(u_1)_{x_2} O13 (u_1)_{x_0} - (u_1)_{x_0} O31 (u_1)_{x_2}
1625 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1626 e) += -2.0 * yt[qx][dy][dz][2][1][0] * wDx * wx;
1628 // (u_1)_{x_0} O33 (u_1)_{x_0})
1629 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1630 e) += yt[qx][dy][dz][5][0][0] * wDx * wDx;
1634 // (u_2)_{x_1} (O11 (u_2)_{x_1} - O12 (u_2)_{x_0}) - (u_2)_{x_0} (O21 (u_2)_{x_1} - O22 (u_2)_{x_0})
1636 // (u_2)_{x_1} O11 (u_2)_{x_1}
1637 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1638 e) += yt[qx][dy][dz][0][0][2] * wx * wx;
1640 // -(u_2)_{x_1} O12 (u_2)_{x_0} - (u_2)_{x_0} O21 (u_2)_{x_1}
1641 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1642 e) += -2.0 * yt[qx][dy][dz][1][0][1] * wDx * wx;
1644 // (u_2)_{x_0} O22 (u_2)_{x_0}
1645 diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1646 e) += yt[qx][dy][dz][3][0][0] * wDx * wDx;
1651 } // end of x contraction
1653 osc += D1Dx * D1Dy * D1Dz;
1655 }); // end of element loop
1658 void CurlCurlIntegrator::AssembleDiagonalPA(Vector& diag)
1662 // Reduce HCURL_MAX_D1D/Q1D to avoid using too much memory
1663 constexpr int MAX_D1D = 4;
1664 constexpr int MAX_Q1D = 5;
1665 PACurlCurlAssembleDiagonal3D<MAX_D1D,MAX_Q1D>(dofs1D, quad1D, ne,
1672 PACurlCurlAssembleDiagonal2D(dofs1D, quad1D, ne,
1673 mapsO->B, mapsC->G, pa_data, diag);
1677 MFEM_ABORT("Unsupported dimension!
");
1681 void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace
1683 const FiniteElementSpace &test_fes)
1685 // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1686 Mesh *mesh = trial_fes.GetMesh();
1687 const FiniteElement *trial_fel = trial_fes.GetFE(0);
1688 const FiniteElement *test_fel = test_fes.GetFE(0);
1690 const NodalTensorFiniteElement *trial_el =
1691 dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1694 const VectorTensorFiniteElement *test_el =
1695 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1698 const IntegrationRule *ir
1699 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
1700 *mesh->GetElementTransformation(0));
1701 const int dims = trial_el->GetDim();
1702 MFEM_VERIFY(dims == 2 || dims == 3, "");
1704 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
1705 const int nq = ir->GetNPoints();
1706 dim = mesh->Dimension();
1707 MFEM_VERIFY(dim == 2 || dim == 3, "");
1709 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1711 ne = trial_fes.GetNE();
1712 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
1713 mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
1714 mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
1715 dofs1D = mapsC->ndof;
1716 quad1D = mapsC->nqpt;
1718 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
1720 pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
1722 Vector coeff(ne * nq);
1726 for (int e=0; e<ne; ++e)
1728 ElementTransformation *tr = mesh->GetElementTransformation(e);
1729 for (int p=0; p<nq; ++p)
1731 coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
1736 // Use the same setup functions as VectorFEMassIntegrator.
1737 if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
1739 PAHcurlSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
1742 else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
1744 PAHcurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
1749 MFEM_ABORT("Unknown kernel.
");
1753 // Apply to x corresponding to DOF's in H^1 (trial), whose gradients are integrated
1754 // against H(curl) test functions corresponding to y.
1755 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
1756 static void PAHcurlH1Apply3D(const int D1D,
1759 const Array<double> &_Bc,
1760 const Array<double> &_Gc,
1761 const Array<double> &_Bot,
1762 const Array<double> &_Bct,
1767 constexpr static int VDIM = 3;
1769 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1770 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1771 auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
1772 auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
1773 auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
1774 auto x = Reshape(_x.Read(), D1D, D1D, D1D, NE);
1775 auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
1779 double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
1781 for (int qz = 0; qz < Q1D; ++qz)
1783 for (int qy = 0; qy < Q1D; ++qy)
1785 for (int qx = 0; qx < Q1D; ++qx)
1787 for (int c = 0; c < VDIM; ++c)
1789 mass[qz][qy][qx][c] = 0.0;
1795 for (int dz = 0; dz < D1D; ++dz)
1797 double gradXY[MAX_Q1D][MAX_Q1D][3];
1798 for (int qy = 0; qy < Q1D; ++qy)
1800 for (int qx = 0; qx < Q1D; ++qx)
1802 gradXY[qy][qx][0] = 0.0;
1803 gradXY[qy][qx][1] = 0.0;
1804 gradXY[qy][qx][2] = 0.0;
1807 for (int dy = 0; dy < D1D; ++dy)
1809 double gradX[MAX_Q1D][2];
1810 for (int qx = 0; qx < Q1D; ++qx)
1815 for (int dx = 0; dx < D1D; ++dx)
1817 const double s = x(dx,dy,dz,e);
1818 for (int qx = 0; qx < Q1D; ++qx)
1820 gradX[qx][0] += s * Bc(qx,dx);
1821 gradX[qx][1] += s * Gc(qx,dx);
1824 for (int qy = 0; qy < Q1D; ++qy)
1826 const double wy = Bc(qy,dy);
1827 const double wDy = Gc(qy,dy);
1828 for (int qx = 0; qx < Q1D; ++qx)
1830 const double wx = gradX[qx][0];
1831 const double wDx = gradX[qx][1];
1832 gradXY[qy][qx][0] += wDx * wy;
1833 gradXY[qy][qx][1] += wx * wDy;
1834 gradXY[qy][qx][2] += wx * wy;
1838 for (int qz = 0; qz < Q1D; ++qz)
1840 const double wz = Bc(qz,dz);
1841 const double wDz = Gc(qz,dz);
1842 for (int qy = 0; qy < Q1D; ++qy)
1844 for (int qx = 0; qx < Q1D; ++qx)
1846 mass[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
1847 mass[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
1848 mass[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
1854 // Apply D operator.
1855 for (int qz = 0; qz < Q1D; ++qz)
1857 for (int qy = 0; qy < Q1D; ++qy)
1859 for (int qx = 0; qx < Q1D; ++qx)
1861 const double O11 = op(qx,qy,qz,0,e);
1862 const double O12 = op(qx,qy,qz,1,e);
1863 const double O13 = op(qx,qy,qz,2,e);
1864 const double O22 = op(qx,qy,qz,3,e);
1865 const double O23 = op(qx,qy,qz,4,e);
1866 const double O33 = op(qx,qy,qz,5,e);
1867 const double massX = mass[qz][qy][qx][0];
1868 const double massY = mass[qz][qy][qx][1];
1869 const double massZ = mass[qz][qy][qx][2];
1870 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
1871 mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ);
1872 mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ);
1877 for (int qz = 0; qz < Q1D; ++qz)
1879 double massXY[MAX_D1D][MAX_D1D];
1883 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
1885 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
1886 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1887 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1889 for (int dy = 0; dy < D1Dy; ++dy)
1891 for (int dx = 0; dx < D1Dx; ++dx)
1896 for (int qy = 0; qy < Q1D; ++qy)
1898 double massX[MAX_D1D];
1899 for (int dx = 0; dx < D1Dx; ++dx)
1903 for (int qx = 0; qx < Q1D; ++qx)
1905 for (int dx = 0; dx < D1Dx; ++dx)
1907 massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
1910 for (int dy = 0; dy < D1Dy; ++dy)
1912 const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
1913 for (int dx = 0; dx < D1Dx; ++dx)
1915 massXY[dy][dx] += massX[dx] * wy;
1920 for (int dz = 0; dz < D1Dz; ++dz)
1922 const double wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
1923 for (int dy = 0; dy < D1Dy; ++dy)
1925 for (int dx = 0; dx < D1Dx; ++dx)
1927 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
1932 osc += D1Dx * D1Dy * D1Dz;
1935 }); // end of element loop
1938 // Apply to x corresponding to DOF's in H^1 (trial), whose gradients are integrated
1939 // against H(curl) test functions corresponding to y.
1940 static void PAHcurlH1Apply2D(const int D1D,
1943 const Array<double> &_Bc,
1944 const Array<double> &_Gc,
1945 const Array<double> &_Bot,
1946 const Array<double> &_Bct,
1951 constexpr static int VDIM = 2;
1953 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1954 auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1955 auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
1956 auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
1957 auto op = Reshape(_op.Read(), Q1D, Q1D, 3, NE);
1958 auto x = Reshape(_x.Read(), D1D, D1D, NE);
1959 auto y = Reshape(_y.ReadWrite(), 2*(D1D-1)*D1D, NE);
1963 double mass[MAX_Q1D][MAX_Q1D][VDIM];
1965 for (int qy = 0; qy < Q1D; ++qy)
1967 for (int qx = 0; qx < Q1D; ++qx)
1969 for (int c = 0; c < VDIM; ++c)
1971 mass[qy][qx][c] = 0.0;
1976 for (int dy = 0; dy < D1D; ++dy)
1978 double gradX[MAX_Q1D][2];
1979 for (int qx = 0; qx < Q1D; ++qx)
1984 for (int dx = 0; dx < D1D; ++dx)
1986 const double s = x(dx,dy,e);
1987 for (int qx = 0; qx < Q1D; ++qx)
1989 gradX[qx][0] += s * Bc(qx,dx);
1990 gradX[qx][1] += s * Gc(qx,dx);
1993 for (int qy = 0; qy < Q1D; ++qy)
1995 const double wy = Bc(qy,dy);
1996 const double wDy = Gc(qy,dy);
1997 for (int qx = 0; qx < Q1D; ++qx)
1999 const double wx = gradX[qx][0];
2000 const double wDx = gradX[qx][1];
2001 mass[qy][qx][0] += wDx * wy;
2002 mass[qy][qx][1] += wx * wDy;
2007 // Apply D operator.
2008 for (int qy = 0; qy < Q1D; ++qy)
2010 for (int qx = 0; qx < Q1D; ++qx)
2012 const double O11 = op(qx,qy,0,e);
2013 const double O12 = op(qx,qy,1,e);
2014 const double O22 = op(qx,qy,2,e);
2015 const double massX = mass[qy][qx][0];
2016 const double massY = mass[qy][qx][1];
2017 mass[qy][qx][0] = (O11*massX)+(O12*massY);
2018 mass[qy][qx][1] = (O12*massX)+(O22*massY);
2022 for (int qy = 0; qy < Q1D; ++qy)
2026 for (int c = 0; c < VDIM; ++c) // loop over x, y components
2028 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
2029 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
2031 double massX[MAX_D1D];
2032 for (int dx = 0; dx < D1Dx; ++dx)
2036 for (int qx = 0; qx < Q1D; ++qx)
2038 for (int dx = 0; dx < D1Dx; ++dx)
2040 massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
2044 for (int dy = 0; dy < D1Dy; ++dy)
2046 const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
2048 for (int dx = 0; dx < D1Dx; ++dx)
2050 y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
2057 }); // end of element loop
2060 void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const
2063 PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
2064 mapsO->Bt, mapsC->Bt, pa_data, x, y);
2066 PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
2067 mapsO->Bt, mapsC->Bt, pa_data, x, y);
2070 MFEM_ABORT("Unsupported dimension!
");
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for Finite Elements.
int GetDim() const
Returns the reference space dimension for the finite element.
Class for an integration rule - an Array of IntegrationPoint.
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
int GetNE() const
Returns number of elements in the mesh.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Implements CalcCurlShape methods.
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags)
Return the mesh geometric factors corresponding to the given integration rule.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Mesh * GetMesh() const
Returns the mesh.
constexpr int HCURL_MAX_D1D
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i'th element.
constexpr int HCURL_MAX_Q1D