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