1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Mass operator example using MFEM 10 11 #include <ceed.h> 12 13 #include <mfem.hpp> 14 15 #include "bp1.h" 16 17 /// Wrapper for a mass CeedOperator as an mfem::Operator 18 class CeedMassOperator : public mfem::Operator { 19 protected: 20 const mfem::FiniteElementSpace *fes; 21 CeedOperator build_oper, oper; 22 CeedBasis basis, mesh_basis; 23 CeedElemRestriction restr, mesh_restr, restr_i, mesh_restr_i; 24 CeedQFunction apply_qfunc, build_qfunc; 25 CeedQFunctionContext build_ctx; 26 CeedVector node_coords, qdata; 27 CeedVector u, v; 28 29 BuildContext build_ctx_data; 30 31 static void FESpace2Ceed(const mfem::FiniteElementSpace *fes, const mfem::IntegrationRule &ir, Ceed ceed, CeedBasis *basis, 32 CeedElemRestriction *restr) { 33 mfem::Mesh *mesh = fes->GetMesh(); 34 const mfem::FiniteElement *fe = fes->GetFE(0); 35 const int order = fes->GetOrder(0); 36 mfem::Array<int> dof_map; 37 switch (mesh->Dimension()) { 38 case 1: { 39 const mfem::H1_SegmentElement *h1_fe = dynamic_cast<const mfem::H1_SegmentElement *>(fe); 40 MFEM_VERIFY(h1_fe, "invalid FE"); 41 h1_fe->GetDofMap().Copy(dof_map); 42 break; 43 } 44 case 2: { 45 const mfem::H1_QuadrilateralElement *h1_fe = dynamic_cast<const mfem::H1_QuadrilateralElement *>(fe); 46 MFEM_VERIFY(h1_fe, "invalid FE"); 47 h1_fe->GetDofMap().Copy(dof_map); 48 break; 49 } 50 case 3: { 51 const mfem::H1_HexahedronElement *h1_fe = dynamic_cast<const mfem::H1_HexahedronElement *>(fe); 52 MFEM_VERIFY(h1_fe, "invalid FE"); 53 h1_fe->GetDofMap().Copy(dof_map); 54 break; 55 } 56 } 57 const mfem::FiniteElement *fe1d = fes->FEColl()->FiniteElementForGeometry(mfem::Geometry::SEGMENT); 58 mfem::DenseMatrix shape1d(fe1d->GetDof(), ir.GetNPoints()); 59 mfem::DenseMatrix grad_1d(fe1d->GetDof(), ir.GetNPoints()); 60 mfem::Vector q_ref_1d(ir.GetNPoints()), q_weight_1d(ir.GetNPoints()); 61 mfem::Vector shape_i(shape1d.Height()); 62 mfem::DenseMatrix grad_i(grad_1d.Height(), 1); 63 const mfem::H1_SegmentElement *h1_fe1d = dynamic_cast<const mfem::H1_SegmentElement *>(fe1d); 64 MFEM_VERIFY(h1_fe1d, "invalid FE"); 65 const mfem::Array<int> &dof_map_1d = h1_fe1d->GetDofMap(); 66 for (int i = 0; i < ir.GetNPoints(); i++) { 67 const mfem::IntegrationPoint &ip = ir.IntPoint(i); 68 q_ref_1d(i) = ip.x; 69 q_weight_1d(i) = ip.weight; 70 fe1d->CalcShape(ip, shape_i); 71 fe1d->CalcDShape(ip, grad_i); 72 for (int j = 0; j < shape1d.Height(); j++) { 73 shape1d(j, i) = shape_i(dof_map_1d[j]); 74 grad_1d(j, i) = grad_i(dof_map_1d[j], 0); 75 } 76 } 77 CeedBasisCreateTensorH1(ceed, mesh->Dimension(), fes->GetVDim(), order + 1, ir.GetNPoints(), shape1d.GetData(), grad_1d.GetData(), 78 q_ref_1d.GetData(), q_weight_1d.GetData(), basis); 79 80 const mfem::Table &el_dof = fes->GetElementToDofTable(); 81 mfem::Array<int> tp_el_dof(el_dof.Size_of_connections()); 82 for (int i = 0; i < mesh->GetNE(); i++) { 83 const int el_offset = fe->GetDof() * i; 84 for (int j = 0; j < fe->GetDof(); j++) { 85 tp_el_dof[j + el_offset] = el_dof.GetJ()[dof_map[j] + el_offset]; 86 } 87 } 88 CeedElemRestrictionCreate(ceed, mesh->GetNE(), fe->GetDof(), fes->GetVDim(), fes->GetNDofs(), (fes->GetVDim()) * (fes->GetNDofs()), CEED_MEM_HOST, 89 CEED_COPY_VALUES, tp_el_dof.GetData(), restr); 90 } 91 92 public: 93 /// Constructor. Assumes @a fes is a scalar FE space. 94 CeedMassOperator(Ceed ceed, const mfem::FiniteElementSpace *fes) : Operator(fes->GetNDofs()), fes(fes) { 95 mfem::Mesh *mesh = fes->GetMesh(); 96 const int order = fes->GetOrder(0); 97 const int ir_order = 2 * (order + 2) - 1; // <----- 98 const mfem::IntegrationRule &ir = mfem::IntRules.Get(mfem::Geometry::SEGMENT, ir_order); 99 CeedInt num_elem = mesh->GetNE(), dim = mesh->SpaceDimension(), ncompx = dim, nqpts; 100 101 FESpace2Ceed(fes, ir, ceed, &basis, &restr); 102 103 const mfem::FiniteElementSpace *mesh_fes = mesh->GetNodalFESpace(); 104 MFEM_VERIFY(mesh_fes, "the Mesh has no nodal FE space"); 105 FESpace2Ceed(mesh_fes, ir, ceed, &mesh_basis, &mesh_restr); 106 CeedBasisGetNumQuadraturePoints(basis, &nqpts); 107 108 CeedInt strides[3] = {1, nqpts, nqpts}; 109 CeedElemRestrictionCreateStrided(ceed, num_elem, nqpts, 1, nqpts * num_elem, strides, &restr_i); 110 111 CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords); 112 CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, mesh->GetNodes()->GetData()); 113 114 CeedVectorCreate(ceed, num_elem * nqpts, &qdata); 115 116 // Context data to be passed to the 'f_build_mass' Q-function. 117 build_ctx_data.dim = mesh->Dimension(); 118 build_ctx_data.space_dim = dim; 119 CeedQFunctionContextCreate(ceed, &build_ctx); 120 CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data); 121 122 // Create the Q-function that builds the mass operator (i.e. computes its quadrature data) and set its context data. 123 CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &build_qfunc); 124 CeedQFunctionAddInput(build_qfunc, "dx", ncompx * dim, CEED_EVAL_GRAD); 125 CeedQFunctionAddInput(build_qfunc, "weights", 1, CEED_EVAL_WEIGHT); 126 CeedQFunctionAddOutput(build_qfunc, "qdata", 1, CEED_EVAL_NONE); 127 CeedQFunctionSetContext(build_qfunc, build_ctx); 128 129 // Create the operator that builds the quadrature data for the mass operator. 130 CeedOperatorCreate(ceed, build_qfunc, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &build_oper); 131 CeedOperatorSetField(build_oper, "dx", mesh_restr, mesh_basis, CEED_VECTOR_ACTIVE); 132 CeedOperatorSetField(build_oper, "weights", CEED_ELEMRESTRICTION_NONE, mesh_basis, CEED_VECTOR_NONE); 133 CeedOperatorSetField(build_oper, "qdata", restr_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 134 135 // Compute the quadrature data for the mass operator. 136 CeedOperatorApply(build_oper, node_coords, qdata, CEED_REQUEST_IMMEDIATE); 137 138 // Create the Q-function that defines the action of the mass operator. 139 CeedQFunctionCreateInterior(ceed, 1, f_apply_mass, f_apply_mass_loc, &apply_qfunc); 140 CeedQFunctionAddInput(apply_qfunc, "u", 1, CEED_EVAL_INTERP); 141 CeedQFunctionAddInput(apply_qfunc, "qdata", 1, CEED_EVAL_NONE); 142 CeedQFunctionAddOutput(apply_qfunc, "v", 1, CEED_EVAL_INTERP); 143 144 // Create the mass operator. 145 CeedOperatorCreate(ceed, apply_qfunc, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &oper); 146 CeedOperatorSetField(oper, "u", restr, basis, CEED_VECTOR_ACTIVE); 147 CeedOperatorSetField(oper, "qdata", restr_i, CEED_BASIS_NONE, qdata); 148 CeedOperatorSetField(oper, "v", restr, basis, CEED_VECTOR_ACTIVE); 149 150 CeedVectorCreate(ceed, fes->GetNDofs(), &u); 151 CeedVectorCreate(ceed, fes->GetNDofs(), &v); 152 } 153 154 /// Destructor 155 ~CeedMassOperator() { 156 CeedVectorDestroy(&u); 157 CeedVectorDestroy(&v); 158 CeedVectorDestroy(&qdata); 159 CeedVectorDestroy(&node_coords); 160 CeedOperatorDestroy(&build_oper); 161 CeedQFunctionDestroy(&build_qfunc); 162 CeedOperatorDestroy(&oper); 163 CeedQFunctionDestroy(&apply_qfunc); 164 CeedQFunctionContextDestroy(&build_ctx); 165 CeedBasisDestroy(&basis); 166 CeedBasisDestroy(&mesh_basis); 167 CeedElemRestrictionDestroy(&restr); 168 CeedElemRestrictionDestroy(&mesh_restr); 169 CeedElemRestrictionDestroy(&restr_i); 170 } 171 172 /// Operator action 173 virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const { 174 CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData()); 175 CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData()); 176 CeedOperatorApply(oper, u, v, CEED_REQUEST_IMMEDIATE); 177 CeedVectorSyncArray(v, CEED_MEM_HOST); 178 } 179 }; 180