xref: /libCEED/examples/mfem/bp3.hpp (revision 21617c043cf2be06120e43f272b68c5ebd5f09fa)
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 /// MFEM diffusion operator based on libCEED
19 
20 #include <ceed.h>
21 #include <mfem.hpp>
22 
23 /// A structure used to pass additional data to f_build_diff and f_apply_diff
24 struct DiffContext { CeedInt dim, space_dim; };
25 
26 /// libCEED Q-function for building quadrature data for a diffusion operator
27 static int f_build_diff(void *ctx, void *qdata, CeedInt Q,
28                         const CeedScalar *const *u, CeedScalar *const *v) {
29   // u[1] is Jacobians, size (Q x nc x dim) with column-major layout
30   // u[4] is quadrature weights, size (Q)
31   //
32   // At every quadrature point, compute qw/det(J).adj(J).adj(J)^T and store
33   // the symmetric part of the result.
34   DiffContext *dc = (DiffContext*)ctx;
35   CeedScalar *qd = (CeedScalar*)qdata;
36   const CeedScalar *J = u[1], *qw = u[4];
37   switch (dc->dim + 10*dc->space_dim) {
38   case 11:
39     for (CeedInt i=0; i<Q; i++) {
40       qd[i] = qw[i] / J[i];
41     }
42     break;
43   case 22:
44     for (CeedInt i=0; i<Q; i++) {
45       // J: 0 2   qd: 0 1   adj(J):  J22 -J12
46       //    1 3       1 2           -J21  J11
47       const CeedScalar J11 = J[i+Q*0];
48       const CeedScalar J21 = J[i+Q*1];
49       const CeedScalar J12 = J[i+Q*2];
50       const CeedScalar J22 = J[i+Q*3];
51       const CeedScalar w = qw[i] / (J11*J22 - J21*J12);
52       qd[i+Q*0] =   w * (J12*J12 + J22*J22);
53       qd[i+Q*1] = - w * (J11*J12 + J21*J22);
54       qd[i+Q*2] =   w * (J11*J11 + J21*J21);
55     }
56     break;
57   case 33:
58     for (CeedInt i=0; i<Q; i++) {
59       // J: 0 3 6   qd: 0 1 2
60       //    1 4 7       1 3 4
61       //    2 5 8       2 4 5
62       const CeedScalar J11 = J[i+Q*0];
63       const CeedScalar J21 = J[i+Q*1];
64       const CeedScalar J31 = J[i+Q*2];
65       const CeedScalar J12 = J[i+Q*3];
66       const CeedScalar J22 = J[i+Q*4];
67       const CeedScalar J32 = J[i+Q*5];
68       const CeedScalar J13 = J[i+Q*6];
69       const CeedScalar J23 = J[i+Q*7];
70       const CeedScalar J33 = J[i+Q*8];
71       const CeedScalar A11 = J22*J33 - J23*J32;
72       const CeedScalar A12 = J13*J32 - J12*J33;
73       const CeedScalar A13 = J12*J23 - J13*J22;
74       const CeedScalar A21 = J23*J31 - J21*J33;
75       const CeedScalar A22 = J11*J33 - J13*J31;
76       const CeedScalar A23 = J13*J21 - J11*J23;
77       const CeedScalar A31 = J21*J32 - J22*J31;
78       const CeedScalar A32 = J12*J31 - J11*J32;
79       const CeedScalar A33 = J11*J22 - J12*J21;
80       const CeedScalar w = qw[i] / (J11*A11 + J21*A12 + J31*A13);
81       qd[i+Q*0] = w * (A11*A11 + A12*A12 + A13*A13);
82       qd[i+Q*1] = w * (A11*A21 + A12*A22 + A13*A23);
83       qd[i+Q*2] = w * (A11*A31 + A12*A32 + A13*A33);
84       qd[i+Q*3] = w * (A21*A21 + A22*A22 + A23*A23);
85       qd[i+Q*4] = w * (A21*A31 + A22*A32 + A23*A33);
86       qd[i+Q*5] = w * (A31*A31 + A32*A32 + A33*A33);
87     }
88     break;
89   default:
90     return CeedError(NULL, 1, "dim=%d, space_dim=%d is not supported",
91                      dc->dim, dc->space_dim);
92   }
93   return 0;
94 }
95 
96 /// libCEED Q-function for applying a diff operator
97 static int f_apply_diff(void *ctx, void *qdata, CeedInt Q,
98                         const CeedScalar *const *u, CeedScalar *const *v) {
99   DiffContext *dc = (DiffContext*)ctx;
100   const CeedScalar *qd = (const CeedScalar*)qdata;
101   // u[1], v[1]: size: (Q x nc x dim) with column-major layout (nc == 1)
102   const CeedScalar *ug = u[1];
103   CeedScalar *vg = v[1];
104   switch (dc->dim) {
105   case 1:
106     for (CeedInt i=0; i<Q; i++) {
107       vg[i] = ug[i] * qd[i];
108     }
109     break;
110   case 2:
111     for (CeedInt i=0; i<Q; i++) {
112       const CeedScalar ug0 = ug[i+Q*0];
113       const CeedScalar ug1 = ug[i+Q*1];
114       vg[i+Q*0] = qd[i+Q*0]*ug0 + qd[i+Q*1]*ug1;
115       vg[i+Q*1] = qd[i+Q*1]*ug0 + qd[i+Q*2]*ug1;
116     }
117     break;
118   case 3:
119     for (CeedInt i=0; i<Q; i++) {
120       const CeedScalar ug0 = ug[i+Q*0];
121       const CeedScalar ug1 = ug[i+Q*1];
122       const CeedScalar ug2 = ug[i+Q*2];
123       vg[i+Q*0] = qd[i+Q*0]*ug0 + qd[i+Q*1]*ug1 + qd[i+Q*2]*ug2;
124       vg[i+Q*1] = qd[i+Q*1]*ug0 + qd[i+Q*3]*ug1 + qd[i+Q*4]*ug2;
125       vg[i+Q*2] = qd[i+Q*2]*ug0 + qd[i+Q*4]*ug1 + qd[i+Q*5]*ug2;
126     }
127     break;
128   default:
129     return CeedError(NULL, 1, "topo_dim=%d is not supported", dc->dim);
130   }
131   return 0;
132 }
133 
134 /// Wrapper for a diffusion CeedOperator as an mfem::Operator
135 class CeedDiffusionOperator : public mfem::Operator {
136  protected:
137   const mfem::FiniteElementSpace *fes;
138   CeedOperator build_oper, oper;
139   CeedBasis basis, mesh_basis;
140   CeedElemRestriction restr, mesh_restr;
141   CeedQFunction apply_qfunc, build_qfunc;
142   CeedVector node_coords, qdata;
143 
144   DiffContext diff_ctx;
145 
146   CeedVector u, v;
147 
148   static void FESpace2Ceed(const mfem::FiniteElementSpace *fes,
149                            const mfem::IntegrationRule &ir,
150                            Ceed ceed, CeedBasis *basis,
151                            CeedElemRestriction *restr) {
152     mfem::Mesh *mesh = fes->GetMesh();
153     const mfem::FiniteElement *fe = fes->GetFE(0);
154     const int order = fes->GetOrder(0);
155     mfem::Array<int> dof_map;
156     switch (mesh->Dimension()) {
157     case 1: {
158       const mfem::H1_SegmentElement *h1_fe =
159         dynamic_cast<const mfem::H1_SegmentElement*>(fe);
160       MFEM_VERIFY(h1_fe, "invalid FE");
161       h1_fe->GetDofMap().Copy(dof_map);
162       break;
163     }
164     case 2: {
165       const mfem::H1_QuadrilateralElement *h1_fe =
166         dynamic_cast<const mfem::H1_QuadrilateralElement*>(fe);
167       MFEM_VERIFY(h1_fe, "invalid FE");
168       h1_fe->GetDofMap().Copy(dof_map);
169       break;
170     }
171     case 3: {
172       const mfem::H1_HexahedronElement *h1_fe =
173         dynamic_cast<const mfem::H1_HexahedronElement*>(fe);
174       MFEM_VERIFY(h1_fe, "invalid FE");
175       h1_fe->GetDofMap().Copy(dof_map);
176       break;
177     }
178     }
179     const mfem::FiniteElement *fe1d =
180       fes->FEColl()->FiniteElementForGeometry(mfem::Geometry::SEGMENT);
181     mfem::DenseMatrix shape1d(fe1d->GetDof(), ir.GetNPoints());
182     mfem::DenseMatrix grad1d(fe1d->GetDof(), ir.GetNPoints());
183     mfem::Vector qref1d(ir.GetNPoints()), qweight1d(ir.GetNPoints());
184     mfem::Vector shape_i(shape1d.Height());
185     mfem::DenseMatrix grad_i(grad1d.Height(), 1);
186     const mfem::H1_SegmentElement *h1_fe1d =
187       dynamic_cast<const mfem::H1_SegmentElement*>(fe1d);
188     MFEM_VERIFY(h1_fe1d, "invalid FE");
189     const mfem::Array<int> &dof_map_1d = h1_fe1d->GetDofMap();
190     for (int i = 0; i < ir.GetNPoints(); i++) {
191       const mfem::IntegrationPoint &ip = ir.IntPoint(i);
192       qref1d(i) = ip.x;
193       qweight1d(i) = ip.weight;
194       fe1d->CalcShape(ip, shape_i);
195       fe1d->CalcDShape(ip, grad_i);
196       for (int j = 0; j < shape1d.Height(); j++) {
197         shape1d(j,i) = shape_i(dof_map_1d[j]);
198         grad1d(j,i) = grad_i(dof_map_1d[j],0);
199       }
200     }
201     CeedBasisCreateTensorH1(ceed, mesh->Dimension(), fes->GetVDim(), order+1,
202                             ir.GetNPoints(), shape1d.GetData(),
203                             grad1d.GetData(), qref1d.GetData(),
204                             qweight1d.GetData(), basis);
205 
206     const mfem::Table &el_dof = fes->GetElementToDofTable();
207     mfem::Array<int> tp_el_dof(el_dof.Size_of_connections());
208     for (int i = 0; i < mesh->GetNE(); i++) {
209       const int el_offset = fe->GetDof()*i;
210       for (int j = 0; j < fe->GetDof(); j++) {
211         tp_el_dof[j + el_offset] = el_dof.GetJ()[dof_map[j] + el_offset];
212       }
213     }
214     CeedElemRestrictionCreate(ceed, mesh->GetNE(), fe->GetDof(),
215                               fes->GetNDofs(), CEED_MEM_HOST, CEED_COPY_VALUES,
216                               tp_el_dof.GetData(), restr);
217   }
218 
219  public:
220   /// Constructor. Assumes @a fes is a scalar FE space.
221   CeedDiffusionOperator(Ceed ceed, const mfem::FiniteElementSpace *fes)
222     : Operator(fes->GetNDofs()),
223       fes(fes) {
224     mfem::Mesh *mesh = fes->GetMesh();
225     const int order = fes->GetOrder(0);
226     const int ir_order = 2*(order + 2) - 1; // <-----
227     const mfem::IntegrationRule &ir =
228       mfem::IntRules.Get(mfem::Geometry::SEGMENT, ir_order);
229 
230     FESpace2Ceed(fes, ir, ceed, &basis, &restr);
231 
232     const mfem::FiniteElementSpace *mesh_fes = mesh->GetNodalFESpace();
233     MFEM_VERIFY(mesh_fes, "the Mesh has no nodal FE space");
234     FESpace2Ceed(mesh_fes, ir, ceed, &mesh_basis, &mesh_restr);
235 
236     CeedVectorCreate(ceed, mesh->GetNodes()->Size(), &node_coords);
237     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER,
238                        mesh->GetNodes()->GetData());
239 
240     const int dim = mesh->Dimension();
241     diff_ctx.dim = dim;
242     diff_ctx.space_dim = mesh->SpaceDimension();
243 
244     const int qsize = dim*(dim+1)/2;
245     CeedQFunctionCreateInterior(ceed, 1, 1, qsize*sizeof(CeedScalar),
246                                 (CeedEvalMode)(CEED_EVAL_GRAD|CEED_EVAL_WEIGHT),
247                                 CEED_EVAL_NONE, f_build_diff,
248                                 __FILE__":f_build_diff", &build_qfunc);
249     CeedQFunctionSetContext(build_qfunc, &diff_ctx, sizeof(diff_ctx));
250     CeedOperatorCreate(ceed, mesh_restr, mesh_basis, build_qfunc, NULL, NULL,
251                        &build_oper);
252     CeedOperatorGetQData(build_oper, &qdata);
253     CeedOperatorApply(build_oper, qdata, node_coords, NULL,
254                       CEED_REQUEST_IMMEDIATE);
255 
256     CeedQFunctionCreateInterior(ceed, 1, 1, qsize*sizeof(CeedScalar),
257                                 CEED_EVAL_GRAD, CEED_EVAL_GRAD, f_apply_diff,
258                                 __FILE__":f_apply_diff", &apply_qfunc);
259     CeedQFunctionSetContext(apply_qfunc, &diff_ctx, sizeof(diff_ctx));
260     CeedOperatorCreate(ceed, restr, basis, apply_qfunc, NULL, NULL, &oper);
261 
262     CeedVectorCreate(ceed, fes->GetNDofs(), &u);
263     CeedVectorCreate(ceed, fes->GetNDofs(), &v);
264   }
265 
266   /// Destructor
267   ~CeedDiffusionOperator() {
268     CeedVectorDestroy(&v);
269     CeedVectorDestroy(&u);
270     CeedOperatorDestroy(&oper);
271     CeedQFunctionDestroy(&apply_qfunc);
272     // qdata is owned by build_oper
273     CeedOperatorDestroy(&build_oper);
274     CeedQFunctionDestroy(&build_qfunc);
275     CeedVectorDestroy(&node_coords);
276     CeedElemRestrictionDestroy(&mesh_restr);
277     CeedBasisDestroy(&mesh_basis);
278     CeedElemRestrictionDestroy(&restr);
279     CeedBasisDestroy(&basis);
280   }
281 
282   /// Operator action
283   virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const {
284     CeedVectorSetArray(u, CEED_MEM_HOST, CEED_USE_POINTER, x.GetData());
285     CeedVectorSetArray(v, CEED_MEM_HOST, CEED_USE_POINTER, y.GetData());
286 
287     CeedOperatorApply(oper, qdata, u, v, CEED_REQUEST_IMMEDIATE);
288   }
289 };
290