xref: /libCEED/examples/petsc/area.c (revision 66087c0803b66a861339d30698fcb9988ebff34d)
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, CeedInt P, CeedInt ncomp,
64                                  CeedElemRestriction *Erestrict, DM dm) {
65   PetscInt ierr;
66   PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset;
67   PetscSection section;
68   Vec Uloc;
69 
70   PetscFunctionBegin;
71 
72   // Get Nelem
73   ierr = DMGetSection(dm, &section); CHKERRQ(ierr);
74   ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr);
75   nelem = cEnd - cStart;
76 
77   // Get indices
78   ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr);
79   for (c=cStart, eoffset = 0; c<cEnd; c++) {
80     PetscInt numindices, *indices, i;
81     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
82                                    &indices, NULL); CHKERRQ(ierr);
83     for (i=0; i<numindices; i+=ncomp) {
84       for (PetscInt j=0; j<ncomp; j++) {
85         if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i])))
86           SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
87                    "Cell %D closure indices not interlaced", c);
88       }
89       // NO BC on closed surfaces
90       PetscInt loc = indices[i];
91       erestrict[eoffset++] = loc/ncomp;
92     }
93     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
94                                        &indices, NULL); CHKERRQ(ierr);
95   }
96 
97   // Setup CEED restriction
98   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
99   ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr);
100 
101   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
102   CeedElemRestrictionCreate(ceed, nelem, P*P, nnodes/ncomp, ncomp,
103                             CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
104                             Erestrict);
105   ierr = PetscFree(erestrict); CHKERRQ(ierr);
106 
107   PetscFunctionReturn(0);
108 }
109 
110 // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
111 static PetscErrorCode ProjectToUnitSphere(DM dm) {
112   Vec            coordinates;
113   PetscScalar   *coords;
114   PetscInt       Nv, v, dim, d;
115   PetscErrorCode ierr;
116 
117   PetscFunctionBeginUser;
118   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
119   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
120   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
121   Nv  /= dim;
122   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
123   for (v = 0; v < Nv; ++v) {
124     PetscReal r = 0.0;
125 
126     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
127     r = PetscSqrtReal(r);
128     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
129   }
130   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
131   PetscFunctionReturn(0);
132 }
133 
134 int main(int argc, char **argv) {
135   PetscInt ierr;
136   MPI_Comm comm;
137   char filename[PETSC_MAX_PATH_LEN],
138        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
139   PetscInt lsize, gsize, xlsize,
140            qextra  = 1, // default number of extra quadrature points
141            ncompx  = 3, // number of components of 3D physical coordinates
142            ncompu  = 1, // dimension of field to which apply mass operator
143            topodim = 2, // topological dimension of manifold
144            degree  = 3; // default degree for finite element bases
145   PetscBool read_mesh = PETSC_FALSE,
146             test_mode = PETSC_FALSE;
147   PetscSpace sp;
148   PetscFE fe;
149   Vec X, Xloc, V, Vloc;
150   DM  dm, dmcoord;
151   Ceed ceed;
152   CeedInt P, Q;
153   CeedOperator op_setupgeo, op_apply;
154   CeedQFunction qf_setupgeo, qf_apply;
155   CeedBasis basisx, basisu;
156   CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi,
157                       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     return CeedError(ceed, 1, "Function not found in the list");
203 
204   // Setup DM
205   if (read_mesh) {
206     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
207     CHKERRQ(ierr);
208   } else {
209     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box.
210     PetscBool simplex = PETSC_FALSE;
211     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm);
212     CHKERRQ(ierr);
213     // Set the object name
214     ierr = PetscObjectSetName((PetscObject)dm, problemtype); CHKERRQ(ierr);
215     // Distribute mesh over processes
216     {
217       DM dmDist = NULL;
218       PetscPartitioner part;
219 
220       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
221       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
222       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
223       if (dmDist) {
224         ierr = DMDestroy(&dm); CHKERRQ(ierr);
225         dm  = dmDist;
226       }
227     }
228     // Refine DMPlex with uniform refinement using runtime option -dm_refine
229     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
230     ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
231     if (!strcmp(problemtype, "sphere"))
232       ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);
233     // View DMPlex via runtime option
234     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
235   }
236 
237   // Create FE
238   ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL,
239                               PETSC_DETERMINE, &fe);
240   CHKERRQ(ierr);
241   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
242   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
243   ierr = DMCreateDS(dm); CHKERRQ(ierr);
244   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
245   CHKERRQ(ierr);
246 
247   // Get basis space degree
248   ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr);
249   ierr = PetscSpaceGetDegree(sp, &degree, NULL); CHKERRQ(ierr);
250   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
251   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
252                              "-petscspace_degree %D must be at least 1", degree);
253 
254   // Create vectors
255   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
256   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
257   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
258   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
259   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
260   ierr = VecDuplicate(X, &V); CHKERRQ(ierr);
261   ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr);
262 
263   // Set up libCEED
264   CeedInit(ceedresource, &ceed);
265 
266   // Print summary
267   P = degree + 1;
268   Q = P + qextra;
269   const char *usedresource;
270   CeedGetResource(ceed, &usedresource);
271   if (!test_mode) {
272     ierr = PetscPrintf(comm,
273                        "\n-- libCEED + PETSc Surface Area of a Manifold --\n"
274                        "  libCEED:\n"
275                        "    libCEED Backend                    : %s\n"
276                        "  Mesh:\n"
277                        "    Number of 1D Basis Nodes (p)       : %d\n"
278                        "    Number of 1D Quadrature Points (q) : %d\n"
279                        "    Global nodes                       : %D\n",
280                        usedresource, P, Q,  gsize/ncompu);
281     CHKERRQ(ierr);
282   }
283 
284   // Setup libCEED's objects:
285   // Create bases
286   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q,
287                                   CEED_GAUSS, &basisu);
288   CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q,
289                                   CEED_GAUSS, &basisx);
290 
291   // CEED restrictions
292   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
293   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
294   CHKERRQ(ierr);
295 
296   CreateRestrictionPlex(ceed, 2, ncompx, &Erestrictx, dmcoord); CHKERRQ(ierr);
297   CreateRestrictionPlex(ceed, P, ncompu, &Erestrictu, dm); CHKERRQ(ierr);
298 
299   CeedInt cStart, cEnd;
300   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
301   const CeedInt nelem = cEnd - cStart;
302 
303   // CEED identity restrictions
304   const CeedInt qdatasize = 1;
305   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q,
306                                     qdatasize, &Erestrictqdi);
307   CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q, 1,
308                                     &Erestrictxi);
309 
310   // Element coordinates
311   Vec coords;
312   const PetscScalar *coordArray;
313   PetscSection section;
314   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
315   ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr);
316   ierr = DMGetSection(dmcoord, &section); CHKERRQ(ierr);
317 
318   CeedVector xcoord;
319   CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL);
320   CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES,
321                      (PetscScalar *)coordArray);
322   ierr = VecRestoreArrayRead(coords, &coordArray);
323 
324   // Create the vectors that will be needed in setup and apply
325   CeedVector uceed, vceed, qdata;
326   CeedInt nqpts;
327   CeedBasisGetNumQuadraturePoints(basisu, &nqpts);
328   CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata);
329   CeedVectorCreate(ceed, xlsize, &uceed);
330   CeedVectorCreate(ceed, xlsize, &vceed);
331 
332   // Create the Q-function that builds the operator for the geomteric factors
333   //   (i.e., the quadrature data)
334   CeedQFunctionCreateInterior(ceed, 1, geomfp, str, &qf_setupgeo);
335   CeedQFunctionAddInput(qf_setupgeo, "x", ncompx, CEED_EVAL_INTERP);
336   CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD);
337   CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT);
338   CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE);
339 
340   // Set up the mass operator
341   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply);
342   CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP);
343   CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE);
344   CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP);
345 
346   // Create the operator that builds the quadrature data for the operator
347   CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo);
348   CeedOperatorSetField(op_setupgeo, "x", Erestrictx, CEED_TRANSPOSE,
349                        basisx, CEED_VECTOR_ACTIVE);
350   CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_TRANSPOSE,
351                        basisx, CEED_VECTOR_ACTIVE);
352   CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE,
353                        basisx, CEED_VECTOR_NONE);
354   CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
355                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
356 
357   // Create the mass operator
358   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
359   CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE,
360                        basisu, CEED_VECTOR_ACTIVE);
361   CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE,
362                        CEED_BASIS_COLLOCATED, qdata);
363   CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE,
364                        basisu, CEED_VECTOR_ACTIVE);
365 
366   // Compute the quadrature data for the mass operator
367   CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE);
368 
369   PetscScalar *v;
370   ierr = VecZeroEntries(Vloc); CHKERRQ(ierr);
371   ierr = VecGetArray(Vloc, &v);
372   CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v);
373 
374   // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1
375   if (!test_mode) {
376     ierr = PetscPrintf(comm,
377                        "Computing the mesh area using the formula: area = 1^T M 1\n");
378     CHKERRQ(ierr);
379   }
380 
381   // Initialize u and v with ones
382   CeedVectorSetValue(uceed, 1.0);
383 
384   // Apply the mass operator: 'u' -> 'v'
385   CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE);
386   CeedVectorSyncArray(vceed, CEED_MEM_HOST);
387 
388   // Gather output vector
389   ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr);
390   ierr = VecZeroEntries(V); CHKERRQ(ierr);
391   ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
392   ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
393 
394   // Compute and print the sum of the entries of 'v' giving the mesh surface area
395   PetscScalar area;
396   ierr = VecSum(V, &area); CHKERRQ(ierr);
397 
398   // Compute the exact surface area and print the result
399   CeedScalar exact_surfarea = 4 * M_PI;
400   if (!strcmp(problemtype, "cube")) {
401     PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
402     exact_surfarea = 6 * (2*l) * (2*l);
403   }
404 
405   if (!test_mode) {
406     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
407                        exact_surfarea); CHKERRQ(ierr);
408     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
409     CHKERRQ(ierr);
410     ierr = PetscPrintf(comm, "Area error                 : % .14g\n",
411                        fabs(area - exact_surfarea)); CHKERRQ(ierr);
412   }
413 
414   // PETSc cleanup
415   ierr = DMDestroy(&dm); CHKERRQ(ierr);
416   ierr = VecDestroy(&X); CHKERRQ(ierr);
417   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
418   ierr = VecDestroy(&V); CHKERRQ(ierr);
419   ierr = VecDestroy(&Vloc); CHKERRQ(ierr);
420 
421   // libCEED cleanup
422   CeedQFunctionDestroy(&qf_setupgeo);
423   CeedOperatorDestroy(&op_setupgeo);
424   CeedVectorDestroy(&xcoord);
425   CeedVectorDestroy(&uceed);
426   CeedVectorDestroy(&vceed);
427   CeedVectorDestroy(&qdata);
428   CeedBasisDestroy(&basisx);
429   CeedBasisDestroy(&basisu);
430   CeedElemRestrictionDestroy(&Erestrictu);
431   CeedElemRestrictionDestroy(&Erestrictx);
432   CeedElemRestrictionDestroy(&Erestrictxi);
433   CeedElemRestrictionDestroy(&Erestrictqdi);
434   CeedQFunctionDestroy(&qf_apply);
435   CeedOperatorDestroy(&op_apply);
436   CeedDestroy(&ceed);
437   return PetscFinalize();
438 }
439