xref: /libCEED/examples/mfem/bp1.hpp (revision dc8efd83546faf0200bf0bfcfb1678fae1874cc5)
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 
20 #include <ceed.h>
21 #include <mfem.hpp>
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   CeedQFunctionContext build_ctx;
33   CeedVector node_coords, qdata;
34   CeedVector u, v;
35 
36   BuildContext build_ctx_data;
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->GetVDim(), fes->GetNDofs(),
106                               (fes->GetVDim())*(fes->GetNDofs()),
107                               CEED_MEM_HOST, CEED_COPY_VALUES,
108                               tp_el_dof.GetData(), restr);
109   }
110 
111  public:
112   /// Constructor. Assumes @a fes is a scalar FE space.
113   CeedMassOperator(Ceed ceed, const mfem::FiniteElementSpace *fes)
114     : Operator(fes->GetNDofs()),
115       fes(fes) {
116     mfem::Mesh *mesh = fes->GetMesh();
117     const int order = fes->GetOrder(0);
118     const int ir_order = 2*(order + 2) - 1; // <-----
119     const mfem::IntegrationRule &ir =
120       mfem::IntRules.Get(mfem::Geometry::SEGMENT, ir_order);
121     CeedInt nelem = mesh->GetNE(), dim = mesh->SpaceDimension(),
122             ncompx = dim, nqpts;
123 
124     FESpace2Ceed(fes, ir, ceed, &basis, &restr);
125 
126     const mfem::FiniteElementSpace *mesh_fes = mesh->GetNodalFESpace();
127     MFEM_VERIFY(mesh_fes, "the Mesh has no nodal FE space");
128     FESpace2Ceed(mesh_fes, ir, ceed, &mesh_basis, &mesh_restr);
129     CeedBasisGetNumQuadraturePoints(basis, &nqpts);
130 
131     CeedInt strides[3] = {1, nqpts, nqpts};
132     CeedElemRestrictionCreateStrided(ceed, nelem, nqpts, 1, nqpts*nelem,
133                                      strides, &restr_i);
134 
135     CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords);
136     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER,
137                        mesh->GetNodes()->GetData());
138 
139     CeedVectorCreate(ceed, nelem*nqpts, &qdata);
140 
141     // Context data to be passed to the 'f_build_mass' Q-function.
142     build_ctx_data.dim = mesh->Dimension();
143     build_ctx_data.space_dim = dim;
144     CeedQFunctionContextCreate(ceed, &build_ctx);
145     CeedQFunctionContextSetData(build_ctx, CEED_MEM_HOST, CEED_USE_POINTER,
146                                 sizeof(build_ctx_data), &build_ctx_data);
147 
148     // Create the Q-function that builds the mass operator (i.e. computes its
149     // quadrature data) and set its context data.
150     CeedQFunctionCreateInterior(ceed, 1, f_build_mass,
151                                 f_build_mass_loc, &build_qfunc);
152     CeedQFunctionAddInput(build_qfunc, "dx", ncompx*dim,
153                           CEED_EVAL_GRAD);
154     CeedQFunctionAddInput(build_qfunc, "weights", 1, CEED_EVAL_WEIGHT);
155     CeedQFunctionAddOutput(build_qfunc, "qdata", 1, CEED_EVAL_NONE);
156     CeedQFunctionSetContext(build_qfunc, build_ctx);
157 
158     // Create the operator that builds the quadrature data for the mass operator.
159     CeedOperatorCreate(ceed, build_qfunc, CEED_QFUNCTION_NONE,
160                        CEED_QFUNCTION_NONE, &build_oper);
161     CeedOperatorSetField(build_oper, "dx", mesh_restr, mesh_basis,
162                          CEED_VECTOR_ACTIVE);
163     CeedOperatorSetField(build_oper, "weights", CEED_ELEMRESTRICTION_NONE,
164                          mesh_basis, CEED_VECTOR_NONE);
165     CeedOperatorSetField(build_oper, "qdata", restr_i, CEED_BASIS_COLLOCATED,
166                          CEED_VECTOR_ACTIVE);
167 
168     // Compute the quadrature data for the mass operator.
169     CeedOperatorApply(build_oper, node_coords, qdata,
170                       CEED_REQUEST_IMMEDIATE);
171 
172     // Create the Q-function that defines the action of the mass operator.
173     CeedQFunctionCreateInterior(ceed, 1, f_apply_mass,
174                                 f_apply_mass_loc, &apply_qfunc);
175     CeedQFunctionAddInput(apply_qfunc, "u", 1, CEED_EVAL_INTERP);
176     CeedQFunctionAddInput(apply_qfunc, "qdata", 1, CEED_EVAL_NONE);
177     CeedQFunctionAddOutput(apply_qfunc, "v", 1, CEED_EVAL_INTERP);
178 
179     // Create the mass operator.
180     CeedOperatorCreate(ceed, apply_qfunc, CEED_QFUNCTION_NONE,
181                        CEED_QFUNCTION_NONE, &oper);
182     CeedOperatorSetField(oper, "u", restr, basis, CEED_VECTOR_ACTIVE);
183     CeedOperatorSetField(oper, "qdata", restr_i, CEED_BASIS_COLLOCATED, qdata);
184     CeedOperatorSetField(oper, "v", restr, basis, CEED_VECTOR_ACTIVE);
185 
186     CeedVectorCreate(ceed, fes->GetNDofs(), &u);
187     CeedVectorCreate(ceed, fes->GetNDofs(), &v);
188   }
189 
190   /// Destructor
191   ~CeedMassOperator() {
192     CeedVectorDestroy(&u);
193     CeedVectorDestroy(&v);
194     CeedVectorDestroy(&qdata);
195     CeedVectorDestroy(&node_coords);
196     CeedOperatorDestroy(&build_oper);
197     CeedQFunctionDestroy(&build_qfunc);
198     CeedOperatorDestroy(&oper);
199     CeedQFunctionDestroy(&apply_qfunc);
200     CeedQFunctionContextDestroy(&build_ctx);
201     CeedBasisDestroy(&basis);
202     CeedBasisDestroy(&mesh_basis);
203     CeedElemRestrictionDestroy(&restr);
204     CeedElemRestrictionDestroy(&mesh_restr);
205     CeedElemRestrictionDestroy(&restr_i);
206   }
207 
208   /// Operator action
209   virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const {
210     CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData());
211     CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData());
212     CeedOperatorApply(oper, u, v, CEED_REQUEST_IMMEDIATE);
213     CeedVectorSyncArray(v, CEED_MEM_HOST);
214   }
215 };
216