1182fbe45STzanio // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2182fbe45STzanio // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3182fbe45STzanio // reserved. See files LICENSE and NOTICE for details. 4182fbe45STzanio // 5182fbe45STzanio // This file is part of CEED, a collection of benchmarks, miniapps, software 6182fbe45STzanio // libraries and APIs for efficient high-order finite element and spectral 7182fbe45STzanio // element discretizations for exascale applications. For more information and 8182fbe45STzanio // source code availability see http://github.com/ceed. 9182fbe45STzanio // 10182fbe45STzanio // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11182fbe45STzanio // a collaborative effort of two U.S. Department of Energy organizations (Office 12182fbe45STzanio // of Science and the National Nuclear Security Administration) responsible for 13182fbe45STzanio // the planning and preparation of a capable exascale ecosystem, including 14182fbe45STzanio // software, applications, hardware, advanced system engineering and early 15182fbe45STzanio // testbed platforms, in support of the nation's exascale computing imperative. 16182fbe45STzanio 17182fbe45STzanio /// @file 18182fbe45STzanio /// MFEM diffusion operator based on libCEED 19182fbe45STzanio 20182fbe45STzanio #include <ceed.h> 21182fbe45STzanio #include <mfem.hpp> 22182fbe45STzanio 23182fbe45STzanio /// A structure used to pass additional data to f_build_diff and f_apply_diff 24a48d94bfSjeremylt struct BuildContext { CeedInt dim, space_dim; }; 25182fbe45STzanio 26182fbe45STzanio /// libCEED Q-function for building quadrature data for a diffusion operator 2754251743Sjeremylt static int f_build_diff(void *ctx, CeedInt Q, 2854251743Sjeremylt const CeedScalar *const *in, CeedScalar *const *out) { 29a48d94bfSjeremylt BuildContext *bc = (BuildContext*)ctx; 30*ecf6354eSJed Brown // in[0] is Jacobians with shape [dim, nc=dim, Q] 3154251743Sjeremylt // in[1] is quadrature weights, size (Q) 32182fbe45STzanio // 33182fbe45STzanio // At every quadrature point, compute qw/det(J).adj(J).adj(J)^T and store 34182fbe45STzanio // the symmetric part of the result. 357ca8db16Sjeremylt const CeedScalar *J = in[0], *qw = in[1]; 367ca8db16Sjeremylt CeedScalar *qd = out[0]; 37a48d94bfSjeremylt switch (bc->dim + 10*bc->space_dim) { 38182fbe45STzanio case 11: 39182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 40182fbe45STzanio qd[i] = qw[i] / J[i]; 41182fbe45STzanio } 42182fbe45STzanio break; 43182fbe45STzanio case 22: 44182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 45182fbe45STzanio // J: 0 2 qd: 0 1 adj(J): J22 -J12 46182fbe45STzanio // 1 3 1 2 -J21 J11 47182fbe45STzanio const CeedScalar J11 = J[i+Q*0]; 48182fbe45STzanio const CeedScalar J21 = J[i+Q*1]; 49182fbe45STzanio const CeedScalar J12 = J[i+Q*2]; 50182fbe45STzanio const CeedScalar J22 = J[i+Q*3]; 51182fbe45STzanio const CeedScalar w = qw[i] / (J11*J22 - J21*J12); 52182fbe45STzanio qd[i+Q*0] = w * (J12*J12 + J22*J22); 53182fbe45STzanio qd[i+Q*1] = - w * (J11*J12 + J21*J22); 54182fbe45STzanio qd[i+Q*2] = w * (J11*J11 + J21*J21); 55182fbe45STzanio } 56182fbe45STzanio break; 57182fbe45STzanio case 33: 58182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 59182fbe45STzanio // J: 0 3 6 qd: 0 1 2 60182fbe45STzanio // 1 4 7 1 3 4 61182fbe45STzanio // 2 5 8 2 4 5 62182fbe45STzanio const CeedScalar J11 = J[i+Q*0]; 63182fbe45STzanio const CeedScalar J21 = J[i+Q*1]; 64182fbe45STzanio const CeedScalar J31 = J[i+Q*2]; 65182fbe45STzanio const CeedScalar J12 = J[i+Q*3]; 66182fbe45STzanio const CeedScalar J22 = J[i+Q*4]; 67182fbe45STzanio const CeedScalar J32 = J[i+Q*5]; 68182fbe45STzanio const CeedScalar J13 = J[i+Q*6]; 69182fbe45STzanio const CeedScalar J23 = J[i+Q*7]; 70182fbe45STzanio const CeedScalar J33 = J[i+Q*8]; 71182fbe45STzanio const CeedScalar A11 = J22*J33 - J23*J32; 72182fbe45STzanio const CeedScalar A12 = J13*J32 - J12*J33; 73182fbe45STzanio const CeedScalar A13 = J12*J23 - J13*J22; 74182fbe45STzanio const CeedScalar A21 = J23*J31 - J21*J33; 75182fbe45STzanio const CeedScalar A22 = J11*J33 - J13*J31; 76182fbe45STzanio const CeedScalar A23 = J13*J21 - J11*J23; 77182fbe45STzanio const CeedScalar A31 = J21*J32 - J22*J31; 78182fbe45STzanio const CeedScalar A32 = J12*J31 - J11*J32; 79182fbe45STzanio const CeedScalar A33 = J11*J22 - J12*J21; 80182fbe45STzanio const CeedScalar w = qw[i] / (J11*A11 + J21*A12 + J31*A13); 81182fbe45STzanio qd[i+Q*0] = w * (A11*A11 + A12*A12 + A13*A13); 82182fbe45STzanio qd[i+Q*1] = w * (A11*A21 + A12*A22 + A13*A23); 83182fbe45STzanio qd[i+Q*2] = w * (A11*A31 + A12*A32 + A13*A33); 84182fbe45STzanio qd[i+Q*3] = w * (A21*A21 + A22*A22 + A23*A23); 85182fbe45STzanio qd[i+Q*4] = w * (A21*A31 + A22*A32 + A23*A33); 86182fbe45STzanio qd[i+Q*5] = w * (A31*A31 + A32*A32 + A33*A33); 87182fbe45STzanio } 88182fbe45STzanio break; 89182fbe45STzanio default: 90182fbe45STzanio return CeedError(NULL, 1, "dim=%d, space_dim=%d is not supported", 91a48d94bfSjeremylt bc->dim, bc->space_dim); 92182fbe45STzanio } 93182fbe45STzanio return 0; 94182fbe45STzanio } 95182fbe45STzanio 96182fbe45STzanio /// libCEED Q-function for applying a diff operator 9754251743Sjeremylt static int f_apply_diff(void *ctx, CeedInt Q, 9854251743Sjeremylt const CeedScalar *const *in, CeedScalar *const *out) { 99a48d94bfSjeremylt BuildContext *bc = (BuildContext*)ctx; 100*ecf6354eSJed Brown // in[0], out[0] have shape [dim, nc=1, Q] 1017ca8db16Sjeremylt const CeedScalar *ug = in[0], *qd = in[1]; 1027ca8db16Sjeremylt CeedScalar *vg = out[0]; 103a48d94bfSjeremylt switch (bc->dim) { 104182fbe45STzanio case 1: 105182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 106182fbe45STzanio vg[i] = ug[i] * qd[i]; 107182fbe45STzanio } 108182fbe45STzanio break; 109182fbe45STzanio case 2: 110182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 111182fbe45STzanio const CeedScalar ug0 = ug[i+Q*0]; 112182fbe45STzanio const CeedScalar ug1 = ug[i+Q*1]; 113182fbe45STzanio vg[i+Q*0] = qd[i+Q*0]*ug0 + qd[i+Q*1]*ug1; 114182fbe45STzanio vg[i+Q*1] = qd[i+Q*1]*ug0 + qd[i+Q*2]*ug1; 115182fbe45STzanio } 116182fbe45STzanio break; 117182fbe45STzanio case 3: 118182fbe45STzanio for (CeedInt i=0; i<Q; i++) { 119182fbe45STzanio const CeedScalar ug0 = ug[i+Q*0]; 120182fbe45STzanio const CeedScalar ug1 = ug[i+Q*1]; 121182fbe45STzanio const CeedScalar ug2 = ug[i+Q*2]; 122182fbe45STzanio vg[i+Q*0] = qd[i+Q*0]*ug0 + qd[i+Q*1]*ug1 + qd[i+Q*2]*ug2; 123182fbe45STzanio vg[i+Q*1] = qd[i+Q*1]*ug0 + qd[i+Q*3]*ug1 + qd[i+Q*4]*ug2; 124182fbe45STzanio vg[i+Q*2] = qd[i+Q*2]*ug0 + qd[i+Q*4]*ug1 + qd[i+Q*5]*ug2; 125182fbe45STzanio } 126182fbe45STzanio break; 127182fbe45STzanio default: 128a48d94bfSjeremylt return CeedError(NULL, 1, "topo_dim=%d is not supported", bc->dim); 129182fbe45STzanio } 130182fbe45STzanio return 0; 131182fbe45STzanio } 132182fbe45STzanio 133182fbe45STzanio /// Wrapper for a diffusion CeedOperator as an mfem::Operator 134182fbe45STzanio class CeedDiffusionOperator : public mfem::Operator { 135182fbe45STzanio protected: 136182fbe45STzanio const mfem::FiniteElementSpace *fes; 137182fbe45STzanio CeedOperator build_oper, oper; 138182fbe45STzanio CeedBasis basis, mesh_basis; 139135a076eSjeremylt CeedElemRestriction restr, mesh_restr, restr_i, mesh_restr_i; 140182fbe45STzanio CeedQFunction apply_qfunc, build_qfunc; 1417ca8db16Sjeremylt CeedVector node_coords, rho; 142182fbe45STzanio 143a48d94bfSjeremylt BuildContext build_ctx; 144182fbe45STzanio 145182fbe45STzanio CeedVector u, v; 146182fbe45STzanio 147182fbe45STzanio static void FESpace2Ceed(const mfem::FiniteElementSpace *fes, 148182fbe45STzanio const mfem::IntegrationRule &ir, 149182fbe45STzanio Ceed ceed, CeedBasis *basis, 150182fbe45STzanio CeedElemRestriction *restr) { 151182fbe45STzanio mfem::Mesh *mesh = fes->GetMesh(); 152182fbe45STzanio const mfem::FiniteElement *fe = fes->GetFE(0); 153182fbe45STzanio const int order = fes->GetOrder(0); 154182fbe45STzanio mfem::Array<int> dof_map; 155182fbe45STzanio switch (mesh->Dimension()) { 156182fbe45STzanio case 1: { 157182fbe45STzanio const mfem::H1_SegmentElement *h1_fe = 158182fbe45STzanio dynamic_cast<const mfem::H1_SegmentElement*>(fe); 159182fbe45STzanio MFEM_VERIFY(h1_fe, "invalid FE"); 160182fbe45STzanio h1_fe->GetDofMap().Copy(dof_map); 161182fbe45STzanio break; 162182fbe45STzanio } 163182fbe45STzanio case 2: { 164182fbe45STzanio const mfem::H1_QuadrilateralElement *h1_fe = 165182fbe45STzanio dynamic_cast<const mfem::H1_QuadrilateralElement*>(fe); 166182fbe45STzanio MFEM_VERIFY(h1_fe, "invalid FE"); 167182fbe45STzanio h1_fe->GetDofMap().Copy(dof_map); 168182fbe45STzanio break; 169182fbe45STzanio } 170182fbe45STzanio case 3: { 171182fbe45STzanio const mfem::H1_HexahedronElement *h1_fe = 172182fbe45STzanio dynamic_cast<const mfem::H1_HexahedronElement*>(fe); 173182fbe45STzanio MFEM_VERIFY(h1_fe, "invalid FE"); 174182fbe45STzanio h1_fe->GetDofMap().Copy(dof_map); 175182fbe45STzanio break; 176182fbe45STzanio } 177182fbe45STzanio } 178182fbe45STzanio const mfem::FiniteElement *fe1d = 179182fbe45STzanio fes->FEColl()->FiniteElementForGeometry(mfem::Geometry::SEGMENT); 180182fbe45STzanio mfem::DenseMatrix shape1d(fe1d->GetDof(), ir.GetNPoints()); 181182fbe45STzanio mfem::DenseMatrix grad1d(fe1d->GetDof(), ir.GetNPoints()); 182182fbe45STzanio mfem::Vector qref1d(ir.GetNPoints()), qweight1d(ir.GetNPoints()); 183182fbe45STzanio mfem::Vector shape_i(shape1d.Height()); 184182fbe45STzanio mfem::DenseMatrix grad_i(grad1d.Height(), 1); 185182fbe45STzanio const mfem::H1_SegmentElement *h1_fe1d = 186182fbe45STzanio dynamic_cast<const mfem::H1_SegmentElement*>(fe1d); 187182fbe45STzanio MFEM_VERIFY(h1_fe1d, "invalid FE"); 188182fbe45STzanio const mfem::Array<int> &dof_map_1d = h1_fe1d->GetDofMap(); 189182fbe45STzanio for (int i = 0; i < ir.GetNPoints(); i++) { 190182fbe45STzanio const mfem::IntegrationPoint &ip = ir.IntPoint(i); 191182fbe45STzanio qref1d(i) = ip.x; 192182fbe45STzanio qweight1d(i) = ip.weight; 193182fbe45STzanio fe1d->CalcShape(ip, shape_i); 194182fbe45STzanio fe1d->CalcDShape(ip, grad_i); 195182fbe45STzanio for (int j = 0; j < shape1d.Height(); j++) { 196182fbe45STzanio shape1d(j,i) = shape_i(dof_map_1d[j]); 197182fbe45STzanio grad1d(j,i) = grad_i(dof_map_1d[j],0); 198182fbe45STzanio } 199182fbe45STzanio } 200182fbe45STzanio CeedBasisCreateTensorH1(ceed, mesh->Dimension(), fes->GetVDim(), order+1, 201182fbe45STzanio ir.GetNPoints(), shape1d.GetData(), 202182fbe45STzanio grad1d.GetData(), qref1d.GetData(), 203182fbe45STzanio qweight1d.GetData(), basis); 204182fbe45STzanio 205182fbe45STzanio const mfem::Table &el_dof = fes->GetElementToDofTable(); 206182fbe45STzanio mfem::Array<int> tp_el_dof(el_dof.Size_of_connections()); 207182fbe45STzanio for (int i = 0; i < mesh->GetNE(); i++) { 208182fbe45STzanio const int el_offset = fe->GetDof()*i; 209182fbe45STzanio for (int j = 0; j < fe->GetDof(); j++) { 210182fbe45STzanio tp_el_dof[j + el_offset] = el_dof.GetJ()[dof_map[j] + el_offset]; 211182fbe45STzanio } 212182fbe45STzanio } 213182fbe45STzanio CeedElemRestrictionCreate(ceed, mesh->GetNE(), fe->GetDof(), 2147ca8db16Sjeremylt fes->GetNDofs(), fes->GetVDim(), CEED_MEM_HOST, CEED_COPY_VALUES, 215182fbe45STzanio tp_el_dof.GetData(), restr); 216182fbe45STzanio } 217182fbe45STzanio 218182fbe45STzanio public: 219182fbe45STzanio /// Constructor. Assumes @a fes is a scalar FE space. 220182fbe45STzanio CeedDiffusionOperator(Ceed ceed, const mfem::FiniteElementSpace *fes) 221182fbe45STzanio : Operator(fes->GetNDofs()), 222182fbe45STzanio fes(fes) { 223182fbe45STzanio mfem::Mesh *mesh = fes->GetMesh(); 224182fbe45STzanio const int order = fes->GetOrder(0); 225182fbe45STzanio const int ir_order = 2*(order + 2) - 1; // <----- 226182fbe45STzanio const mfem::IntegrationRule &ir = 227182fbe45STzanio mfem::IntRules.Get(mfem::Geometry::SEGMENT, ir_order); 228a48d94bfSjeremylt CeedInt nqpts, nelem = mesh->GetNE(), dim = mesh->SpaceDimension(); 229182fbe45STzanio 230182fbe45STzanio FESpace2Ceed(fes, ir, ceed, &basis, &restr); 231182fbe45STzanio 232182fbe45STzanio const mfem::FiniteElementSpace *mesh_fes = mesh->GetNodalFESpace(); 233182fbe45STzanio MFEM_VERIFY(mesh_fes, "the Mesh has no nodal FE space"); 234182fbe45STzanio FESpace2Ceed(mesh_fes, ir, ceed, &mesh_basis, &mesh_restr); 235a48d94bfSjeremylt CeedBasisGetNumQuadraturePoints(basis, &nqpts); 236182fbe45STzanio 237135a076eSjeremylt CeedElemRestrictionCreateIdentity(ceed, nelem, nqpts*dim*(dim+1)/2, 238135a076eSjeremylt nqpts*nelem*dim*(dim+1)/2, 1, &restr_i); 239135a076eSjeremylt CeedElemRestrictionCreateIdentity(ceed, nelem, nqpts, 240135a076eSjeremylt nqpts*nelem, 1, &mesh_restr_i); 241135a076eSjeremylt 242182fbe45STzanio CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords); 243182fbe45STzanio CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, 244182fbe45STzanio mesh->GetNodes()->GetData()); 245182fbe45STzanio 246a48d94bfSjeremylt CeedVectorCreate(ceed, nelem*nqpts*dim*(dim+1)/2, &rho); 2477ca8db16Sjeremylt 2487ca8db16Sjeremylt // Context data to be passed to the 'f_build_diff' Q-function. 2497ca8db16Sjeremylt build_ctx.dim = mesh->Dimension(); 2507ca8db16Sjeremylt build_ctx.space_dim = mesh->SpaceDimension(); 2517ca8db16Sjeremylt 2527ca8db16Sjeremylt // Create the Q-function that builds the diff operator (i.e. computes its 2537ca8db16Sjeremylt // quadrature data) and set its context data. 25454251743Sjeremylt CeedQFunctionCreateInterior(ceed, 1, f_build_diff, 255182fbe45STzanio __FILE__":f_build_diff", &build_qfunc); 256a48d94bfSjeremylt CeedQFunctionAddInput(build_qfunc, "dx", dim, CEED_EVAL_GRAD); 2577ca8db16Sjeremylt CeedQFunctionAddInput(build_qfunc, "weights", 1, CEED_EVAL_WEIGHT); 258a48d94bfSjeremylt CeedQFunctionAddOutput(build_qfunc, "rho", dim*(dim+1)/2, CEED_EVAL_NONE); 2597ca8db16Sjeremylt CeedQFunctionSetContext(build_qfunc, &build_ctx, sizeof(build_ctx)); 26054251743Sjeremylt 2617ca8db16Sjeremylt // Create the operator that builds the quadrature data for the diff operator. 26254251743Sjeremylt CeedOperatorCreate(ceed, build_qfunc, NULL, NULL, &build_oper); 2637ca8db16Sjeremylt CeedOperatorSetField(build_oper, "dx", mesh_restr, mesh_basis, 2647ca8db16Sjeremylt CEED_VECTOR_ACTIVE); 265135a076eSjeremylt CeedOperatorSetField(build_oper, "weights", mesh_restr_i, 2667ca8db16Sjeremylt mesh_basis, CEED_VECTOR_NONE); 267135a076eSjeremylt CeedOperatorSetField(build_oper, "rho", restr_i, 2687ca8db16Sjeremylt CEED_BASIS_COLOCATED, CEED_VECTOR_ACTIVE); 26954251743Sjeremylt 270a48d94bfSjeremylt // Compute the quadrature data for the diff operator. 2717ca8db16Sjeremylt CeedOperatorApply(build_oper, node_coords, rho, 272182fbe45STzanio CEED_REQUEST_IMMEDIATE); 273182fbe45STzanio 2747ca8db16Sjeremylt // Create the Q-function that defines the action of the diff operator. 27554251743Sjeremylt CeedQFunctionCreateInterior(ceed, 1, f_apply_diff, 276a48d94bfSjeremylt __FILE__":f_apply_diff", &apply_qfunc); 277a48d94bfSjeremylt CeedQFunctionAddInput(apply_qfunc, "u", 1, CEED_EVAL_GRAD); 278a48d94bfSjeremylt CeedQFunctionAddInput(apply_qfunc, "rho", dim*(dim+1)/2, CEED_EVAL_NONE); 279a48d94bfSjeremylt CeedQFunctionAddOutput(apply_qfunc, "v", 1, CEED_EVAL_GRAD); 280a48d94bfSjeremylt CeedQFunctionSetContext(apply_qfunc, &build_ctx, sizeof(build_ctx)); 28154251743Sjeremylt 282a48d94bfSjeremylt // Create the diff operator. 28354251743Sjeremylt CeedOperatorCreate(ceed, apply_qfunc, NULL, NULL, &oper); 284a48d94bfSjeremylt CeedOperatorSetField(oper, "u", restr, basis, CEED_VECTOR_ACTIVE); 285135a076eSjeremylt CeedOperatorSetField(oper, "rho", restr_i, 2867ca8db16Sjeremylt CEED_BASIS_COLOCATED, rho); 287a48d94bfSjeremylt CeedOperatorSetField(oper, "v", restr, basis, CEED_VECTOR_ACTIVE); 288182fbe45STzanio 289182fbe45STzanio CeedVectorCreate(ceed, fes->GetNDofs(), &u); 290182fbe45STzanio CeedVectorCreate(ceed, fes->GetNDofs(), &v); 291182fbe45STzanio } 292182fbe45STzanio 293182fbe45STzanio /// Destructor 294182fbe45STzanio ~CeedDiffusionOperator() { 295182fbe45STzanio CeedVectorDestroy(&u); 2967ca8db16Sjeremylt CeedVectorDestroy(&v); 2977ca8db16Sjeremylt CeedVectorDestroy(&rho); 298182fbe45STzanio CeedVectorDestroy(&node_coords); 299182fbe45STzanio CeedElemRestrictionDestroy(&restr); 3007ca8db16Sjeremylt CeedElemRestrictionDestroy(&mesh_restr); 301135a076eSjeremylt CeedElemRestrictionDestroy(&restr_i); 302135a076eSjeremylt CeedElemRestrictionDestroy(&mesh_restr_i); 303182fbe45STzanio CeedBasisDestroy(&basis); 3047ca8db16Sjeremylt CeedBasisDestroy(&mesh_basis); 3057ca8db16Sjeremylt CeedQFunctionDestroy(&build_qfunc); 3067ca8db16Sjeremylt CeedOperatorDestroy(&build_oper); 3077ca8db16Sjeremylt CeedQFunctionDestroy(&apply_qfunc); 3087ca8db16Sjeremylt CeedOperatorDestroy(&oper); 309182fbe45STzanio } 310182fbe45STzanio 311182fbe45STzanio /// Operator action 312182fbe45STzanio virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const { 313182fbe45STzanio CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData()); 314182fbe45STzanio CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData()); 315182fbe45STzanio 31654251743Sjeremylt CeedOperatorApply(oper, u, v, CEED_REQUEST_IMMEDIATE); 317182fbe45STzanio } 318182fbe45STzanio }; 319