xref: /libCEED/examples/mfem/bp1.hpp (revision 56d8cfc241e0df885b47de7ad2d12edb4b8ae247)
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 /// Mass operator example using MFEM
19 #include <ceed.h>
20 #include <mfem.hpp>
21 
22 #include "bp1.h"
23 
24 /// Wrapper for a mass CeedOperator as an mfem::Operator
25 class CeedMassOperator : 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   CeedVector u, v;
34 
35   BuildContext build_ctx;
36 
37   static void FESpace2Ceed(const mfem::FiniteElementSpace *fes,
38                            const mfem::IntegrationRule &ir,
39                            Ceed ceed, CeedBasis *basis,
40                            CeedElemRestriction *restr) {
41     mfem::Mesh *mesh = fes->GetMesh();
42     const mfem::FiniteElement *fe = fes->GetFE(0);
43     const int order = fes->GetOrder(0);
44     mfem::Array<int> dof_map;
45     switch (mesh->Dimension()) {
46     case 1: {
47       const mfem::H1_SegmentElement *h1_fe =
48         dynamic_cast<const mfem::H1_SegmentElement *>(fe);
49       MFEM_VERIFY(h1_fe, "invalid FE");
50       h1_fe->GetDofMap().Copy(dof_map);
51       break;
52     }
53     case 2: {
54       const mfem::H1_QuadrilateralElement *h1_fe =
55         dynamic_cast<const mfem::H1_QuadrilateralElement *>(fe);
56       MFEM_VERIFY(h1_fe, "invalid FE");
57       h1_fe->GetDofMap().Copy(dof_map);
58       break;
59     }
60     case 3: {
61       const mfem::H1_HexahedronElement *h1_fe =
62         dynamic_cast<const mfem::H1_HexahedronElement *>(fe);
63       MFEM_VERIFY(h1_fe, "invalid FE");
64       h1_fe->GetDofMap().Copy(dof_map);
65       break;
66     }
67     }
68     const mfem::FiniteElement *fe1d =
69       fes->FEColl()->FiniteElementForGeometry(mfem::Geometry::SEGMENT);
70     mfem::DenseMatrix shape1d(fe1d->GetDof(), ir.GetNPoints());
71     mfem::DenseMatrix grad1d(fe1d->GetDof(), ir.GetNPoints());
72     mfem::Vector qref1d(ir.GetNPoints()), qweight1d(ir.GetNPoints());
73     mfem::Vector shape_i(shape1d.Height());
74     mfem::DenseMatrix grad_i(grad1d.Height(), 1);
75     const mfem::H1_SegmentElement *h1_fe1d =
76       dynamic_cast<const mfem::H1_SegmentElement *>(fe1d);
77     MFEM_VERIFY(h1_fe1d, "invalid FE");
78     const mfem::Array<int> &dof_map_1d = h1_fe1d->GetDofMap();
79     for (int i = 0; i < ir.GetNPoints(); i++) {
80       const mfem::IntegrationPoint &ip = ir.IntPoint(i);
81       qref1d(i) = ip.x;
82       qweight1d(i) = ip.weight;
83       fe1d->CalcShape(ip, shape_i);
84       fe1d->CalcDShape(ip, grad_i);
85       for (int j = 0; j < shape1d.Height(); j++) {
86         shape1d(j,i) = shape_i(dof_map_1d[j]);
87         grad1d(j,i) = grad_i(dof_map_1d[j],0);
88       }
89     }
90     CeedBasisCreateTensorH1(ceed, mesh->Dimension(), fes->GetVDim(), order+1,
91                             ir.GetNPoints(), shape1d.GetData(),
92                             grad1d.GetData(), qref1d.GetData(),
93                             qweight1d.GetData(), basis);
94 
95     const mfem::Table &el_dof = fes->GetElementToDofTable();
96     mfem::Array<int> tp_el_dof(el_dof.Size_of_connections());
97     for (int i = 0; i < mesh->GetNE(); i++) {
98       const int el_offset = fe->GetDof()*i;
99       for (int j = 0; j < fe->GetDof(); j++) {
100         tp_el_dof[j + el_offset] = el_dof.GetJ()[dof_map[j] + el_offset];
101       }
102     }
103     CeedInterlaceMode imode = CEED_NONINTERLACED;
104     CeedElemRestrictionCreate(ceed, imode, 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   CeedMassOperator(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 nelem = mesh->GetNE(), dim = mesh->SpaceDimension(),
120             ncompx = dim, nqpts;
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     CeedInt strides[3] = {1, nqpts, nqpts};
130     CeedElemRestrictionCreateStrided(ceed, nelem, nqpts, nqpts*nelem, 1,
131                                      strides, &restr_i);
132 
133     CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords);
134     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER,
135                        mesh->GetNodes()->GetData());
136 
137     CeedVectorCreate(ceed, nelem*nqpts, &qdata);
138 
139     // Context data to be passed to the 'f_build_mass' Q-function.
140     build_ctx.dim = mesh->Dimension();
141     build_ctx.space_dim = dim;
142 
143     // Create the Q-function that builds the mass operator (i.e. computes its
144     // quadrature data) and set its context data.
145     CeedQFunctionCreateInterior(ceed, 1, f_build_mass,
146                                 f_build_mass_loc, &build_qfunc);
147     CeedQFunctionAddInput(build_qfunc, "dx", ncompx*dim,
148                           CEED_EVAL_GRAD);
149     CeedQFunctionAddInput(build_qfunc, "weights", 1, CEED_EVAL_WEIGHT);
150     CeedQFunctionAddOutput(build_qfunc, "qdata", 1, CEED_EVAL_NONE);
151     CeedQFunctionSetContext(build_qfunc, &build_ctx, sizeof(build_ctx));
152 
153     // Create the operator that builds the quadrature data for the mass operator.
154     CeedOperatorCreate(ceed, build_qfunc, CEED_QFUNCTION_NONE,
155                        CEED_QFUNCTION_NONE, &build_oper);
156     CeedOperatorSetField(build_oper, "dx", mesh_restr, mesh_basis,
157                          CEED_VECTOR_ACTIVE);
158     CeedOperatorSetField(build_oper, "weights", CEED_ELEMRESTRICTION_NONE,
159                          mesh_basis, CEED_VECTOR_NONE);
160     CeedOperatorSetField(build_oper, "qdata", restr_i, CEED_BASIS_COLLOCATED,
161                          CEED_VECTOR_ACTIVE);
162 
163     // Compute the quadrature data for the mass operator.
164     CeedOperatorApply(build_oper, node_coords, qdata,
165                       CEED_REQUEST_IMMEDIATE);
166 
167     // Create the Q-function that defines the action of the mass operator.
168     CeedQFunctionCreateInterior(ceed, 1, f_apply_mass,
169                                 f_apply_mass_loc, &apply_qfunc);
170     CeedQFunctionAddInput(apply_qfunc, "u", 1, CEED_EVAL_INTERP);
171     CeedQFunctionAddInput(apply_qfunc, "qdata", 1, CEED_EVAL_NONE);
172     CeedQFunctionAddOutput(apply_qfunc, "v", 1, CEED_EVAL_INTERP);
173 
174     // Create the mass operator.
175     CeedOperatorCreate(ceed, apply_qfunc, CEED_QFUNCTION_NONE,
176                        CEED_QFUNCTION_NONE, &oper);
177     CeedOperatorSetField(oper, "u", restr, basis, CEED_VECTOR_ACTIVE);
178     CeedOperatorSetField(oper, "qdata", restr_i, CEED_BASIS_COLLOCATED, qdata);
179     CeedOperatorSetField(oper, "v", restr, basis, CEED_VECTOR_ACTIVE);
180 
181     CeedVectorCreate(ceed, fes->GetNDofs(), &u);
182     CeedVectorCreate(ceed, fes->GetNDofs(), &v);
183   }
184 
185   /// Destructor
186   ~CeedMassOperator() {
187     CeedVectorDestroy(&u);
188     CeedVectorDestroy(&v);
189     CeedVectorDestroy(&qdata);
190     CeedVectorDestroy(&node_coords);
191     CeedOperatorDestroy(&build_oper);
192     CeedQFunctionDestroy(&build_qfunc);
193     CeedOperatorDestroy(&oper);
194     CeedQFunctionDestroy(&apply_qfunc);
195     CeedBasisDestroy(&basis);
196     CeedBasisDestroy(&mesh_basis);
197     CeedElemRestrictionDestroy(&restr);
198     CeedElemRestrictionDestroy(&mesh_restr);
199     CeedElemRestrictionDestroy(&restr_i);
200   }
201 
202   /// Operator action
203   virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const {
204     CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData());
205     CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData());
206     CeedOperatorApply(oper, u, v, CEED_REQUEST_IMMEDIATE);
207     CeedVectorSyncArray(v, CEED_MEM_HOST);
208   }
209 };
210