xref: /libCEED/examples/petsc/bpssphere.c (revision 7d8d0e25636a94a27ff75b3dec09737e24cdb0fe)
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: CEED BPs
18 //
19 // This example demonstrates a simple usage of libCEED with PETSc to solve the
20 // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps,
21 // on a closed surface, such as the one of a discrete sphere.
22 //
23 // The code uses higher level communication protocols in DMPlex.
24 //
25 // Build with:
26 //
27 //     make bpssphere [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
28 //
29 // Sample runs:
30 //
31 //     bpssphere -problem bp1 -degree 3
32 //     bpssphere -problem bp2 -degree 3
33 //     bpssphere -problem bp3 -degree 3
34 //     bpssphere -problem bp4 -degree 3
35 //     bpssphere -problem bp5 -degree 3 -ceed /cpu/self
36 //     bpssphere -problem bp6 -degree 3 -ceed /gpu/cuda
37 //
38 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -dm_refine 2
39 
40 /// @file
41 /// CEED BPs example using PETSc with DMPlex
42 /// See bps.c for a "raw" implementation using a structured grid.
43 /// and bpsdmplex.c for an implementation using an unstructured grid.
44 static const char help[] = "Solve CEED BPs on a sphere using DMPlex in PETSc\n";
45 
46 #include <stdbool.h>
47 #include <string.h>
48 #include <petscksp.h>
49 #include <petscdmplex.h>
50 #include <ceed.h>
51 #include "setupsphere.h"
52 
53 int main(int argc, char **argv) {
54   PetscInt ierr;
55   MPI_Comm comm;
56   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self",
57                                           filename[PETSC_MAX_PATH_LEN];
58   double my_rt_start, my_rt, rt_min, rt_max;
59   PetscInt degree = 3, qextra, lsize, gsize, topodim = 2, ncompx = 3,
60            ncompu = 1, xlsize;
61   PetscScalar *r;
62   PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex;
63   PetscLogStage solvestage;
64   Vec X, Xloc, rhs, rhsloc;
65   Mat matO;
66   KSP ksp;
67   DM  dm;
68   UserO userO;
69   Ceed ceed;
70   CeedData ceeddata;
71   CeedQFunction qf_error;
72   CeedOperator op_error;
73   CeedVector rhsceed, target;
74   bpType bpChoice;
75 
76   ierr = PetscInitialize(&argc, &argv, NULL, help);
77   if (ierr) return ierr;
78   comm = PETSC_COMM_WORLD;
79 
80   // Read command line options
81   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
82   bpChoice = CEED_BP1;
83   ierr = PetscOptionsEnum("-problem",
84                           "CEED benchmark problem to solve", NULL,
85                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
86                           NULL); CHKERRQ(ierr);
87   ncompu = bpOptions[bpChoice].ncompu;
88   test_mode = PETSC_FALSE;
89   ierr = PetscOptionsBool("-test",
90                           "Testing mode (do not print unless error is large)",
91                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
92   benchmark_mode = PETSC_FALSE;
93   ierr = PetscOptionsBool("-benchmark",
94                           "Benchmarking mode (prints benchmark statistics)",
95                           NULL, benchmark_mode, &benchmark_mode, NULL);
96   CHKERRQ(ierr);
97   write_solution = PETSC_FALSE;
98   ierr = PetscOptionsBool("-write_solution",
99                           "Write solution for visualization",
100                           NULL, write_solution, &write_solution, NULL);
101   CHKERRQ(ierr);
102   degree = test_mode ? 3 : 2;
103   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
104                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
105   qextra = bpOptions[bpChoice].qextra;
106   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
107                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
108   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
109                             NULL, ceedresource, ceedresource,
110                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
111   read_mesh = PETSC_FALSE;
112   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
113                             filename, filename, sizeof(filename), &read_mesh);
114   CHKERRQ(ierr);
115   simplex = PETSC_FALSE;
116   ierr = PetscOptionsBool("-simplex", "Use simplices, or tensor product cells",
117                           NULL, simplex, &simplex, NULL); CHKERRQ(ierr);
118   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
119 
120   // Setup DM
121   if (read_mesh) {
122     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
123     CHKERRQ(ierr);
124   } else {
125     // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box
126     ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, 1., &dm);
127     CHKERRQ(ierr);
128     // Set the object name
129     ierr = PetscObjectSetName((PetscObject)dm, "Sphere"); CHKERRQ(ierr);
130     // Distribute mesh over processes
131     {
132       DM dmDist = NULL;
133       PetscPartitioner part;
134 
135       ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
136       ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
137       ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
138       if (dmDist) {
139         ierr = DMDestroy(&dm); CHKERRQ(ierr);
140         dm  = dmDist;
141       }
142     }
143     // Refine DMPlex with uniform refinement using runtime option -dm_refine
144     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
145     ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
146     ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);
147     // View DMPlex via runtime option
148     ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
149   }
150 
151   // Create DM
152   ierr = SetupDMByDegree(dm, degree, ncompu, topodim); CHKERRQ(ierr);
153 
154   // Create vectors
155   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
156   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
157   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
158   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
159   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
160   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
161 
162   // Operator
163   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
164   ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize,
165                         userO, &matO); CHKERRQ(ierr);
166   ierr = MatShellSetOperation(matO, MATOP_MULT,
167                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
168 
169   // Set up libCEED
170   CeedInit(ceedresource, &ceed);
171 
172   // Print summary
173   if (!test_mode) {
174     PetscInt P = degree + 1, Q = P + qextra;
175     const char *usedresource;
176     CeedGetResource(ceed, &usedresource);
177     ierr = PetscPrintf(comm,
178                        "\n-- CEED Benchmark Problem %d on the Sphere -- libCEED + PETSc --\n"
179                        "  libCEED:\n"
180                        "    libCEED Backend                    : %s\n"
181                        "  Mesh:\n"
182                        "    Number of 1D Basis Nodes (p)       : %d\n"
183                        "    Number of 1D Quadrature Points (q) : %d\n"
184                        "    Global nodes                       : %D\n",
185                        bpChoice+1, ceedresource, P, Q,  gsize/ncompu);
186     CHKERRQ(ierr);
187   }
188 
189   // Create RHS vector
190   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
191   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
192   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
193   CeedVectorCreate(ceed, xlsize, &rhsceed);
194   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
195 
196   // Setup libCEED's objects
197   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
198   ierr = SetupLibceedByDegree(dm, ceed, degree, topodim, qextra,
199                               ncompx, ncompu, gsize, xlsize, bpChoice,
200                               ceeddata, true, rhsceed, &target); CHKERRQ(ierr);
201 
202   // Gather RHS
203   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
204   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
205   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
206   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
207   CeedVectorDestroy(&rhsceed);
208 
209   // Create the error Q-function
210   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
211                               bpOptions[bpChoice].errorfname, &qf_error);
212   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
213   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
214   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
215 
216   // Create the error operator
217   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
218   CeedOperatorSetField(op_error, "u", ceeddata->Erestrictu,
219                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
220   CeedOperatorSetField(op_error, "true_soln", ceeddata->Erestrictui,
221                        CEED_BASIS_COLLOCATED, target);
222   CeedOperatorSetField(op_error, "error", ceeddata->Erestrictui,
223                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
224 
225   // Set up Mat
226   userO->comm = comm;
227   userO->dm = dm;
228   userO->Xloc = Xloc;
229   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
230   userO->xceed = ceeddata->xceed;
231   userO->yceed = ceeddata->yceed;
232   userO->op = ceeddata->op_apply;
233   userO->ceed = ceed;
234 
235   // Setup solver
236   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
237   {
238     PC pc;
239     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
240     if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
241       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
242       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
243     } else {
244       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
245       MatNullSpace nullspace;
246 
247       ierr = MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace);
248       CHKERRQ(ierr);
249       ierr = MatSetNullSpace(matO, nullspace); CHKERRQ(ierr);
250       ierr = MatNullSpaceDestroy(&nullspace); CHKERRQ(ierr);
251     }
252     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
253     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
254     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
255                             PETSC_DEFAULT); CHKERRQ(ierr);
256   }
257   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
258   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
259 
260   // First run, if benchmarking
261   if (benchmark_mode) {
262     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
263     CHKERRQ(ierr);
264     my_rt_start = MPI_Wtime();
265     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
266     my_rt = MPI_Wtime() - my_rt_start;
267     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
268     CHKERRQ(ierr);
269     // Set maxits based on first iteration timing
270     if (my_rt > 0.02) {
271       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
272       CHKERRQ(ierr);
273     } else {
274       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
275       CHKERRQ(ierr);
276     }
277   }
278 
279   // Timed solve
280   ierr = VecZeroEntries(X); CHKERRQ(ierr);
281   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
282 
283   // -- Performance logging
284   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
285   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
286 
287   // -- Solve
288   my_rt_start = MPI_Wtime();
289   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
290   my_rt = MPI_Wtime() - my_rt_start;
291 
292   // -- Performance logging
293   ierr = PetscLogStagePop();
294 
295   // Output results
296   {
297     KSPType ksptype;
298     KSPConvergedReason reason;
299     PetscReal rnorm;
300     PetscInt its;
301     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
302     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
303     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
304     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
305     if (!test_mode || reason < 0 || rnorm > 1e-8) {
306       ierr = PetscPrintf(comm,
307                          "  KSP:\n"
308                          "    KSP Type                           : %s\n"
309                          "    KSP Convergence                    : %s\n"
310                          "    Total KSP Iterations               : %D\n"
311                          "    Final rnorm                        : %e\n",
312                          ksptype, KSPConvergedReasons[reason], its,
313                          (double)rnorm); CHKERRQ(ierr);
314     }
315     if (!test_mode) {
316       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
317     }
318     {
319       PetscReal maxerror;
320       ierr = ComputeErrorMax(userO, op_error, X, target, &maxerror);
321       CHKERRQ(ierr);
322       PetscReal tol = 5e-4;
323       if (!test_mode || maxerror > tol) {
324         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
325         CHKERRQ(ierr);
326         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
327         CHKERRQ(ierr);
328         ierr = PetscPrintf(comm,
329                            "    Pointwise Error (max)              : %e\n"
330                            "    CG Solve Time                      : %g (%g) sec\n",
331                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
332       }
333     }
334     if (benchmark_mode && (!test_mode)) {
335       ierr = PetscPrintf(comm,
336                          "    DoFs/Sec in CG                     : %g (%g) million\n",
337                          1e-6*gsize*its/rt_max,
338                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
339     }
340   }
341 
342   // Output solution
343   if (write_solution) {
344     PetscViewer vtkviewersoln;
345 
346     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
347     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
348     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr);
349     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
350     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
351   }
352 
353   // Cleanup
354   ierr = VecDestroy(&X); CHKERRQ(ierr);
355   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
356   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
357   ierr = MatDestroy(&matO); CHKERRQ(ierr);
358   ierr = PetscFree(userO); CHKERRQ(ierr);
359   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
360   ierr = DMDestroy(&dm); CHKERRQ(ierr);
361 
362   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
363   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
364   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
365   CeedVectorDestroy(&target);
366   CeedQFunctionDestroy(&qf_error);
367   CeedOperatorDestroy(&op_error);
368   CeedDestroy(&ceed);
369   return PetscFinalize();
370 }
371