xref: /libCEED/examples/petsc/area.c (revision 1d83af800c39d52e8338743310d4d813bfba9d04)
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 //                        libCEED + PETSc Example: Surface Area
18 //
19 // This example demonstrates a simple usage of libCEED with PETSc to calculate
20 // the surface area of a simple closed surface, such as the one of a cube or a
21 // tensor-product discrete sphere via the mass operator.
22 //
23 // The code uses higher level communication protocols in DMPlex.
24 //
25 // Build with:
26 //
27 //     make area [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
28 //
29 // Sample runs:
30 //   Sequential:
31 //
32 //     ./area -problem cube -petscspace_degree 3 -dm_refine 2
33 //
34 //     ./area -problem sphere -petscspace_degree 3 -dm_refine 2
35 //
36 //   In parallel:
37 //
38 //     mpiexec -n 4 ./area -probelm cube -petscspace_degree 3 -dm_refine 2
39 //
40 //     mpiexec -n 4 ./area -problem sphere -petscspace_degree 3 -dm_refine 2
41 //
42 //   The above example runs use 2 levels of refinement for the mesh.
43 //   Use -dm_refine k, for k levels of uniform refinement.
44 //
45 //TESTARGS -ceed {ceed_resource} -test -petscspace_degree 3
46 
47 /// @file
48 /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex
49 static const char help[] =
50   "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n";
51 
52 #include <string.h>
53 #include <petscdmplex.h>
54 #include <ceed.h>
55 #include "qfunctions/area/areacube.h"
56 #include "qfunctions/area/areasphere.h"
57 
58 #ifndef M_PI
59 #  define M_PI    3.14159265358979323846
60 #endif
61 
62 // Auxiliary function to define CEED restrictions from DMPlex data
63 static int CreateRestrictionPlex(Ceed ceed, CeedInterlaceMode imode, CeedInt P,
64                                  CeedInt ncomp, CeedElemRestriction *Erestrict,
65                                  DM dm) {
66   PetscInt ierr;
67   PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset;
68   PetscSection section;
69   Vec Uloc;
70 
71   PetscFunctionBegin;
72 
73   // Get Nelem
74   ierr = DMGetSection(dm, &section); CHKERRQ(ierr);
75   ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr);
76   nelem = cEnd - cStart;
77 
78   // Get indices
79   ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr);
80   for (c=cStart, eoffset = 0; c<cEnd; c++) {
81     PetscInt numindices, *indices, i;
82     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
83                                    &indices, NULL); CHKERRQ(ierr);
84     for (i=0; i<numindices; i+=ncomp) {
85       for (PetscInt j=0; j<ncomp; j++) {
86         if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i])))
87           SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
88                    "Cell %D closure indices not interlaced", c);
89       }
90       // NO BC on closed surfaces
91       PetscInt loc = indices[i];
92       erestrict[eoffset++] = loc/ncomp;
93     }
94     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
95                                        &indices, NULL); CHKERRQ(ierr);
96   }
97 
98   // Setup CEED restriction
99   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
100   ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr);
101 
102   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
103   CeedElemRestrictionCreate(ceed, imode, nelem, P*P, nnodes/ncomp, ncomp,
104                             CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
105                             Erestrict);
106   ierr = PetscFree(erestrict); CHKERRQ(ierr);
107 
108   PetscFunctionReturn(0);
109 }
110 
111 // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
112 static PetscErrorCode ProjectToUnitSphere(DM dm) {
113   Vec            coordinates;
114   PetscScalar   *coords;
115   PetscInt       Nv, v, dim, d;
116   PetscErrorCode ierr;
117 
118   PetscFunctionBeginUser;
119   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
120   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
121   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
122   Nv  /= dim;
123   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
124   for (v = 0; v < Nv; ++v) {
125     PetscReal r = 0.0;
126 
127     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
128     r = PetscSqrtReal(r);
129     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
130   }
131   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 int main(int argc, char **argv) {
136   PetscInt ierr;
137   MPI_Comm comm;
138   char filename[PETSC_MAX_PATH_LEN],
139        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
140   PetscInt lsize, gsize, xlsize,
141            qextra  = 1, // default number of extra quadrature points
142            ncompx  = 3, // number of components of 3D physical coordinates
143            ncompu  = 1, // dimension of field to which apply mass operator
144            topodim = 2, // topological dimension of manifold
145            degree  = 3; // default degree for finite element bases
146   PetscBool read_mesh = PETSC_FALSE,
147             test_mode = PETSC_FALSE;
148   PetscSpace sp;
149   PetscFE fe;
150   Vec X, Xloc, V, Vloc;
151   DM  dm, dmcoord;
152   Ceed ceed;
153   CeedInt P, Q;
154   CeedOperator op_setupgeo, op_apply;
155   CeedQFunction qf_setupgeo, qf_apply;
156   CeedBasis basisx, basisu;
157   CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi;
158   PetscFunctionList geomfactorlist = NULL;
159   char problemtype[PETSC_MAX_PATH_LEN] = "sphere";
160 
161   ierr = PetscInitialize(&argc, &argv, NULL, help);
162   if (ierr) return ierr;
163   comm = PETSC_COMM_WORLD;
164 
165   // Set up problem type command line option
166   ierr = PetscFunctionListAdd(&geomfactorlist, "cube", &SetupMassGeoCube);
167   CHKERRQ(ierr);
168   ierr = PetscFunctionListAdd(&geomfactorlist, "sphere", &SetupMassGeoSphere);
169   CHKERRQ(ierr);
170 
171   // Read command line options
172   ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc",
173                            NULL);
174   CHKERRQ(ierr);
175   ierr = PetscOptionsFList("-problem", "Problem to solve", NULL, geomfactorlist,
176                            problemtype, problemtype, sizeof problemtype, NULL);
177   CHKERRQ(ierr);
178   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
179                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
180   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
181                             NULL, ceedresource, ceedresource,
182                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
183   ierr = PetscOptionsBool("-test",
184                           "Testing mode (do not print unless error is large)",
185                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
186   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
187                             filename, filename, sizeof(filename), &read_mesh);
188   CHKERRQ(ierr);
189   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
190 
191   // Setup function pointer for geometric factors
192   int (*geomfp)(void *ctx, const CeedInt Q, const CeedScalar *const *in,
193                 CeedScalar *const *out);
194   ierr = PetscFunctionListFind(geomfactorlist, problemtype,
195                                (void(* *)(void))&geomfp); CHKERRQ(ierr);
196   const char *str;
197   if (geomfp == SetupMassGeoCube)
198     str = SetupMassGeoCube_loc;
199   else if (geomfp == SetupMassGeoSphere)
200     str = SetupMassGeoSphere_loc;
201   else {
202     SETERRQ(comm, PETSC_ERR_USER, "Function not found in the list");
203     return PETSC_ERR_USER;
204   }
205 
206   // Setup DM
207   if (read_mesh) {
208     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
209     CHKERRQ(ierr);
210   } else {
211     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box.
212     PetscBool simplex = PETSC_FALSE;
213     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm);
214     CHKERRQ(ierr);
215     // Set the object name
216     ierr = PetscObjectSetName((PetscObject)dm, problemtype); CHKERRQ(ierr);
217     // Distribute mesh over processes
218     {
219       DM dmDist = NULL;
220       PetscPartitioner part;
221 
222       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
223       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
224       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
225       if (dmDist) {
226         ierr = DMDestroy(&dm); CHKERRQ(ierr);
227         dm  = dmDist;
228       }
229     }
230     // Refine DMPlex with uniform refinement using runtime option -dm_refine
231     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
232     ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
233     if (!strcmp(problemtype, "sphere"))
234       ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);
235     // View DMPlex via runtime option
236     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
237   }
238 
239   // Create FE
240   ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL,
241                               PETSC_DETERMINE, &fe);
242   CHKERRQ(ierr);
243   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
244   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
245   ierr = DMCreateDS(dm); CHKERRQ(ierr);
246   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
247   CHKERRQ(ierr);
248 
249   // Get basis space degree
250   ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr);
251   ierr = PetscSpaceGetDegree(sp, &degree, NULL); CHKERRQ(ierr);
252   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
253   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
254                              "-petscspace_degree %D must be at least 1", degree);
255 
256   // Create vectors
257   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
258   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
259   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
260   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
261   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
262   ierr = VecDuplicate(X, &V); CHKERRQ(ierr);
263   ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr);
264 
265   // Set up libCEED
266   CeedInit(ceedresource, &ceed);
267 
268   // Print summary
269   P = degree + 1;
270   Q = P + qextra;
271   const char *usedresource;
272   CeedGetResource(ceed, &usedresource);
273   if (!test_mode) {
274     ierr = PetscPrintf(comm,
275                        "\n-- libCEED + PETSc Surface Area of a Manifold --\n"
276                        "  libCEED:\n"
277                        "    libCEED Backend                    : %s\n"
278                        "  Mesh:\n"
279                        "    Number of 1D Basis Nodes (p)       : %d\n"
280                        "    Number of 1D Quadrature Points (q) : %d\n"
281                        "    Global nodes                       : %D\n"
282                        "    DoF per node                       : %D\n",
283                        usedresource, P, Q,  gsize/ncompu, ncompu);
284     CHKERRQ(ierr);
285   }
286 
287   // Setup libCEED's objects:
288   // Create bases
289   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q,
290                                   CEED_GAUSS, &basisu);
291   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q,
292                                   CEED_GAUSS, &basisx);
293 
294   // CEED restrictions
295   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
296   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
297   CHKERRQ(ierr);
298 
299   CreateRestrictionPlex(ceed, CEED_INTERLACED, 2, ncompx, &Erestrictx, dmcoord);
300   CHKERRQ(ierr);
301   CreateRestrictionPlex(ceed, CEED_INTERLACED, P, ncompu, &Erestrictu, dm);
302   CHKERRQ(ierr);
303 
304   CeedInt cStart, cEnd;
305   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
306   const CeedInt nelem = cEnd - cStart;
307 
308   // CEED strided restrictions
309   const CeedInt qdatasize = 1;
310   CeedElemRestrictionCreateStrided(ceed, nelem, Q*Q, nelem*Q*Q, qdatasize,
311                                    CEED_STRIDES_BACKEND, &Erestrictqdi);
312 
313   // Element coordinates
314   Vec coords;
315   const PetscScalar *coordArray;
316   PetscSection section;
317   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
318   ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr);
319   ierr = DMGetSection(dmcoord, &section); CHKERRQ(ierr);
320 
321   CeedVector xcoord;
322   CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL);
323   CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES,
324                      (PetscScalar *)coordArray);
325   ierr = VecRestoreArrayRead(coords, &coordArray);
326 
327   // Create the vectors that will be needed in setup and apply
328   CeedVector uceed, vceed, qdata;
329   CeedInt nqpts;
330   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
331   CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata);
332   CeedVectorCreate(ceed, xlsize, &uceed);
333   CeedVectorCreate(ceed, xlsize, &vceed);
334 
335   // Create the Q-function that builds the operator for the geomteric factors
336   //   (i.e., the quadrature data)
337   CeedQFunctionCreateInterior(ceed, 1, geomfp, str, &qf_setupgeo);
338   CeedQFunctionAddInput(qf_setupgeo, "x", ncompx, CEED_EVAL_INTERP);
339   CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD);
340   CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT);
341   CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE);
342 
343   // Set up the mass operator
344   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply);
345   CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP);
346   CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE);
347   CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP);
348 
349   // Create the operator that builds the quadrature data for the operator
350   CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo);
351   CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx,
352                        CEED_VECTOR_ACTIVE);
353   CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx,
354                        CEED_VECTOR_ACTIVE);
355   CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx,
356                        CEED_VECTOR_NONE);
357   CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi,
358                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
359 
360   // Create the mass operator
361   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
362   CeedOperatorSetField(op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
363   CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED,
364                        qdata);
365   CeedOperatorSetField(op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE);
366 
367   // Compute the quadrature data for the mass operator
368   CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
369 
370   PetscScalar *v;
371   ierr = VecZeroEntries(Vloc); CHKERRQ(ierr);
372   ierr = VecGetArray(Vloc, &v);
373   CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v);
374 
375   // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1
376   if (!test_mode) {
377     ierr = PetscPrintf(comm,
378                        "Computing the mesh area using the formula: area = 1^T M 1\n");
379     CHKERRQ(ierr);
380   }
381 
382   // Initialize u with ones
383   CeedVectorSetValue(uceed, 1.0);
384 
385   // Apply the mass operator: 'u' -> 'v'
386   CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE);
387   CeedVectorSyncArray(vceed, CEED_MEM_HOST);
388 
389   // Gather output vector
390   ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr);
391   ierr = VecZeroEntries(V); CHKERRQ(ierr);
392   ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
393   ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
394 
395   // Compute and print the sum of the entries of 'v' giving the mesh surface area
396   PetscScalar area;
397   ierr = VecSum(V, &area); CHKERRQ(ierr);
398 
399   // Compute the exact surface area and print the result
400   CeedScalar exact_surfarea = 4 * M_PI;
401   if (!strcmp(problemtype, "cube")) {
402     PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
403     exact_surfarea = 6 * (2*l) * (2*l);
404   }
405 
406   if (!test_mode) {
407     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
408                        exact_surfarea); CHKERRQ(ierr);
409     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
410     CHKERRQ(ierr);
411     ierr = PetscPrintf(comm, "Area error                 : % .14g\n",
412                        fabs(area - exact_surfarea)); CHKERRQ(ierr);
413   }
414 
415   // PETSc cleanup
416   ierr = DMDestroy(&dm); CHKERRQ(ierr);
417   ierr = VecDestroy(&X); CHKERRQ(ierr);
418   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
419   ierr = VecDestroy(&V); CHKERRQ(ierr);
420   ierr = VecDestroy(&Vloc); CHKERRQ(ierr);
421 
422   // libCEED cleanup
423   CeedQFunctionDestroy(&qf_setupgeo);
424   CeedOperatorDestroy(&op_setupgeo);
425   CeedVectorDestroy(&xcoord);
426   CeedVectorDestroy(&uceed);
427   CeedVectorDestroy(&vceed);
428   CeedVectorDestroy(&qdata);
429   CeedBasisDestroy(&basisx);
430   CeedBasisDestroy(&basisu);
431   CeedElemRestrictionDestroy(&Erestrictu);
432   CeedElemRestrictionDestroy(&Erestrictx);
433   CeedElemRestrictionDestroy(&Erestrictqdi);
434   CeedQFunctionDestroy(&qf_apply);
435   CeedOperatorDestroy(&op_apply);
436   CeedDestroy(&ceed);
437   return PetscFinalize();
438 }
439