xref: /libCEED/examples/petsc/area.c (revision 22785960d60f45e0c65c47523fd8ab76bf9537ed)
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        ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
64   PetscInt l_size, g_size, xl_size,
65            q_extra     = 1, // default number of extra quadrature points
66            num_comp_x  = 3, // number of components of 3D physical coordinates
67            num_comp_u  = 1, // dimension of field to which apply mass operator
68            topo_dim    = 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, U_loc, V, V_loc;
74   DM  dm;
75   UserO user;
76   Ceed ceed;
77   CeedData ceed_data;
78   ProblemType problem_choice;
79   VecType vec_type;
80   PetscMemType mem_type;
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   problem_choice = SPHERE;
91   ierr = PetscOptionsEnum("-problem",
92                           "Problem to solve", NULL,
93                           problem_types, (PetscEnum)problem_choice,
94                           (PetscEnum *)&problem_choice,
95                           NULL); CHKERRQ(ierr);
96   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
97                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
98   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
99                             NULL, ceed_resource, ceed_resource,
100                             sizeof(ceed_resource), 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, topo_dim, simplex, 1., &dm);
120     CHKERRQ(ierr);
121     // Set the object name
122     ierr = PetscObjectSetName((PetscObject)dm, problem_types[problem_choice]);
123     CHKERRQ(ierr);
124     // Distribute mesh over processes
125     {
126       DM dm_dist = NULL;
127       PetscPartitioner part;
128 
129       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
130       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
131       ierr = DMPlexDistribute(dm, 0, NULL, &dm_dist); CHKERRQ(ierr);
132       if (dm_dist) {
133         ierr = DMDestroy(&dm); CHKERRQ(ierr);
134         dm  = dm_dist;
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 (problem_choice == 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, num_comp_u, topo_dim, false,
149                          (BCFunction)NULL);
150   CHKERRQ(ierr);
151 
152   // Create vectors
153   ierr = DMCreateGlobalVector(dm, &U); CHKERRQ(ierr);
154   ierr = VecGetLocalSize(U, &l_size); CHKERRQ(ierr);
155   ierr = VecGetSize(U, &g_size); CHKERRQ(ierr);
156   ierr = DMCreateLocalVector(dm, &U_loc); CHKERRQ(ierr);
157   ierr = VecGetSize(U_loc, &xl_size); CHKERRQ(ierr);
158   ierr = VecDuplicate(U, &V); CHKERRQ(ierr);
159   ierr = VecDuplicate(U_loc, &V_loc); CHKERRQ(ierr);
160 
161   // Setup user structure
162   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
163 
164   // Set up libCEED
165   CeedInit(ceed_resource, &ceed);
166   CeedMemType mem_type_backend;
167   CeedGetPreferredMemType(ceed, &mem_type_backend);
168 
169   ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
170   if (!vec_type) { // Not yet set by user -dm_vec_type
171     switch (mem_type_backend) {
172     case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
173     case CEED_MEM_DEVICE: {
174       const char *resolved;
175       CeedGetResource(ceed, &resolved);
176       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
177       else if (strstr(resolved, "/gpu/hip/occa"))
178         vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
179       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
180       else vec_type = VECSTANDARD;
181     }
182     }
183     ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr);
184   }
185 
186   // Print summary
187   if (!test_mode) {
188     PetscInt P = degree + 1, Q = P + q_extra;
189     const char *used_resource;
190     CeedGetResource(ceed, &used_resource);
191     ierr = PetscPrintf(comm,
192                        "\n-- libCEED + PETSc Surface Area of a Manifold --\n"
193                        "  libCEED:\n"
194                        "    libCEED Backend                    : %s\n"
195                        "    libCEED Backend MemType            : %s\n"
196                        "  Mesh:\n"
197                        "    Number of 1D Basis Nodes (p)       : %d\n"
198                        "    Number of 1D Quadrature Points (q) : %d\n"
199                        "    Global nodes                       : %D\n"
200                        "    DoF per node                       : %D\n"
201                        "    Global DoFs                        : %D\n",
202                        used_resource, CeedMemTypes[mem_type_backend], P, Q,
203                        g_size/num_comp_u, num_comp_u, g_size); CHKERRQ(ierr);
204   }
205 
206   // Setup libCEED's objects and apply setup operator
207   ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr);
208   ierr = SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x,
209                               num_comp_u,
210                               g_size, xl_size, problem_options[problem_choice], ceed_data,
211                               false, (CeedVector)NULL, (CeedVector *)NULL);
212   CHKERRQ(ierr);
213 
214   // Setup output vector
215   PetscScalar *v;
216   ierr = VecZeroEntries(V_loc); CHKERRQ(ierr);
217   ierr = VecGetArrayAndMemType(V_loc, &v, &mem_type); CHKERRQ(ierr);
218   CeedVectorSetArray(ceed_data->y_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER,
219                      v);
220 
221   // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1
222   if (!test_mode) {
223     ierr = PetscPrintf(comm,
224                        "Computing the mesh area using the formula: area = 1^T M 1\n");
225     CHKERRQ(ierr);
226   }
227 
228   // Initialize u with ones
229   CeedVectorSetValue(ceed_data->x_ceed, 1.0);
230 
231   // Apply the mass operator: 'u' -> 'v'
232   CeedOperatorApply(ceed_data->op_apply, ceed_data->x_ceed, ceed_data->y_ceed,
233                     CEED_REQUEST_IMMEDIATE);
234 
235   // Gather output vector
236   CeedVectorTakeArray(ceed_data->y_ceed, CEED_MEM_HOST, NULL);
237   ierr = VecRestoreArrayAndMemType(V_loc, &v); CHKERRQ(ierr);
238   ierr = VecZeroEntries(V); CHKERRQ(ierr);
239   ierr = DMLocalToGlobalBegin(dm, V_loc, ADD_VALUES, V); CHKERRQ(ierr);
240   ierr = DMLocalToGlobalEnd(dm, V_loc, ADD_VALUES, V); CHKERRQ(ierr);
241 
242   // Compute and print the sum of the entries of 'v' giving the mesh surface area
243   PetscScalar area;
244   ierr = VecSum(V, &area); CHKERRQ(ierr);
245 
246   // Compute the exact surface area and print the result
247   CeedScalar exact_surface_area = 4 * M_PI;
248   if (problem_choice == CUBE) {
249     PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube
250     exact_surface_area = 6 * (2*l) * (2*l);
251   }
252 
253   PetscReal error = fabs(area - exact_surface_area);
254   PetscReal tol = 5e-6;
255   if (!test_mode || error > tol) {
256     ierr = PetscPrintf(comm, "Exact mesh surface area    : % .14g\n",
257                        exact_surface_area);
258     CHKERRQ(ierr);
259     ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area);
260     CHKERRQ(ierr);
261     ierr = PetscPrintf(comm, "Area error                 : % .14g\n", error);
262     CHKERRQ(ierr);
263   }
264 
265   // Cleanup
266   ierr = DMDestroy(&dm); CHKERRQ(ierr);
267   ierr = VecDestroy(&U); CHKERRQ(ierr);
268   ierr = VecDestroy(&U_loc); CHKERRQ(ierr);
269   ierr = VecDestroy(&V); CHKERRQ(ierr);
270   ierr = VecDestroy(&V_loc); CHKERRQ(ierr);
271   ierr = PetscFree(user); CHKERRQ(ierr);
272   ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr);
273   CeedDestroy(&ceed);
274   return PetscFinalize();
275 }
276