// Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced // at the Lawrence Livermore National Laboratory. All Rights reserved. See files // LICENSE and NOTICE for details. LLNL-CODE-806117. // // This file is part of the MFEM library. For more information and source code // availability visit https://mfem.org. // // MFEM is free software; you can redistribute it and/or modify it under the // terms of the BSD-3 license. We welcome feedback and contributions, see file // CONTRIBUTING.md for details. #include "../../general/forall.hpp" #include "../bilininteg.hpp" #include "../gridfunc.hpp" #include "../qfunction.hpp" #include "bilininteg_diffusion_kernels.hpp" namespace mfem { // Apply to x corresponding to DOFs in H^1 (trial), whose gradients are // integrated against H(curl) test functions corresponding to y. static void PAHcurlH1Apply2D(const int D1D, const int Q1D, const int NE, const Array &bc, const Array &gc, const Array &bot, const Array &bct, const Vector &pa_data, const Vector &x, Vector &y) { auto Bc = Reshape(bc.Read(), Q1D, D1D); auto Gc = Reshape(gc.Read(), Q1D, D1D); auto Bot = Reshape(bot.Read(), D1D-1, Q1D); auto Bct = Reshape(bct.Read(), D1D, Q1D); auto op = Reshape(pa_data.Read(), Q1D, Q1D, 3, NE); auto X = Reshape(x.Read(), D1D, D1D, NE); auto Y = Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE); mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e) { constexpr static int VDIM = 2; constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D; constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D; double mass[MAX_Q1D][MAX_Q1D][VDIM]; for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { for (int c = 0; c < VDIM; ++c) { mass[qy][qx][c] = 0.0; } } } for (int dy = 0; dy < D1D; ++dy) { double gradX[MAX_Q1D][2]; for (int qx = 0; qx < Q1D; ++qx) { gradX[qx][0] = 0.0; gradX[qx][1] = 0.0; } for (int dx = 0; dx < D1D; ++dx) { const double s = X(dx,dy,e); for (int qx = 0; qx < Q1D; ++qx) { gradX[qx][0] += s * Bc(qx,dx); gradX[qx][1] += s * Gc(qx,dx); } } for (int qy = 0; qy < Q1D; ++qy) { const double wy = Bc(qy,dy); const double wDy = Gc(qy,dy); for (int qx = 0; qx < Q1D; ++qx) { const double wx = gradX[qx][0]; const double wDx = gradX[qx][1]; mass[qy][qx][0] += wDx * wy; mass[qy][qx][1] += wx * wDy; } } } // Apply D operator. for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { const double O11 = op(qx,qy,0,e); const double O12 = op(qx,qy,1,e); const double O22 = op(qx,qy,2,e); const double massX = mass[qy][qx][0]; const double massY = mass[qy][qx][1]; mass[qy][qx][0] = (O11*massX)+(O12*massY); mass[qy][qx][1] = (O12*massX)+(O22*massY); } } for (int qy = 0; qy < Q1D; ++qy) { int osc = 0; for (int c = 0; c < VDIM; ++c) // loop over x, y components { const int D1Dy = (c == 1) ? D1D - 1 : D1D; const int D1Dx = (c == 0) ? D1D - 1 : D1D; double massX[MAX_D1D]; for (int dx = 0; dx < D1Dx; ++dx) { massX[dx] = 0; } for (int qx = 0; qx < Q1D; ++qx) { for (int dx = 0; dx < D1Dx; ++dx) { massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)); } } for (int dy = 0; dy < D1Dy; ++dy) { const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy); for (int dx = 0; dx < D1Dx; ++dx) { Y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy; } } osc += D1Dx * D1Dy; } // loop c } }); // end of element loop } // Apply to x corresponding to DOFs in H(curl), integrated // against gradients of H^1 functions corresponding to y. static void PAHcurlH1ApplyTranspose2D(const int D1D, const int Q1D, const int NE, const Array &bc, const Array &bo, const Array &bct, const Array &gct, const Vector &pa_data, const Vector &x, Vector &y) { auto Bc = Reshape(bc.Read(), Q1D, D1D); auto Bo = Reshape(bo.Read(), Q1D, D1D-1); auto Bt = Reshape(bct.Read(), D1D, Q1D); auto Gt = Reshape(gct.Read(), D1D, Q1D); auto op = Reshape(pa_data.Read(), Q1D, Q1D, 3, NE); auto X = Reshape(x.Read(), 2*(D1D-1)*D1D, NE); auto Y = Reshape(y.ReadWrite(), D1D, D1D, NE); mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e) { constexpr static int VDIM = 2; constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D; constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D; double mass[MAX_Q1D][MAX_Q1D][VDIM]; for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { for (int c = 0; c < VDIM; ++c) { mass[qy][qx][c] = 0.0; } } } int osc = 0; for (int c = 0; c < VDIM; ++c) // loop over x, y components { const int D1Dy = (c == 1) ? D1D - 1 : D1D; const int D1Dx = (c == 0) ? D1D - 1 : D1D; for (int dy = 0; dy < D1Dy; ++dy) { double massX[MAX_Q1D]; for (int qx = 0; qx < Q1D; ++qx) { massX[qx] = 0.0; } for (int dx = 0; dx < D1Dx; ++dx) { const double t = X(dx + (dy * D1Dx) + osc, e); for (int qx = 0; qx < Q1D; ++qx) { massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)); } } for (int qy = 0; qy < Q1D; ++qy) { const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy); for (int qx = 0; qx < Q1D; ++qx) { mass[qy][qx][c] += massX[qx] * wy; } } } osc += D1Dx * D1Dy; } // loop (c) over components // Apply D operator. for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { const double O11 = op(qx,qy,0,e); const double O12 = op(qx,qy,1,e); const double O22 = op(qx,qy,2,e); const double massX = mass[qy][qx][0]; const double massY = mass[qy][qx][1]; mass[qy][qx][0] = (O11*massX)+(O12*massY); mass[qy][qx][1] = (O12*massX)+(O22*massY); } } for (int qy = 0; qy < Q1D; ++qy) { double gradX[MAX_D1D][2]; for (int dx = 0; dx < D1D; ++dx) { gradX[dx][0] = 0; gradX[dx][1] = 0; } for (int qx = 0; qx < Q1D; ++qx) { const double gX = mass[qy][qx][0]; const double gY = mass[qy][qx][1]; for (int dx = 0; dx < D1D; ++dx) { const double wx = Bt(dx,qx); const double wDx = Gt(dx,qx); gradX[dx][0] += gX * wDx; gradX[dx][1] += gY * wx; } } for (int dy = 0; dy < D1D; ++dy) { const double wy = Bt(dy,qy); const double wDy = Gt(dy,qy); for (int dx = 0; dx < D1D; ++dx) { Y(dx,dy,e) += ((gradX[dx][0] * wy) + (gradX[dx][1] * wDy)); } } } }); // end of element loop } // Apply to x corresponding to DOFs in H^1 (trial), whose gradients are // integrated against H(curl) test functions corresponding to y. static void PAHcurlH1Apply3D(const int D1D, const int Q1D, const int NE, const Array &bc, const Array &gc, const Array &bot, const Array &bct, const Vector &pa_data, const Vector &x, Vector &y) { MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D, "Error: D1D > MAX_D1D"); MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D, "Error: Q1D > MAX_Q1D"); constexpr static int VDIM = 3; auto Bc = Reshape(bc.Read(), Q1D, D1D); auto Gc = Reshape(gc.Read(), Q1D, D1D); auto Bot = Reshape(bot.Read(), D1D-1, Q1D); auto Bct = Reshape(bct.Read(), D1D, Q1D); auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, 6, NE); auto X = Reshape(x.Read(), D1D, D1D, D1D, NE); auto Y = Reshape(y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE); mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e) { constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D; constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D; double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM]; for (int qz = 0; qz < Q1D; ++qz) { for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { for (int c = 0; c < VDIM; ++c) { mass[qz][qy][qx][c] = 0.0; } } } } for (int dz = 0; dz < D1D; ++dz) { double gradXY[MAX_Q1D][MAX_Q1D][3]; for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { gradXY[qy][qx][0] = 0.0; gradXY[qy][qx][1] = 0.0; gradXY[qy][qx][2] = 0.0; } } for (int dy = 0; dy < D1D; ++dy) { double gradX[MAX_Q1D][2]; for (int qx = 0; qx < Q1D; ++qx) { gradX[qx][0] = 0.0; gradX[qx][1] = 0.0; } for (int dx = 0; dx < D1D; ++dx) { const double s = X(dx,dy,dz,e); for (int qx = 0; qx < Q1D; ++qx) { gradX[qx][0] += s * Bc(qx,dx); gradX[qx][1] += s * Gc(qx,dx); } } for (int qy = 0; qy < Q1D; ++qy) { const double wy = Bc(qy,dy); const double wDy = Gc(qy,dy); for (int qx = 0; qx < Q1D; ++qx) { const double wx = gradX[qx][0]; const double wDx = gradX[qx][1]; gradXY[qy][qx][0] += wDx * wy; gradXY[qy][qx][1] += wx * wDy; gradXY[qy][qx][2] += wx * wy; } } } for (int qz = 0; qz < Q1D; ++qz) { const double wz = Bc(qz,dz); const double wDz = Gc(qz,dz); for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { mass[qz][qy][qx][0] += gradXY[qy][qx][0] * wz; mass[qz][qy][qx][1] += gradXY[qy][qx][1] * wz; mass[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz; } } } } // Apply D operator. for (int qz = 0; qz < Q1D; ++qz) { for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { const double O11 = op(qx,qy,qz,0,e); const double O12 = op(qx,qy,qz,1,e); const double O13 = op(qx,qy,qz,2,e); const double O22 = op(qx,qy,qz,3,e); const double O23 = op(qx,qy,qz,4,e); const double O33 = op(qx,qy,qz,5,e); const double massX = mass[qz][qy][qx][0]; const double massY = mass[qz][qy][qx][1]; const double massZ = mass[qz][qy][qx][2]; mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ); mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ); mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ); } } } for (int qz = 0; qz < Q1D; ++qz) { double massXY[MAX_D1D][MAX_D1D]; int osc = 0; for (int c = 0; c < VDIM; ++c) // loop over x, y, z components { const int D1Dz = (c == 2) ? D1D - 1 : D1D; const int D1Dy = (c == 1) ? D1D - 1 : D1D; const int D1Dx = (c == 0) ? D1D - 1 : D1D; for (int dy = 0; dy < D1Dy; ++dy) { for (int dx = 0; dx < D1Dx; ++dx) { massXY[dy][dx] = 0.0; } } for (int qy = 0; qy < Q1D; ++qy) { double massX[MAX_D1D]; for (int dx = 0; dx < D1Dx; ++dx) { massX[dx] = 0; } for (int qx = 0; qx < Q1D; ++qx) { for (int dx = 0; dx < D1Dx; ++dx) { massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)); } } for (int dy = 0; dy < D1Dy; ++dy) { const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy); for (int dx = 0; dx < D1Dx; ++dx) { massXY[dy][dx] += massX[dx] * wy; } } } for (int dz = 0; dz < D1Dz; ++dz) { const double wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz); for (int dy = 0; dy < D1Dy; ++dy) { for (int dx = 0; dx < D1Dx; ++dx) { Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz; } } } osc += D1Dx * D1Dy * D1Dz; } // loop c } // loop qz }); // end of element loop } // Apply to x corresponding to DOFs in H(curl), integrated // against gradients of H^1 functions corresponding to y. static void PAHcurlH1ApplyTranspose3D(const int D1D, const int Q1D, const int NE, const Array &bc, const Array &bo, const Array &bct, const Array &gct, const Vector &pa_data, const Vector &x, Vector &y) { MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D, "Error: D1D > MAX_D1D"); MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D, "Error: Q1D > MAX_Q1D"); constexpr static int VDIM = 3; auto Bc = Reshape(bc.Read(), Q1D, D1D); auto Bo = Reshape(bo.Read(), Q1D, D1D-1); auto Bt = Reshape(bct.Read(), D1D, Q1D); auto Gt = Reshape(gct.Read(), D1D, Q1D); auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, 6, NE); auto X = Reshape(x.Read(), 3*(D1D-1)*D1D*D1D, NE); auto Y = Reshape(y.ReadWrite(), D1D, D1D, D1D, NE); mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e) { constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D; constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D; double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM]; for (int qz = 0; qz < Q1D; ++qz) { for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { for (int c = 0; c < VDIM; ++c) { mass[qz][qy][qx][c] = 0.0; } } } } int osc = 0; for (int c = 0; c < VDIM; ++c) // loop over x, y, z components { const int D1Dz = (c == 2) ? D1D - 1 : D1D; const int D1Dy = (c == 1) ? D1D - 1 : D1D; const int D1Dx = (c == 0) ? D1D - 1 : D1D; for (int dz = 0; dz < D1Dz; ++dz) { double massXY[MAX_Q1D][MAX_Q1D]; for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { massXY[qy][qx] = 0.0; } } for (int dy = 0; dy < D1Dy; ++dy) { double massX[MAX_Q1D]; for (int qx = 0; qx < Q1D; ++qx) { massX[qx] = 0.0; } for (int dx = 0; dx < D1Dx; ++dx) { const double t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e); for (int qx = 0; qx < Q1D; ++qx) { massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)); } } for (int qy = 0; qy < Q1D; ++qy) { const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy); for (int qx = 0; qx < Q1D; ++qx) { const double wx = massX[qx]; massXY[qy][qx] += wx * wy; } } } for (int qz = 0; qz < Q1D; ++qz) { const double wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz); for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { mass[qz][qy][qx][c] += massXY[qy][qx] * wz; } } } } osc += D1Dx * D1Dy * D1Dz; } // loop (c) over components // Apply D operator. for (int qz = 0; qz < Q1D; ++qz) { for (int qy = 0; qy < Q1D; ++qy) { for (int qx = 0; qx < Q1D; ++qx) { const double O11 = op(qx,qy,qz,0,e); const double O12 = op(qx,qy,qz,1,e); const double O13 = op(qx,qy,qz,2,e); const double O22 = op(qx,qy,qz,3,e); const double O23 = op(qx,qy,qz,4,e); const double O33 = op(qx,qy,qz,5,e); const double massX = mass[qz][qy][qx][0]; const double massY = mass[qz][qy][qx][1]; const double massZ = mass[qz][qy][qx][2]; mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ); mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ); mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ); } } } for (int qz = 0; qz < Q1D; ++qz) { double gradXY[MAX_D1D][MAX_D1D][3]; for (int dy = 0; dy < D1D; ++dy) { for (int dx = 0; dx < D1D; ++dx) { gradXY[dy][dx][0] = 0; gradXY[dy][dx][1] = 0; gradXY[dy][dx][2] = 0; } } for (int qy = 0; qy < Q1D; ++qy) { double gradX[MAX_D1D][3]; for (int dx = 0; dx < D1D; ++dx) { gradX[dx][0] = 0; gradX[dx][1] = 0; gradX[dx][2] = 0; } for (int qx = 0; qx < Q1D; ++qx) { const double gX = mass[qz][qy][qx][0]; const double gY = mass[qz][qy][qx][1]; const double gZ = mass[qz][qy][qx][2]; for (int dx = 0; dx < D1D; ++dx) { const double wx = Bt(dx,qx); const double wDx = Gt(dx,qx); gradX[dx][0] += gX * wDx; gradX[dx][1] += gY * wx; gradX[dx][2] += gZ * wx; } } for (int dy = 0; dy < D1D; ++dy) { const double wy = Bt(dy,qy); const double wDy = Gt(dy,qy); for (int dx = 0; dx < D1D; ++dx) { gradXY[dy][dx][0] += gradX[dx][0] * wy; gradXY[dy][dx][1] += gradX[dx][1] * wDy; gradXY[dy][dx][2] += gradX[dx][2] * wy; } } } for (int dz = 0; dz < D1D; ++dz) { const double wz = Bt(dz,qz); const double wDz = Gt(dz,qz); for (int dy = 0; dy < D1D; ++dy) { for (int dx = 0; dx < D1D; ++dx) { Y(dx,dy,dz,e) += ((gradXY[dy][dx][0] * wz) + (gradXY[dy][dx][1] * wz) + (gradXY[dy][dx][2] * wDz)); } } } } // loop qz }); // end of element loop } void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) { // Assumes tensor-product elements, with a vector test space and H^1 trial space. Mesh *mesh = trial_fes.GetMesh(); const FiniteElement *trial_fel = trial_fes.GetFE(0); const FiniteElement *test_fel = test_fes.GetFE(0); const NodalTensorFiniteElement *trial_el = dynamic_cast(trial_fel); MFEM_VERIFY(trial_el != NULL, "Only NodalTensorFiniteElement is supported!"); const VectorTensorFiniteElement *test_el = dynamic_cast(test_fel); MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!"); const IntegrationRule *ir = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el, *mesh->GetElementTransformation(0)); const int dims = trial_el->GetDim(); MFEM_VERIFY(dims == 2 || dims == 3, ""); const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6 const int nq = ir->GetNPoints(); dim = mesh->Dimension(); MFEM_VERIFY(dim == 2 || dim == 3, ""); MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), ""); ne = trial_fes.GetNE(); geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS); mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR); mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR); dofs1D = mapsC->ndof; quad1D = mapsC->nqpt; MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, ""); pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType()); QuadratureSpace qs(*mesh, *ir); CoefficientVector coeff(Q, qs, CoefficientStorage::FULL); // Use the same setup functions as VectorFEMassIntegrator. if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3) { internal::PADiffusionSetup3D(quad1D, 1, ne, ir->GetWeights(), geom->J, coeff, pa_data); } else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2) { internal::PADiffusionSetup2D<2>(quad1D, 1, ne, ir->GetWeights(), geom->J, coeff, pa_data); } else { MFEM_ABORT("Unknown kernel."); } } void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const { if (dim == 3) { PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G, mapsO->Bt, mapsC->Bt, pa_data, x, y); } else if (dim == 2) { PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G, mapsO->Bt, mapsC->Bt, pa_data, x, y); } else { MFEM_ABORT("Unsupported dimension!"); } } void MixedVectorGradientIntegrator::AddMultTransposePA(const Vector &x, Vector &y) const { if (dim == 3) { PAHcurlH1ApplyTranspose3D(dofs1D, quad1D, ne, mapsC->B, mapsO->B, mapsC->Bt, mapsC->Gt, pa_data, x, y); } else if (dim == 2) { PAHcurlH1ApplyTranspose2D(dofs1D, quad1D, ne, mapsC->B, mapsO->B, mapsC->Bt, mapsC->Gt, pa_data, x, y); } else { MFEM_ABORT("Unsupported dimension!"); } } } // namespace mfem