xref: /libCEED/examples/solids/elasticity.c (revision ccaff0309dc399f656ea11018b919b30feb8b0fa)
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: Elasticity
18 //
19 // This example demonstrates a simple usage of libCEED with PETSc to solve
20 //   elasticity problems.
21 //
22 // The code uses higher level communication protocols in DMPlex.
23 //
24 // Build with:
25 //
26 //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27 //
28 // Sample runs:
29 //
30 //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem linElas -forcing mms
31 //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperSS -forcing none -ceed /cpu/self
32 //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperFS -forcing none -ceed /gpu/occa
33 //
34 // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
35 //
36 //TESTARGS -ceed {ceed_resource} -test -degree 2 -nu 0.3 -E 1
37 
38 /// @file
39 /// CEED elasticity example using PETSc with DMPlex
40 
41 const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
42 
43 #include "elasticity.h"
44 
45 int main(int argc, char **argv) {
46   PetscInt       ierr;
47   MPI_Comm       comm;
48   // Context structs
49   AppCtx         appCtx;                 // Contains problem options
50   Physics        phys;                   // Contains physical constants
51   Units          units;                  // Contains units scaling
52   // PETSc objects
53   PetscLogStage  stageDMSetup, stageLibceedSetup,
54                  stageSnesSetup, stageSnesSolve;
55   DM             dmOrig;                 // Distributed DM to clone
56   DM             *levelDMs;
57   Vec            U, *Ug, *Uloc;          // U: solution, R: residual, F: forcing
58   Vec            R, Rloc, F, Floc;       // g: global, loc: local
59   SNES           snes, snesCoarse;
60   Mat            *jacobMat, jacobMatCoarse, *prolongRestrMat;
61   // PETSc data
62   UserMult       resCtx, jacobCoarseCtx, *jacobCtx;
63   FormJacobCtx   formJacobCtx;
64   UserMultProlongRestr *prolongRestrCtx;
65   PCMGCycleType  pcmgCycleType = PC_MG_CYCLE_V;
66   // libCEED objects
67   Ceed           ceed, ceedFine = NULL;
68   CeedData       *ceedData;
69   CeedQFunction  qfRestrict, qfProlong;
70   // Parameters
71   PetscInt       ncompu = 3;             // 3 DoFs in 3D
72   PetscInt       numLevels = 1, fineLevel = 0;
73   PetscInt       *Ugsz, *Ulsz, *Ulocsz;  // sz: size
74   PetscInt       snesIts = 0;
75   // Timing
76   double         startTime, elapsedTime, minTime, maxTime;
77 
78   ierr = PetscInitialize(&argc, &argv, NULL, help);
79   if (ierr)
80     return ierr;
81 
82   // ---------------------------------------------------------------------------
83   // Process command line options
84   // ---------------------------------------------------------------------------
85   comm = PETSC_COMM_WORLD;
86 
87   // -- Set mesh file, polynomial degree, problem type
88   ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr);
89   ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr);
90   numLevels = appCtx->numLevels;
91   fineLevel = numLevels - 1;
92 
93   // -- Set Poison's ratio, Young's Modulus
94   ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
95   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
96   ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr);
97 
98   // ---------------------------------------------------------------------------
99   // Setup DM
100   // ---------------------------------------------------------------------------
101   // Performance logging
102   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup);
103   CHKERRQ(ierr);
104   ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr);
105 
106   // -- Create distributed DM from mesh file
107   ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr);
108 
109   // -- Setup DM by polynomial degree
110   ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr);
111   for (int level = 0; level < numLevels; level++) {
112     ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr);
113     ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level],
114                            ncompu); CHKERRQ(ierr);
115   }
116 
117   // ---------------------------------------------------------------------------
118   // Setup solution and work vectors
119   // ---------------------------------------------------------------------------
120   // Allocate arrays
121   ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr);
122   ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr);
123   ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr);
124   ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr);
125   ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr);
126 
127   // -- Setup solution vectors for each level
128   for (int level = 0; level < numLevels; level++) {
129     // -- Create global unknown vector U
130     ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr);
131     ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr);
132     // Note: Local size for matShell
133     ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr);
134 
135     // -- Create local unknown vector Uloc
136     ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr);
137     // Note: local size for libCEED
138     ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr);
139   }
140 
141   // -- Create residual and forcing vectors
142   ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr);
143   ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr);
144   ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr);
145   ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr);
146   ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr);
147 
148   // Performance logging
149   ierr = PetscLogStagePop();
150 
151   // ---------------------------------------------------------------------------
152   // Set up libCEED
153   // ---------------------------------------------------------------------------
154   // Performance logging
155   ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup);
156   CHKERRQ(ierr);
157   ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr);
158 
159   // Initalize backend
160   CeedInit(appCtx->ceedResource, &ceed);
161   if (appCtx->degree > 4 && appCtx->ceedResourceFine[0])
162     CeedInit(appCtx->ceedResourceFine, &ceedFine);
163 
164   // -- Create libCEED local forcing vector
165   CeedVector forceCeed;
166   CeedScalar *f;
167   if (appCtx->forcingChoice != FORCE_NONE) {
168     ierr = VecGetArray(Floc, &f); CHKERRQ(ierr);
169     CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed);
170     CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f);
171   }
172 
173   // -- Restriction and prolongation QFunction
174   if (appCtx->multigridChoice != MULTIGRID_NONE) {
175     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
176                                 &qfRestrict);
177     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
178                                 &qfProlong);
179   }
180 
181   // -- Setup libCEED objects
182   ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr);
183   // ---- Setup residual evaluator and geometric information
184   ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr);
185   {
186     bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine);
187     Ceed levelCeed = highOrder ? ceedFine : ceed;
188     ierr = SetupLibceedFineLevel(levelDMs[fineLevel], levelCeed, appCtx,
189                                  phys, ceedData, fineLevel, ncompu,
190                                  Ugsz[fineLevel], Ulocsz[fineLevel], forceCeed,
191                                  qfRestrict, qfProlong);
192     CHKERRQ(ierr);
193   }
194   // ---- Setup Jacobian evaluator and prolongation/restriction
195   for (int level = 0; level < numLevels; level++) {
196     if (level != fineLevel) {
197       ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr);
198     }
199 
200     // Note: use high order ceed, if specified and degree > 3
201     bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine);
202     Ceed levelCeed = highOrder ? ceedFine : ceed;
203     ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys,
204                              ceedData,  level, ncompu, Ugsz[level],
205                              Ulocsz[level], forceCeed, qfRestrict,
206                              qfProlong); CHKERRQ(ierr);
207   }
208 
209   // Performance logging
210   ierr = PetscLogStagePop();
211 
212   // ---------------------------------------------------------------------------
213   // Setup global forcing vector
214   // ---------------------------------------------------------------------------
215   ierr = VecZeroEntries(F); CHKERRQ(ierr);
216 
217   if (appCtx->forcingChoice != FORCE_NONE) {
218     ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr);
219     ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F);
220     CHKERRQ(ierr);
221     CeedVectorDestroy(&forceCeed);
222   }
223 
224   // ---------------------------------------------------------------------------
225   // Print problem summary
226   // ---------------------------------------------------------------------------
227   if (!appCtx->testMode) {
228     const char *usedresource;
229     CeedGetResource(ceed, &usedresource);
230 
231     ierr = PetscPrintf(comm,
232                        "\n-- Elastisticy Example - libCEED + PETSc --\n"
233                        "  libCEED:\n"
234                        "    libCEED Backend                    : %s\n",
235                        usedresource); CHKERRQ(ierr);
236 
237     if (ceedFine) {
238       ierr = PetscPrintf(comm,
239                          "    libCEED Backend - high order       : %s\n",
240                          appCtx->ceedResourceFine); CHKERRQ(ierr);
241     }
242 
243     ierr = PetscPrintf(comm,
244                        "  Problem:\n"
245                        "    Problem Name                       : %s\n"
246                        "    Forcing Function                   : %s\n"
247                        "  Mesh:\n"
248                        "    File                               : %s\n"
249                        "    Number of 1D Basis Nodes (p)       : %d\n"
250                        "    Number of 1D Quadrature Points (q) : %d\n"
251                        "    Global nodes                       : %D\n"
252                        "    Owned nodes                        : %D\n"
253                        "    DoF per node                       : %D\n"
254                        "  Multigrid:\n"
255                        "    Type                               : %s\n"
256                        "    Number of Levels                   : %d\n",
257                        problemTypesForDisp[appCtx->problemChoice],
258                        forcingTypesForDisp[appCtx->forcingChoice],
259                        appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh",
260                        appCtx->degree + 1, appCtx->degree + 1,
261                        Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu,
262                        multigridTypesForDisp[appCtx->multigridChoice],
263                        numLevels); CHKERRQ(ierr);
264 
265     if (appCtx->multigridChoice != MULTIGRID_NONE) {
266       for (int i = 0; i < 2; i++) {
267         CeedInt level = i ? fineLevel : 0;
268         ierr = PetscPrintf(comm,
269                            "    Level %D (%s):\n"
270                            "      Number of 1D Basis Nodes (p)     : %d\n"
271                            "      Global Nodes                     : %D\n"
272                            "      Owned Nodes                      : %D\n",
273                            level, i ? "fine" : "coarse",
274                            appCtx->levelDegrees[level] + 1,
275                            Ugsz[level]/ncompu, Ulsz[level]/ncompu);
276         CHKERRQ(ierr);
277       }
278     }
279   }
280 
281   // ---------------------------------------------------------------------------
282   // Setup SNES
283   // ---------------------------------------------------------------------------
284   // Performance logging
285   ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup);
286   CHKERRQ(ierr);
287   ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr);
288 
289   // Create SNES
290   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
291   ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr);
292 
293   // -- Jacobian evaluators
294   ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr);
295   ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr);
296   for (int level = 0; level < numLevels; level++) {
297     // -- Jacobian context for level
298     ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr);
299     ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level],
300                             Uloc[level], ceedData[level], ceed,
301                             jacobCtx[level]); CHKERRQ(ierr);
302 
303     // -- Form Action of Jacobian on delta_u
304     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level],
305                           Ugsz[level], jacobCtx[level], &jacobMat[level]);
306     CHKERRQ(ierr);
307     ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT,
308                                 (void (*)(void))ApplyJacobian_Ceed);
309     CHKERRQ(ierr);
310     ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL,
311                                 (void(*)(void))GetDiag_Ceed);
312 
313   }
314   // Note: FormJacobian updates the state count of the Jacobian diagonals
315   //   and assembles the Jpre matrix, if needed
316   ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr);
317   formJacobCtx->jacobCtx = jacobCtx;
318   formJacobCtx->numLevels = numLevels;
319   formJacobCtx->jacobMat = jacobMat;
320   ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel],
321                          FormJacobian, formJacobCtx); CHKERRQ(ierr);
322 
323   // -- Residual evaluation function
324   ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr);
325   ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel],
326                      sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr);
327   resCtx->op = ceedData[fineLevel]->opApply;
328   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr);
329 
330   // -- Prolongation/Restriction evaluation
331   ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr);
332   ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr);
333   for (int level = 1; level < numLevels; level++) {
334     // ---- Prolongation/restriction context for level
335     ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr);
336     ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level],
337                                    Ug[level], Uloc[level-1], Uloc[level],
338                                    ceedData[level-1], ceedData[level], ceed,
339                                    prolongRestrCtx[level]); CHKERRQ(ierr);
340 
341     // ---- Form Action of Jacobian on delta_u
342     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level],
343                           Ugsz[level-1], prolongRestrCtx[level],
344                           &prolongRestrMat[level]); CHKERRQ(ierr);
345     // Note: In PCMG, restriction is the transpose of prolongation
346     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT,
347                                 (void (*)(void))Prolong_Ceed);
348     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE,
349                                 (void (*)(void))Restrict_Ceed);
350     CHKERRQ(ierr);
351   }
352 
353   // ---------------------------------------------------------------------------
354   // Setup dummy SNES for AMG coarse solve
355   // ---------------------------------------------------------------------------
356   ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr);
357   ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr);
358   ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr);
359 
360   // -- Jacobian matrix
361   ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr);
362   ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr);
363   ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL,
364                          NULL); CHKERRQ(ierr);
365 
366   // -- Residual evaluation function
367   ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr);
368   ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0]));
369   CHKERRQ(ierr);
370   ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed,
371                          jacobCoarseCtx); CHKERRQ(ierr);
372 
373   // -- Update formJacobCtx
374   formJacobCtx->Ucoarse = Ug[0];
375   formJacobCtx->snesCoarse = snesCoarse;
376   formJacobCtx->jacobMatCoarse = jacobMatCoarse;
377 
378   // ---------------------------------------------------------------------------
379   // Setup KSP
380   // ---------------------------------------------------------------------------
381   {
382     PC pc;
383     KSP ksp;
384 
385     // -- KSP
386     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
387     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
388     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
389     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
390                             PETSC_DEFAULT); CHKERRQ(ierr);
391     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
392 
393     // -- Preconditioning
394     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
395     ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr);
396     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
397 
398     if (appCtx->multigridChoice == MULTIGRID_NONE) {
399       // ---- No Multigrid
400       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
401       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
402     } else {
403       // ---- PCMG
404       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
405 
406       // ------ PCMG levels
407       ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr);
408       for (int level = 0; level < numLevels; level++) {
409         // -------- Smoother
410         KSP kspSmoother, kspEst;
411         PC pcSmoother;
412 
413         // ---------- Smoother KSP
414         ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr);
415         ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr);
416         ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr);
417 
418         // ---------- Chebyshev options
419         ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr);
420         ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1);
421         CHKERRQ(ierr);
422         ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr);
423         ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr);
424         ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE);
425         CHKERRQ(ierr);
426         ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]);
427         CHKERRQ(ierr);
428 
429         // ---------- Smoother preconditioner
430         ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr);
431         ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr);
432         ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
433 
434         // -------- Work vector
435         if (level != fineLevel) {
436           ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr);
437         }
438 
439         // -------- Level prolongation/restriction operator
440         if (level > 0) {
441           ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]);
442           CHKERRQ(ierr);
443           ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]);
444           CHKERRQ(ierr);
445         }
446       }
447 
448       // ------ PCMG coarse solve
449       KSP kspCoarse;
450       PC pcCoarse;
451 
452       // -------- Coarse KSP
453       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
454       ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr);
455       ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse);
456       CHKERRQ(ierr);
457       ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr);
458 
459       // -------- Coarse preconditioner
460       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
461       ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr);
462       ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr);
463 
464       ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr);
465       ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr);
466 
467       // ------ PCMG options
468       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
469       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
470       ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr);
471     }
472     ierr = KSPSetFromOptions(ksp);
473     ierr = PCSetFromOptions(pc);
474   }
475   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
476 
477   // Performance logging
478   ierr = PetscLogStagePop();
479 
480   // ---------------------------------------------------------------------------
481   // Set initial guess
482   // ---------------------------------------------------------------------------
483   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
484 
485   // View solution
486   if (appCtx->viewSoln) {
487     ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr);
488   }
489 
490   // ---------------------------------------------------------------------------
491   // Solve SNES
492   // ---------------------------------------------------------------------------
493   // Performance logging
494   ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve);
495   CHKERRQ(ierr);
496   ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr);
497 
498   // Timing
499   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
500   startTime = MPI_Wtime();
501 
502   // Solve for each load increment
503   for (PetscInt increment = 1; increment <= appCtx->numIncrements; increment++) {
504     // -- Scale the problem
505     PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements,
506                 scalingFactor = loadIncrement /
507                                 (increment == 1 ? 1 : resCtx->loadIncrement);
508     resCtx->loadIncrement = loadIncrement;
509     if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) {
510       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
511     }
512 
513     // -- Solve
514     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
515 
516     // -- View solution
517     if (appCtx->viewSoln) {
518       ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr);
519     }
520 
521     // -- Update SNES iteration count
522     PetscInt its;
523     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
524     snesIts += its;
525   }
526 
527   // Timing
528   elapsedTime = MPI_Wtime() - startTime;
529 
530   // Performance logging
531   ierr = PetscLogStagePop();
532 
533   // ---------------------------------------------------------------------------
534   // Output summary
535   // ---------------------------------------------------------------------------
536   if (!appCtx->testMode) {
537     // -- SNES
538     SNESType snesType;
539     SNESConvergedReason reason;
540     PetscReal rnorm;
541     ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr);
542     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
543     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
544     ierr = PetscPrintf(comm,
545                        "  SNES:\n"
546                        "    SNES Type                          : %s\n"
547                        "    SNES Convergence                   : %s\n"
548                        "    Number of Load Increments          : %d\n"
549                        "    Total SNES Iterations              : %D\n"
550                        "    Final rnorm                        : %e\n",
551                        snesType, SNESConvergedReasons[reason],
552                        appCtx->numIncrements, snesIts, (double)rnorm);
553     CHKERRQ(ierr);
554 
555     // -- KSP
556     KSP ksp;
557     KSPType kspType;
558     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
559     ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr);
560     ierr = PetscPrintf(comm,
561                        "  Linear Solver:\n"
562                        "    KSP Type                           : %s\n",
563                        kspType); CHKERRQ(ierr);
564 
565     // -- PC
566     if (appCtx->multigridChoice != MULTIGRID_NONE) {
567       PC pc;
568       PCMGType pcmgType;
569       ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
570       ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr);
571       ierr = PetscPrintf(comm,
572                          "  P-Multigrid:\n"
573                          "    PCMG Type                          : %s\n"
574                          "    PCMG Cycle Type                    : %s\n",
575                          PCMGTypes[pcmgType],
576                          PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr);
577 
578       // -- Coarse Solve
579       KSP kspCoarse;
580       PC pcCoarse;
581       PCType pcType;
582 
583       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
584       ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr);
585       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
586       ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr);
587       ierr = PetscPrintf(comm,
588                          "    Coarse Solve:\n"
589                          "      KSP Type                         : %s\n"
590                          "      PC Type                          : %s\n",
591                          kspType, pcType); CHKERRQ(ierr);
592     }
593   }
594 
595   // ---------------------------------------------------------------------------
596   // Compute solve time
597   // ---------------------------------------------------------------------------
598   if (!appCtx->testMode) {
599     ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm);
600     CHKERRQ(ierr);
601     ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm);
602     CHKERRQ(ierr);
603     ierr = PetscPrintf(comm,
604                        "  Performance:\n"
605                        "    SNES Solve Time                    : %g (%g) sec\n",
606                        maxTime, minTime); CHKERRQ(ierr);
607   }
608 
609   // ---------------------------------------------------------------------------
610   // Compute error
611   // ---------------------------------------------------------------------------
612   if (appCtx->forcingChoice == FORCE_MMS) {
613     CeedScalar l2Error = 1., l2Unorm = 1.;
614     const CeedScalar *truearray;
615     Vec errorVec, trueVec;
616 
617     // -- Work vectors
618     ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr);
619     ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr);
620     ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr);
621     ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr);
622 
623     // -- Assemble global true solution vector
624     CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST,
625                            &truearray);
626     ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr);
627     ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec);
628     CHKERRQ(ierr);
629     ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr);
630     CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray);
631 
632     // -- Compute L2 error
633     ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr);
634     ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr);
635     ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr);
636     l2Error /= l2Unorm;
637 
638     // -- Output
639     if (!appCtx->testMode || l2Error > 0.05) {
640       ierr = PetscPrintf(comm,
641                          "    L2 Error                           : %e\n",
642                          l2Error); CHKERRQ(ierr);
643     }
644 
645     // -- Cleanup
646     ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
647     ierr = VecDestroy(&trueVec); CHKERRQ(ierr);
648   }
649 
650   // ---------------------------------------------------------------------------
651   // Free objects
652   // ---------------------------------------------------------------------------
653   // Data in arrays per level
654   for (int level = 0; level < numLevels; level++) {
655     // Vectors
656     ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr);
657     ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr);
658 
659     // Jacobian matrix and data
660     ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr);
661     ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr);
662     ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr);
663 
664     // Prolongation/Restriction matrix and data
665     if (level > 0) {
666       ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr);
667       ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr);
668       ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr);
669     }
670 
671     // DM
672     ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr);
673 
674     // libCEED objects
675     ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr);
676   }
677 
678   // Arrays
679   ierr = PetscFree(Ug); CHKERRQ(ierr);
680   ierr = PetscFree(Uloc); CHKERRQ(ierr);
681   ierr = PetscFree(Ugsz); CHKERRQ(ierr);
682   ierr = PetscFree(Ulsz); CHKERRQ(ierr);
683   ierr = PetscFree(Ulocsz); CHKERRQ(ierr);
684   ierr = PetscFree(jacobCtx); CHKERRQ(ierr);
685   ierr = PetscFree(jacobMat); CHKERRQ(ierr);
686   ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr);
687   ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr);
688   ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr);
689   ierr = PetscFree(ceedData); CHKERRQ(ierr);
690 
691   // libCEED objects
692   CeedQFunctionDestroy(&qfRestrict);
693   CeedQFunctionDestroy(&qfProlong);
694   CeedDestroy(&ceed);
695   CeedDestroy(&ceedFine);
696 
697   // PETSc objects
698   ierr = VecDestroy(&U); CHKERRQ(ierr);
699   ierr = VecDestroy(&R); CHKERRQ(ierr);
700   ierr = VecDestroy(&Rloc); CHKERRQ(ierr);
701   ierr = VecDestroy(&F); CHKERRQ(ierr);
702   ierr = VecDestroy(&Floc); CHKERRQ(ierr);
703   ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr);
704   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
705   ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr);
706   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
707   ierr = PetscFree(levelDMs); CHKERRQ(ierr);
708 
709   // Structs
710   ierr = PetscFree(resCtx); CHKERRQ(ierr);
711   ierr = PetscFree(formJacobCtx); CHKERRQ(ierr);
712   ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr);
713   ierr = PetscFree(appCtx); CHKERRQ(ierr);
714   ierr = PetscFree(phys); CHKERRQ(ierr);
715   ierr = PetscFree(units); CHKERRQ(ierr);
716 
717   return PetscFinalize();
718 }
719