xref: /libCEED/examples/petsc/bps.c (revision f6d735e939aa19e4cf1cd71bcb0ef5410bccc942)
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 //
22 // The code uses higher level communication protocols in DMPlex.
23 //
24 // Build with:
25 //
26 //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27 //
28 // Sample runs:
29 //
30 //     ./bps -problem bp1 -degree 3
31 //     ./bps -problem bp2 -ceed /cpu/self -degree 3
32 //     ./bps -problem bp3 -ceed /gpu/occa -degree 3
33 //     ./bps -problem bp4 -ceed /cpu/occa -degree 3
34 //     ./bps -problem bp5 -ceed /omp/occa -degree 3
35 //     ./bps -problem bp6 -ceed /ocl/occa -degree 3
36 //
37 //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3
38 
39 /// @file
40 /// CEED BPs example using PETSc with DMPlex
41 /// See bpsraw.c for a "raw" implementation using a structured grid.
42 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
43 
44 #include <stdbool.h>
45 #include <string.h>
46 #include <petscksp.h>
47 #include <petscdmplex.h>
48 #include <ceed.h>
49 #include "setup.h"
50 
51 int main(int argc, char **argv) {
52   PetscInt ierr;
53   MPI_Comm comm;
54   char filename[PETSC_MAX_PATH_LEN],
55        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
56   double my_rt_start, my_rt, rt_min, rt_max;
57   PetscInt degree = 3, qextra, lsize, gsize, dim = 3, melem[3] = {3, 3, 3},
58            ncompu = 1, xlsize;
59   PetscScalar *r;
60   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
61   PetscLogStage solvestage;
62   Vec X, Xloc, rhs, rhsloc;
63   Mat matO;
64   KSP ksp;
65   DM  dm;
66   UserO userO;
67   Ceed ceed;
68   CeedData ceeddata;
69   CeedQFunction qferror;
70   CeedOperator operror;
71   CeedVector rhsceed, target;
72   CeedMemType memtyperequested;
73   bpType bpchoice;
74 
75   // Check PETSc CUDA support
76   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
77   // *INDENT-OFF*
78   #ifdef PETSC_HAVE_CUDA
79   petschavecuda = PETSC_TRUE;
80   #else
81   petschavecuda = PETSC_FALSE;
82   #endif
83   // *INDENT-ON*
84 
85   ierr = PetscInitialize(&argc, &argv, NULL, help);
86   if (ierr) return ierr;
87   comm = PETSC_COMM_WORLD;
88 
89   // Read command line options
90   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
91   bpchoice = CEED_BP1;
92   ierr = PetscOptionsEnum("-problem",
93                           "CEED benchmark problem to solve", NULL,
94                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
95                           NULL); CHKERRQ(ierr);
96   ncompu = bpOptions[bpchoice].ncompu;
97   test_mode = PETSC_FALSE;
98   ierr = PetscOptionsBool("-test",
99                           "Testing mode (do not print unless error is large)",
100                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
101   benchmark_mode = PETSC_FALSE;
102   ierr = PetscOptionsBool("-benchmark",
103                           "Benchmarking mode (prints benchmark statistics)",
104                           NULL, benchmark_mode, &benchmark_mode, NULL);
105   CHKERRQ(ierr);
106   write_solution = PETSC_FALSE;
107   ierr = PetscOptionsBool("-write_solution",
108                           "Write solution for visualization",
109                           NULL, write_solution, &write_solution, NULL);
110   CHKERRQ(ierr);
111   degree = test_mode ? 3 : 2;
112   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
113                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
114   qextra = bpOptions[bpchoice].qextra;
115   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
116                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
117   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
118                             NULL, ceedresource, ceedresource,
119                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
120   read_mesh = PETSC_FALSE;
121   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
122                             filename, filename, sizeof(filename), &read_mesh);
123   CHKERRQ(ierr);
124   if (!read_mesh) {
125     PetscInt tmp = dim;
126     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
127                                 melem, &tmp, NULL); CHKERRQ(ierr);
128   }
129 
130   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
131   ierr = PetscOptionsEnum("-memtype",
132                           "CEED MemType requested", NULL,
133                           memTypes, (PetscEnum)memtyperequested,
134                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
135   CHKERRQ(ierr);
136 
137   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
138 
139   // Setup DM
140   if (read_mesh) {
141     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
142     CHKERRQ(ierr);
143   } else {
144     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
145                                NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
146   }
147 
148   {
149     DM dmDist = NULL;
150     PetscPartitioner part;
151 
152     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
153     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
154     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
155     if (dmDist) {
156       ierr = DMDestroy(&dm); CHKERRQ(ierr);
157       dm  = dmDist;
158     }
159   }
160 
161   // Set up libCEED
162   CeedInit(ceedresource, &ceed);
163   CeedMemType memtypebackend;
164   CeedGetPreferredMemType(ceed, &memtypebackend);
165 
166   // Check memtype compatibility
167   if (!setmemtyperequest)
168     memtyperequested = memtypebackend;
169   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
170     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
171              "PETSc was not built with CUDA. "
172              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
173 
174   // Create DM
175   ierr = SetupDMByDegree(dm, degree, ncompu, bpchoice);
176   CHKERRQ(ierr);
177 
178   // Create vectors
179   if (memtyperequested == CEED_MEM_DEVICE) {
180     ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr);
181   }
182   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
183   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
184   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
185   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
186   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
187   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
188 
189   // Operator
190   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
191   ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize,
192                         userO, &matO); CHKERRQ(ierr);
193   ierr = MatShellSetOperation(matO, MATOP_MULT,
194                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
195   ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL,
196                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
197   if (memtyperequested == CEED_MEM_DEVICE) {
198     ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr);
199   }
200 
201   // Print summary
202   if (!test_mode) {
203     PetscInt P = degree + 1, Q = P + qextra;
204 
205     const char *usedresource;
206     CeedGetResource(ceed, &usedresource);
207 
208     VecType vectype;
209     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
210 
211     ierr = PetscPrintf(comm,
212                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
213                        "  PETSc:\n"
214                        "    PETSc Vec Type                     : %s\n"
215                        "  libCEED:\n"
216                        "    libCEED Backend                    : %s\n"
217                        "    libCEED Backend MemType            : %s\n"
218                        "    libCEED User Requested MemType     : %s\n"
219                        "  Mesh:\n"
220                        "    Number of 1D Basis Nodes (p)       : %d\n"
221                        "    Number of 1D Quadrature Points (q) : %d\n"
222                        "    Global nodes                       : %D\n"
223                        "    Owned nodes                        : %D\n"
224                        "    DoF per node                       : %D\n",
225                        bpchoice+1, vectype, usedresource,
226                        CeedMemTypes[memtypebackend],
227                        (setmemtyperequest) ?
228                        CeedMemTypes[memtyperequested] : "none",
229                        P, Q, gsize/ncompu, lsize/ncompu, ncompu); CHKERRQ(ierr);
230   }
231 
232   // Create RHS vector
233   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
234   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
235   if (memtyperequested == CEED_MEM_HOST) {
236     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
237   } else {
238     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
239   }
240   CeedVectorCreate(ceed, xlsize, &rhsceed);
241   CeedVectorSetArray(rhsceed, memtyperequested, CEED_USE_POINTER, r);
242 
243   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
244   ierr = SetupLibceedByDegree(dm, ceed, degree, dim, qextra,
245                               ncompu, gsize, xlsize, bpchoice, ceeddata,
246                               true, rhsceed, &target); CHKERRQ(ierr);
247 
248   // Gather RHS
249   CeedVectorSyncArray(rhsceed, memtyperequested);
250   if (memtyperequested == CEED_MEM_HOST) {
251     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
252   } else {
253     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
254   }
255   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
256   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
257   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
258   CeedVectorDestroy(&rhsceed);
259 
260   // Create the error Q-function
261   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
262                               bpOptions[bpchoice].errorfname, &qferror);
263   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
264   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
265   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
266 
267   // Create the error operator
268   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
269                      &operror);
270   CeedOperatorSetField(operror, "u", ceeddata->Erestrictu,
271                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
272   CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui,
273                        CEED_BASIS_COLLOCATED, target);
274   CeedOperatorSetField(operror, "error", ceeddata->Erestrictui,
275                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
276 
277   // Set up Mat
278   userO->comm = comm;
279   userO->dm = dm;
280   userO->Xloc = Xloc;
281   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
282   userO->xceed = ceeddata->xceed;
283   userO->yceed = ceeddata->yceed;
284   userO->op = ceeddata->opapply;
285   userO->ceed = ceed;
286   userO->memtype = memtyperequested;
287   if (memtyperequested == CEED_MEM_HOST) {
288     userO->VecGetArray = VecGetArray;
289     userO->VecGetArrayRead = VecGetArrayRead;
290     userO->VecRestoreArray = VecRestoreArray;
291     userO->VecRestoreArrayRead = VecRestoreArrayRead;
292   } else {
293     userO->VecGetArray = VecCUDAGetArray;
294     userO->VecGetArrayRead = VecCUDAGetArrayRead;
295     userO->VecRestoreArray = VecCUDARestoreArray;
296     userO->VecRestoreArrayRead = VecCUDARestoreArrayRead;
297   }
298 
299   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
300   {
301     PC pc;
302     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
303     if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
304       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
305       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
306     } else {
307       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
308     }
309     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
310     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
311     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
312                             PETSC_DEFAULT); CHKERRQ(ierr);
313   }
314   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
315 
316   // First run, if benchmarking
317   if (benchmark_mode) {
318     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
319     CHKERRQ(ierr);
320     my_rt_start = MPI_Wtime();
321     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
322     my_rt = MPI_Wtime() - my_rt_start;
323     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
324     CHKERRQ(ierr);
325     // Set maxits based on first iteration timing
326     if (my_rt > 0.02) {
327       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
328       CHKERRQ(ierr);
329     } else {
330       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
331       CHKERRQ(ierr);
332     }
333   }
334   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
335 
336   // Timed solve
337   ierr = VecZeroEntries(X); CHKERRQ(ierr);
338   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
339 
340   // -- Performance logging
341   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
342   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
343 
344   // -- Solve
345   my_rt_start = MPI_Wtime();
346   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
347   my_rt = MPI_Wtime() - my_rt_start;
348 
349   // -- Performance logging
350   ierr = PetscLogStagePop();
351 
352   // Output results
353   {
354     KSPType ksptype;
355     KSPConvergedReason reason;
356     PetscReal rnorm;
357     PetscInt its;
358     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
359     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
360     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
361     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
362     if (!test_mode || reason < 0 || rnorm > 1e-8) {
363       ierr = PetscPrintf(comm,
364                          "  KSP:\n"
365                          "    KSP Type                           : %s\n"
366                          "    KSP Convergence                    : %s\n"
367                          "    Total KSP Iterations               : %D\n"
368                          "    Final rnorm                        : %e\n",
369                          ksptype, KSPConvergedReasons[reason], its,
370                          (double)rnorm); CHKERRQ(ierr);
371     }
372     if (!test_mode) {
373       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
374     }
375     {
376       PetscReal maxerror;
377       ierr = ComputeErrorMax(userO, operror, X, target, &maxerror);
378       CHKERRQ(ierr);
379       PetscReal tol = 5e-2;
380       if (!test_mode || maxerror > tol) {
381         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
382         CHKERRQ(ierr);
383         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
384         CHKERRQ(ierr);
385         ierr = PetscPrintf(comm,
386                            "    Pointwise Error (max)              : %e\n"
387                            "    CG Solve Time                      : %g (%g) sec\n",
388                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
389       }
390     }
391     if (benchmark_mode && (!test_mode)) {
392       ierr = PetscPrintf(comm,
393                          "    DoFs/Sec in CG                     : %g (%g) million\n",
394                          1e-6*gsize*its/rt_max,
395                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
396     }
397   }
398 
399   if (write_solution) {
400     PetscViewer vtkviewersoln;
401 
402     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
403     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
404     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
405     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
406     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
407   }
408 
409   // Cleanup
410   ierr = VecDestroy(&X); CHKERRQ(ierr);
411   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
412   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
413   ierr = MatDestroy(&matO); CHKERRQ(ierr);
414   ierr = PetscFree(userO); CHKERRQ(ierr);
415   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
416   ierr = DMDestroy(&dm); CHKERRQ(ierr);
417 
418   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
419   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
420   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
421   CeedVectorDestroy(&target);
422   CeedQFunctionDestroy(&qferror);
423   CeedOperatorDestroy(&operror);
424   CeedDestroy(&ceed);
425   return PetscFinalize();
426 }
427