xref: /libCEED/examples/petsc/area.c (revision e83e87a50f1b2e8810b376a735430565127e4d25) !
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 -degree 3 -dm_refine 2
33 //     ./area -problem sphere -degree 3 -dm_refine 2
34 //
35 //   In parallel:
36 //
37 //     mpiexec -n 4 ./area -problem cube -degree 3 -dm_refine 2
38 //     mpiexec -n 4 ./area -problem sphere -degree 3 -dm_refine 2
39 //
40 //   The above example runs use 2 levels of refinement for the mesh.
41 //   Use -dm_refine k, for k levels of uniform refinement.
42 //
43 //TESTARGS -ceed {ceed_resource} -test -degree 3 -dm_refine 1
44 
45 /// @file
46 /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex
47 static const char help[] =
48   "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n";
49 
50 #include <ceed.h>
51 #include <petscdmplex.h>
52 #include <string.h>
53 #include "area.h"
54 
55 #ifndef M_PI
56 #  define M_PI 3.14159265358979323846
57 #endif
58 
59 int main(int argc, char **argv) {
60   PetscInt ierr;
61   MPI_Comm comm;
62   char filename[PETSC_MAX_PATH_LEN],
63        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
64   PetscInt lsize, gsize, xlsize,
65            qextra  = 1, // default number of extra quadrature points
66            ncompx  = 3, // number of components of 3D physical coordinates
67            ncompu  = 1, // dimension of field to which apply mass operator
68            topodim = 2, // topological dimension of manifold
69            degree  = 3; // default degree for finite element bases
70   PetscBool read_mesh = PETSC_FALSE,
71             test_mode = PETSC_FALSE,
72             simplex = PETSC_FALSE;
73   Vec U, Uloc, V, Vloc;
74   DM  dm;
75   UserO user;
76   Ceed ceed;
77   CeedData ceeddata;
78   problemType problemChoice;
79   VecType vectype;
80   PetscMemType memtype;
81 
82   ierr = PetscInitialize(&argc, &argv, NULL, help);
83   if (ierr) return ierr;
84   comm = PETSC_COMM_WORLD;
85 
86   // Read command line options
87   ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc",
88                            NULL);
89   CHKERRQ(ierr);
90   problemChoice = SPHERE;
91   ierr = PetscOptionsEnum("-problem",
92                           "Problem to solve", NULL,
93                           problemTypes, (PetscEnum)problemChoice,
94                           (PetscEnum *)&problemChoice,
95                           NULL); CHKERRQ(ierr);
96   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
97                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
98   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
99                             NULL, ceedresource, ceedresource,
100                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
101   ierr = PetscOptionsBool("-test",
102                           "Testing mode (do not print unless error is large)",
103                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
104   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
105                             filename, filename, sizeof(filename), &read_mesh);
106   CHKERRQ(ierr);
107   ierr = PetscOptionsBool("-simplex", "Use simplices, or tensor product cells",
108                           NULL, simplex, &simplex, NULL); CHKERRQ(ierr);
109   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
110                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
111   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
112 
113   // Setup DM
114   if (read_mesh) {
115     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
116     CHKERRQ(ierr);
117   } else {
118     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box
119     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, 1., &dm);
120     CHKERRQ(ierr);
121     // Set the object name
122     ierr = PetscObjectSetName((PetscObject)dm, problemTypes[problemChoice]);
123     CHKERRQ(ierr);
124     // Distribute mesh over processes
125     {
126       DM dmDist = NULL;
127       PetscPartitioner part;
128 
129       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
130       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
131       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
132       if (dmDist) {
133         ierr = DMDestroy(&dm); CHKERRQ(ierr);
134         dm  = dmDist;
135       }
136     }
137     // Refine DMPlex with uniform refinement using runtime option -dm_refine
138     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
139     ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
140     if (problemChoice == SPHERE) {
141       ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);
142     }
143     // View DMPlex via runtime option
144     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
145   }
146 
147   // Create DM
148   ierr = SetupDMByDegree(dm, degree, ncompu, topodim, false, (BCFunction)NULL);
149   CHKERRQ(ierr);
150 
151   // Create vectors
152   ierr = DMCreateGlobalVector(dm, &U); CHKERRQ(ierr);
153   ierr = VecGetLocalSize(U, &lsize); CHKERRQ(ierr);
154   ierr = VecGetSize(U, &gsize); CHKERRQ(ierr);
155   ierr = DMCreateLocalVector(dm, &Uloc); CHKERRQ(ierr);
156   ierr = VecGetSize(Uloc, &xlsize); CHKERRQ(ierr);
157   ierr = VecDuplicate(U, &V); CHKERRQ(ierr);
158   ierr = VecDuplicate(Uloc, &Vloc); CHKERRQ(ierr);
159 
160   // Setup user structure
161   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
162 
163   // Set up libCEED
164   CeedInit(ceedresource, &ceed);
165   CeedMemType memtypebackend;
166   CeedGetPreferredMemType(ceed, &memtypebackend);
167 
168   ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr);
169   if (!vectype) { // Not yet set by user -dm_vec_type
170     switch (memtypebackend) {
171     case CEED_MEM_HOST: vectype = VECSTANDARD; break;
172     case CEED_MEM_DEVICE: {
173       const char *resolved;
174       CeedGetResource(ceed, &resolved);
175       if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
176       else if (strstr(resolved, "/gpu/hip/occa"))
177         vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
178       else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
179       else vectype = VECSTANDARD;
180     }
181     }
182     ierr = DMSetVecType(dm, vectype); CHKERRQ(ierr);
183   }
184 
185   // Print summary
186   if (!test_mode) {
187     PetscInt P = degree + 1, Q = P + qextra;
188     const char *usedresource;
189     CeedGetResource(ceed, &usedresource);
190     ierr = PetscPrintf(comm,
191                        "\n-- libCEED + PETSc Surface Area of a Manifold --\n"
192                        "  libCEED:\n"
193                        "    libCEED Backend                    : %s\n"
194                        "    libCEED Backend MemType            : %s\n"
195                        "  Mesh:\n"
196                        "    Number of 1D Basis Nodes (p)       : %d\n"
197                        "    Number of 1D Quadrature Points (q) : %d\n"
198                        "    Global nodes                       : %D\n"
199                        "    DoF per node                       : %D\n"
200                        "    Global DoFs                        : %D\n",
201                        usedresource, CeedMemTypes[memtypebackend], P, Q,
202                        gsize/ncompu, ncompu, gsize); CHKERRQ(ierr);
203   }
204 
205   // Setup libCEED's objects and apply setup operator
206   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
207   ierr = SetupLibceedByDegree(dm, ceed, degree, topodim, qextra, ncompx, ncompu,
208                               gsize, xlsize, problemOptions[problemChoice], ceeddata,
209                               false, (CeedVector)NULL, (CeedVector *)NULL);
210   CHKERRQ(ierr);
211 
212   // Setup output vector
213   PetscScalar *v;
214   ierr = VecZeroEntries(Vloc); CHKERRQ(ierr);
215   ierr = VecGetArrayAndMemType(Vloc, &v, &memtype); CHKERRQ(ierr);
216   CeedVectorSetArray(ceeddata->Yceed, MemTypeP2C(memtype), CEED_USE_POINTER, v);
217 
218   // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1
219   if (!test_mode) {
220     ierr = PetscPrintf(comm,
221                        "Computing the mesh area using the formula: area = 1^T M 1\n");
222     CHKERRQ(ierr);
223   }
224 
225   // Initialize u with ones
226   CeedVectorSetValue(ceeddata->Xceed, 1.0);
227 
228   // Apply the mass operator: 'u' -> 'v'
229   CeedOperatorApply(ceeddata->opApply, ceeddata->Xceed, ceeddata->Yceed,
230                     CEED_REQUEST_IMMEDIATE);
231 
232   // Gather output vector
233   CeedVectorTakeArray(ceeddata->Yceed, CEED_MEM_HOST, NULL);
234   ierr = VecRestoreArrayAndMemType(Vloc, &v); CHKERRQ(ierr);
235   ierr = VecZeroEntries(V); CHKERRQ(ierr);
236   ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
237   ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr);
238 
239   // Compute and print the sum of the entries of 'v' giving the mesh surface area
240   PetscScalar area;
241   ierr = VecSum(V, &area); CHKERRQ(ierr);
242 
243   // Compute the exact surface area and print the result
244   CeedScalar exact_surfarea = 4 * M_PI;
245   if (problemChoice == CUBE) {
246     PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
247     exact_surfarea = 6 * (2*l) * (2*l);
248   }
249 
250   PetscReal error = fabs(area - exact_surfarea);
251   PetscReal tol = 5e-6;
252   if (!test_mode || error > tol) {
253     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
254                        exact_surfarea);
255     CHKERRQ(ierr);
256     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
257     CHKERRQ(ierr);
258     ierr = PetscPrintf(comm, "Area error                 : % .14g\n", error);
259     CHKERRQ(ierr);
260   }
261 
262   // Cleanup
263   ierr = DMDestroy(&dm); CHKERRQ(ierr);
264   ierr = VecDestroy(&U); CHKERRQ(ierr);
265   ierr = VecDestroy(&Uloc); CHKERRQ(ierr);
266   ierr = VecDestroy(&V); CHKERRQ(ierr);
267   ierr = VecDestroy(&Vloc); CHKERRQ(ierr);
268   ierr = PetscFree(user); CHKERRQ(ierr);
269   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
270   CeedDestroy(&ceed);
271   return PetscFinalize();
272 }
273