1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 /// @file 18 /// Diffusion operator example using MFEM 19 #include <ceed.h> 20 #include <mfem.hpp> 21 22 #include "bp3.h" 23 24 /// Wrapper for a diffusion CeedOperator as an mfem::Operator 25 class CeedDiffusionOperator : public mfem::Operator { 26 protected: 27 const mfem::FiniteElementSpace *fes; 28 CeedOperator build_oper, oper; 29 CeedBasis basis, mesh_basis; 30 CeedElemRestriction restr, mesh_restr, restr_i, mesh_restr_i; 31 CeedQFunction apply_qfunc, build_qfunc; 32 CeedVector node_coords, qdata; 33 34 BuildContext build_ctx; 35 36 CeedVector u, v; 37 38 static void FESpace2Ceed(const mfem::FiniteElementSpace *fes, 39 const mfem::IntegrationRule &ir, 40 Ceed ceed, CeedBasis *basis, 41 CeedElemRestriction *restr) { 42 mfem::Mesh *mesh = fes->GetMesh(); 43 const mfem::FiniteElement *fe = fes->GetFE(0); 44 const int order = fes->GetOrder(0); 45 mfem::Array<int> dof_map; 46 switch (mesh->Dimension()) { 47 case 1: { 48 const mfem::H1_SegmentElement *h1_fe = 49 dynamic_cast<const mfem::H1_SegmentElement *>(fe); 50 MFEM_VERIFY(h1_fe, "invalid FE"); 51 h1_fe->GetDofMap().Copy(dof_map); 52 break; 53 } 54 case 2: { 55 const mfem::H1_QuadrilateralElement *h1_fe = 56 dynamic_cast<const mfem::H1_QuadrilateralElement *>(fe); 57 MFEM_VERIFY(h1_fe, "invalid FE"); 58 h1_fe->GetDofMap().Copy(dof_map); 59 break; 60 } 61 case 3: { 62 const mfem::H1_HexahedronElement *h1_fe = 63 dynamic_cast<const mfem::H1_HexahedronElement *>(fe); 64 MFEM_VERIFY(h1_fe, "invalid FE"); 65 h1_fe->GetDofMap().Copy(dof_map); 66 break; 67 } 68 } 69 const mfem::FiniteElement *fe1d = 70 fes->FEColl()->FiniteElementForGeometry(mfem::Geometry::SEGMENT); 71 mfem::DenseMatrix shape1d(fe1d->GetDof(), ir.GetNPoints()); 72 mfem::DenseMatrix grad1d(fe1d->GetDof(), ir.GetNPoints()); 73 mfem::Vector qref1d(ir.GetNPoints()), qweight1d(ir.GetNPoints()); 74 mfem::Vector shape_i(shape1d.Height()); 75 mfem::DenseMatrix grad_i(grad1d.Height(), 1); 76 const mfem::H1_SegmentElement *h1_fe1d = 77 dynamic_cast<const mfem::H1_SegmentElement *>(fe1d); 78 MFEM_VERIFY(h1_fe1d, "invalid FE"); 79 const mfem::Array<int> &dof_map_1d = h1_fe1d->GetDofMap(); 80 for (int i = 0; i < ir.GetNPoints(); i++) { 81 const mfem::IntegrationPoint &ip = ir.IntPoint(i); 82 qref1d(i) = ip.x; 83 qweight1d(i) = ip.weight; 84 fe1d->CalcShape(ip, shape_i); 85 fe1d->CalcDShape(ip, grad_i); 86 for (int j = 0; j < shape1d.Height(); j++) { 87 shape1d(j,i) = shape_i(dof_map_1d[j]); 88 grad1d(j,i) = grad_i(dof_map_1d[j],0); 89 } 90 } 91 CeedBasisCreateTensorH1(ceed, mesh->Dimension(), fes->GetVDim(), order+1, 92 ir.GetNPoints(), shape1d.GetData(), 93 grad1d.GetData(), qref1d.GetData(), 94 qweight1d.GetData(), basis); 95 96 const mfem::Table &el_dof = fes->GetElementToDofTable(); 97 mfem::Array<int> tp_el_dof(el_dof.Size_of_connections()); 98 for (int i = 0; i < mesh->GetNE(); i++) { 99 const int el_offset = fe->GetDof()*i; 100 for (int j = 0; j < fe->GetDof(); j++) { 101 tp_el_dof[j + el_offset] = el_dof.GetJ()[dof_map[j] + el_offset]; 102 } 103 } 104 CeedElemRestrictionCreate(ceed, mesh->GetNE(), fe->GetDof(), 105 fes->GetNDofs(), fes->GetVDim(), CEED_MEM_HOST, 106 CEED_COPY_VALUES, tp_el_dof.GetData(), restr); 107 } 108 109 public: 110 /// Constructor. Assumes @a fes is a scalar FE space. 111 CeedDiffusionOperator(Ceed ceed, const mfem::FiniteElementSpace *fes) 112 : Operator(fes->GetNDofs()), 113 fes(fes) { 114 mfem::Mesh *mesh = fes->GetMesh(); 115 const int order = fes->GetOrder(0); 116 const int ir_order = 2*(order + 2) - 1; // <----- 117 const mfem::IntegrationRule &ir = 118 mfem::IntRules.Get(mfem::Geometry::SEGMENT, ir_order); 119 CeedInt nqpts, nelem = mesh->GetNE(), dim = mesh->SpaceDimension(), 120 ncompx = dim; 121 122 FESpace2Ceed(fes, ir, ceed, &basis, &restr); 123 124 const mfem::FiniteElementSpace *mesh_fes = mesh->GetNodalFESpace(); 125 MFEM_VERIFY(mesh_fes, "the Mesh has no nodal FE space"); 126 FESpace2Ceed(mesh_fes, ir, ceed, &mesh_basis, &mesh_restr); 127 CeedBasisGetNumQuadraturePoints(basis, &nqpts); 128 129 CeedElemRestrictionCreateIdentity(ceed, nelem, nqpts, 130 nqpts*nelem, dim*(dim+1)/2, &restr_i); 131 CeedElemRestrictionCreateIdentity(ceed, nelem, nqpts, 132 nqpts*nelem, 1, &mesh_restr_i); 133 134 CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords); 135 CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, 136 mesh->GetNodes()->GetData()); 137 138 CeedVectorCreate(ceed, nelem*nqpts*dim*(dim+1)/2, &qdata); 139 140 // Context data to be passed to the 'f_build_diff' Q-function. 141 build_ctx.dim = mesh->Dimension(); 142 build_ctx.space_dim = mesh->SpaceDimension(); 143 144 // Create the Q-function that builds the diff operator (i.e. computes its 145 // quadrature data) and set its context data. 146 CeedQFunctionCreateInterior(ceed, 1, f_build_diff, 147 f_build_diff_loc, &build_qfunc); 148 CeedQFunctionAddInput(build_qfunc, "dx", ncompx*dim, CEED_EVAL_GRAD); 149 CeedQFunctionAddInput(build_qfunc, "weights", 1, CEED_EVAL_WEIGHT); 150 CeedQFunctionAddOutput(build_qfunc, "qdata", dim*(dim+1)/2, CEED_EVAL_NONE); 151 CeedQFunctionSetContext(build_qfunc, &build_ctx, sizeof(build_ctx)); 152 153 // Create the operator that builds the quadrature data for the diff operator. 154 CeedOperatorCreate(ceed, build_qfunc, NULL, NULL, &build_oper); 155 CeedOperatorSetField(build_oper, "dx", mesh_restr, CEED_NOTRANSPOSE, 156 mesh_basis, CEED_VECTOR_ACTIVE); 157 CeedOperatorSetField(build_oper, "weights", mesh_restr_i, CEED_NOTRANSPOSE, 158 mesh_basis, CEED_VECTOR_NONE); 159 CeedOperatorSetField(build_oper, "qdata", restr_i, CEED_NOTRANSPOSE, 160 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 161 162 // Compute the quadrature data for the diff operator. 163 CeedOperatorApply(build_oper, node_coords, qdata, 164 CEED_REQUEST_IMMEDIATE); 165 166 // Create the Q-function that defines the action of the diff operator. 167 CeedQFunctionCreateInterior(ceed, 1, f_apply_diff, 168 f_apply_diff_loc, &apply_qfunc); 169 CeedQFunctionAddInput(apply_qfunc, "u", dim, CEED_EVAL_GRAD); 170 CeedQFunctionAddInput(apply_qfunc, "qdata", dim*(dim+1)/2, CEED_EVAL_NONE); 171 CeedQFunctionAddOutput(apply_qfunc, "v", dim, CEED_EVAL_GRAD); 172 CeedQFunctionSetContext(apply_qfunc, &build_ctx, sizeof(build_ctx)); 173 174 // Create the diff operator. 175 CeedOperatorCreate(ceed, apply_qfunc, NULL, NULL, &oper); 176 CeedOperatorSetField(oper, "u", restr, CEED_NOTRANSPOSE, 177 basis, CEED_VECTOR_ACTIVE); 178 CeedOperatorSetField(oper, "qdata", restr_i, CEED_NOTRANSPOSE, 179 CEED_BASIS_COLLOCATED, qdata); 180 CeedOperatorSetField(oper, "v", restr, CEED_NOTRANSPOSE, 181 basis, CEED_VECTOR_ACTIVE); 182 183 CeedVectorCreate(ceed, fes->GetNDofs(), &u); 184 CeedVectorCreate(ceed, fes->GetNDofs(), &v); 185 } 186 187 /// Destructor 188 ~CeedDiffusionOperator() { 189 CeedVectorDestroy(&u); 190 CeedVectorDestroy(&v); 191 CeedVectorDestroy(&qdata); 192 CeedVectorDestroy(&node_coords); 193 CeedElemRestrictionDestroy(&restr); 194 CeedElemRestrictionDestroy(&mesh_restr); 195 CeedElemRestrictionDestroy(&restr_i); 196 CeedElemRestrictionDestroy(&mesh_restr_i); 197 CeedBasisDestroy(&basis); 198 CeedBasisDestroy(&mesh_basis); 199 CeedQFunctionDestroy(&build_qfunc); 200 CeedOperatorDestroy(&build_oper); 201 CeedQFunctionDestroy(&apply_qfunc); 202 CeedOperatorDestroy(&oper); 203 } 204 205 /// Operator action 206 virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const { 207 CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData()); 208 CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData()); 209 210 CeedOperatorApply(oper, u, v, CEED_REQUEST_IMMEDIATE); 211 CeedVectorSyncArray(v, CEED_MEM_HOST); 212 } 213 }; 214