xref: /libCEED/examples/solids/elasticity.c (revision 8ca58ff396e94dede9101064bc6e758773b6f761)
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_clamp 998,999 -bc_clamp_998_translate 0.1,0.2,0.3 -problem hyperSS -forcing none -ceed /cpu/self
32 //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_clamp 998,999 -bc_clamp_998 rotate 1,0,0,0.2 -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 -dm_plex_box_faces 3,3,3
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   Physics        physSmoother = NULL;    // Separate context if nuSmoother set
52   Units          units;                  // Contains units scaling
53   // PETSc objects
54   PetscLogStage  stageDMSetup, stageLibceedSetup,
55                  stageSnesSetup, stageSnesSolve;
56   DM             dmOrig;                 // Distributed DM to clone
57   DM             dmEnergy, dmDiagnostic; // DMs for postprocessing
58   DM             *levelDMs;
59   Vec            U, *Ug, *Uloc;          // U: solution, R: residual, F: forcing
60   Vec            R, Rloc, F, Floc;       // g: global, loc: local
61   SNES           snes, snesCoarse = NULL;
62   Mat            *jacobMat, jacobMatCoarse, *prolongRestrMat;
63   // PETSc data
64   UserMult       resCtx, jacobCoarseCtx = NULL, *jacobCtx;
65   FormJacobCtx   formJacobCtx;
66   UserMultProlongRestr *prolongRestrCtx;
67   PCMGCycleType  pcmgCycleType = PC_MG_CYCLE_V;
68   // libCEED objects
69   Ceed           ceed, ceedFine = NULL;
70   CeedData       *ceedData;
71   CeedQFunction  qfRestrict = NULL, qfProlong = NULL;
72   // Parameters
73   PetscInt       ncompu = 3;             // 3 DoFs in 3D
74   PetscInt       ncompe = 1, ncompd = 2; // 1 energy output, 2 diagnostic
75   PetscInt       numLevels = 1, fineLevel = 0;
76   PetscInt       *Ugsz, *Ulsz, *Ulocsz;  // sz: size
77   PetscInt       snesIts = 0;
78   // Timing
79   double         startTime, elapsedTime, minTime, maxTime;
80 
81   ierr = PetscInitialize(&argc, &argv, NULL, help);
82   if (ierr)
83     return ierr;
84 
85   // ---------------------------------------------------------------------------
86   // Process command line options
87   // ---------------------------------------------------------------------------
88   comm = PETSC_COMM_WORLD;
89 
90   // -- Set mesh file, polynomial degree, problem type
91   ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr);
92   ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr);
93   numLevels = appCtx->numLevels;
94   fineLevel = numLevels - 1;
95 
96   // -- Set Poison's ratio, Young's Modulus
97   ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
98   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
99   ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr);
100   if (fabs(appCtx->nuSmoother) > 1E-14) {
101     ierr = PetscMalloc1(1, &physSmoother); CHKERRQ(ierr);
102     ierr = PetscMemcpy(physSmoother, phys, sizeof(*phys)); CHKERRQ(ierr);
103     physSmoother->nu = appCtx->nuSmoother;
104   }
105 
106   // ---------------------------------------------------------------------------
107   // Setup DM
108   // ---------------------------------------------------------------------------
109   // Performance logging
110   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup);
111   CHKERRQ(ierr);
112   ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr);
113 
114   // -- Create distributed DM from mesh file
115   ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr);
116 
117   // -- Setup DM by polynomial degree
118   ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr);
119   for (PetscInt level = 0; level < numLevels; level++) {
120     ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr);
121     ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level],
122                            PETSC_TRUE, ncompu); CHKERRQ(ierr);
123     // -- Label field components for viewing
124     // Empty name for conserved field (because there is only one field)
125     PetscSection section;
126     ierr = DMGetLocalSection(levelDMs[level], &section); CHKERRQ(ierr);
127     ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr);
128     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
129     CHKERRQ(ierr);
130     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
131     CHKERRQ(ierr);
132     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
133     CHKERRQ(ierr);
134   }
135 
136   // -- Setup postprocessing DMs
137   ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr);
138   ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel],
139                          PETSC_FALSE, ncompe); CHKERRQ(ierr);
140   ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr);
141   ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel],
142                          PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr);
143   {
144     // -- Label field components for viewing
145     // Empty name for conserved field (because there is only one field)
146     PetscSection section;
147     ierr = DMGetLocalSection(dmDiagnostic, &section); CHKERRQ(ierr);
148     ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr);
149     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
150     CHKERRQ(ierr);
151     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
152     CHKERRQ(ierr);
153     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
154     CHKERRQ(ierr);
155     ierr = PetscSectionSetComponentName(section, 0, 3, "CondensedPressure");
156     CHKERRQ(ierr);
157     ierr = PetscSectionSetComponentName(section, 0, 4, "StrainEnergyDensity");
158     CHKERRQ(ierr);
159   }
160 
161   // ---------------------------------------------------------------------------
162   // Setup solution and work vectors
163   // ---------------------------------------------------------------------------
164   // Allocate arrays
165   ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr);
166   ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr);
167   ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr);
168   ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr);
169   ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr);
170 
171   // -- Setup solution vectors for each level
172   for (PetscInt level = 0; level < numLevels; level++) {
173     // -- Create global unknown vector U
174     ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr);
175     ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr);
176     // Note: Local size for matShell
177     ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr);
178 
179     // -- Create local unknown vector Uloc
180     ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr);
181     // Note: local size for libCEED
182     ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr);
183   }
184 
185   // -- Create residual and forcing vectors
186   ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr);
187   ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr);
188   ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr);
189   ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr);
190   ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr);
191 
192   // Performance logging
193   ierr = PetscLogStagePop();
194 
195   // ---------------------------------------------------------------------------
196   // Set up libCEED
197   // ---------------------------------------------------------------------------
198   // Performance logging
199   ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup);
200   CHKERRQ(ierr);
201   ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr);
202 
203   // Initalize backend
204   CeedInit(appCtx->ceedResource, &ceed);
205   if (appCtx->degree > 4 && appCtx->ceedResourceFine[0])
206     CeedInit(appCtx->ceedResourceFine, &ceedFine);
207 
208   // -- Create libCEED local forcing vector
209   CeedVector forceCeed;
210   CeedScalar *f;
211   if (appCtx->forcingChoice != FORCE_NONE) {
212     ierr = VecGetArray(Floc, &f); CHKERRQ(ierr);
213     CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed);
214     CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f);
215   }
216 
217   // -- Restriction and prolongation QFunction
218   if (appCtx->multigridChoice != MULTIGRID_NONE) {
219     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
220                                 &qfRestrict);
221     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
222                                 &qfProlong);
223   }
224 
225   // -- Setup libCEED objects
226   ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr);
227   // ---- Setup residual evaluator and geometric information
228   ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr);
229   {
230     bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine);
231     Ceed levelCeed = highOrder ? ceedFine : ceed;
232     ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic,
233                                  levelCeed, appCtx, phys, ceedData, fineLevel,
234                                  ncompu, Ugsz[fineLevel], Ulocsz[fineLevel],
235                                  forceCeed, qfRestrict, qfProlong);
236     CHKERRQ(ierr);
237   }
238   // ---- Setup Jacobian evaluator and prolongation/restriction
239   for (PetscInt level = 0; level < numLevels; level++) {
240     if (level != fineLevel) {
241       ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr);
242     }
243 
244     // Note: use high order ceed, if specified and degree > 4
245     bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine);
246     Ceed levelCeed = highOrder ? ceedFine : ceed;
247     ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys,
248                              ceedData,  level, ncompu, Ugsz[level],
249                              Ulocsz[level], forceCeed, qfRestrict,
250                              qfProlong); CHKERRQ(ierr);
251   }
252 
253   // Performance logging
254   ierr = PetscLogStagePop();
255 
256   // ---------------------------------------------------------------------------
257   // Setup global forcing vector
258   // ---------------------------------------------------------------------------
259   ierr = VecZeroEntries(F); CHKERRQ(ierr);
260 
261   if (appCtx->forcingChoice != FORCE_NONE) {
262     ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr);
263     ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F);
264     CHKERRQ(ierr);
265     CeedVectorDestroy(&forceCeed);
266   }
267 
268   // ---------------------------------------------------------------------------
269   // Print problem summary
270   // ---------------------------------------------------------------------------
271   if (!appCtx->testMode) {
272     const char *usedresource;
273     CeedGetResource(ceed, &usedresource);
274 
275     ierr = PetscPrintf(comm,
276                        "\n-- Elasticity Example - libCEED + PETSc --\n"
277                        "  libCEED:\n"
278                        "    libCEED Backend                    : %s\n",
279                        usedresource); CHKERRQ(ierr);
280 
281     if (ceedFine) {
282       ierr = PetscPrintf(comm,
283                          "    libCEED Backend - high order       : %s\n",
284                          appCtx->ceedResourceFine); CHKERRQ(ierr);
285     }
286 
287     ierr = PetscPrintf(comm,
288                        "  Problem:\n"
289                        "    Problem Name                       : %s\n"
290                        "    Forcing Function                   : %s\n"
291                        "  Mesh:\n"
292                        "    File                               : %s\n"
293                        "    Number of 1D Basis Nodes (p)       : %d\n"
294                        "    Number of 1D Quadrature Points (q) : %d\n"
295                        "    Global nodes                       : %D\n"
296                        "    Owned nodes                        : %D\n"
297                        "    DoF per node                       : %D\n"
298                        "  Multigrid:\n"
299                        "    Type                               : %s\n"
300                        "    Number of Levels                   : %d\n",
301                        problemTypesForDisp[appCtx->problemChoice],
302                        forcingTypesForDisp[appCtx->forcingChoice],
303                        appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh",
304                        appCtx->degree + 1, appCtx->degree + 1,
305                        Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu,
306                        (appCtx->degree == 1 &&
307                         appCtx->multigridChoice != MULTIGRID_NONE) ?
308                        "Algebraic multigrid" :
309                        multigridTypesForDisp[appCtx->multigridChoice],
310                        (appCtx->degree == 1 ||
311                         appCtx->multigridChoice == MULTIGRID_NONE) ?
312                        0 : numLevels); CHKERRQ(ierr);
313 
314     if (appCtx->multigridChoice != MULTIGRID_NONE) {
315       for (PetscInt i = 0; i < 2; i++) {
316         CeedInt level = i ? fineLevel : 0;
317         ierr = PetscPrintf(comm,
318                            "    Level %D (%s):\n"
319                            "      Number of 1D Basis Nodes (p)     : %d\n"
320                            "      Global Nodes                     : %D\n"
321                            "      Owned Nodes                      : %D\n",
322                            level, i ? "fine" : "coarse",
323                            appCtx->levelDegrees[level] + 1,
324                            Ugsz[level]/ncompu, Ulsz[level]/ncompu);
325         CHKERRQ(ierr);
326       }
327     }
328   }
329 
330   // ---------------------------------------------------------------------------
331   // Setup SNES
332   // ---------------------------------------------------------------------------
333   // Performance logging
334   ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup);
335   CHKERRQ(ierr);
336   ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr);
337 
338   // Create SNES
339   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
340   ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr);
341 
342   // -- Jacobian evaluators
343   ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr);
344   ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr);
345   for (PetscInt level = 0; level < numLevels; level++) {
346     // -- Jacobian context for level
347     ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr);
348     ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level],
349                             Uloc[level], ceedData[level], ceed, phys,
350                             physSmoother, jacobCtx[level]); CHKERRQ(ierr);
351 
352     // -- Form Action of Jacobian on delta_u
353     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level],
354                           Ugsz[level], jacobCtx[level], &jacobMat[level]);
355     CHKERRQ(ierr);
356     ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT,
357                                 (void (*)(void))ApplyJacobian_Ceed);
358     CHKERRQ(ierr);
359     ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL,
360                                 (void(*)(void))GetDiag_Ceed);
361 
362   }
363   // Note: FormJacobian updates Jacobian matrices on each level
364   //   and assembles the Jpre matrix, if needed
365   ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr);
366   formJacobCtx->jacobCtx = jacobCtx;
367   formJacobCtx->numLevels = numLevels;
368   formJacobCtx->jacobMat = jacobMat;
369 
370   // -- Residual evaluation function
371   ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr);
372   ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel],
373                      sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr);
374   resCtx->op = ceedData[fineLevel]->opApply;
375   resCtx->qf = ceedData[fineLevel]->qfApply;
376   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr);
377 
378   // -- Prolongation/Restriction evaluation
379   ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr);
380   ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr);
381   for (PetscInt level = 1; level < numLevels; level++) {
382     // ---- Prolongation/restriction context for level
383     ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr);
384     ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level],
385                                    Ug[level], Uloc[level-1], Uloc[level],
386                                    ceedData[level-1], ceedData[level], ceed,
387                                    prolongRestrCtx[level]); CHKERRQ(ierr);
388 
389     // ---- Form Action of Jacobian on delta_u
390     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level],
391                           Ugsz[level-1], prolongRestrCtx[level],
392                           &prolongRestrMat[level]); CHKERRQ(ierr);
393     // Note: In PCMG, restriction is the transpose of prolongation
394     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT,
395                                 (void (*)(void))Prolong_Ceed);
396     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE,
397                                 (void (*)(void))Restrict_Ceed);
398     CHKERRQ(ierr);
399   }
400 
401   // ---------------------------------------------------------------------------
402   // Setup dummy SNES for AMG coarse solve
403   // ---------------------------------------------------------------------------
404   if (appCtx->multigridChoice != MULTIGRID_NONE) {
405     // -- Jacobian Matrix
406     ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr);
407     ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr);
408 
409     if (appCtx->degree > 1) {
410       ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr);
411       ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr);
412       ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr);
413 
414       // -- Jacobian function
415       ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL,
416                              NULL); CHKERRQ(ierr);
417 
418       // -- Residual evaluation function
419       ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr);
420       ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0]));
421       CHKERRQ(ierr);
422       ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed,
423                              jacobCoarseCtx); CHKERRQ(ierr);
424 
425       // -- Update formJacobCtx
426       formJacobCtx->Ucoarse = Ug[0];
427       formJacobCtx->snesCoarse = snesCoarse;
428       formJacobCtx->jacobMatCoarse = jacobMatCoarse;
429     }
430   }
431 
432   // Set Jacobian function
433   if (appCtx->degree > 1) {
434     ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel],
435                            FormJacobian, formJacobCtx); CHKERRQ(ierr);
436   } else {
437     ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse,
438                            SNESComputeJacobianDefaultColor, NULL);
439     CHKERRQ(ierr);
440   }
441 
442   // ---------------------------------------------------------------------------
443   // Setup KSP
444   // ---------------------------------------------------------------------------
445   {
446     PC pc;
447     KSP ksp;
448 
449     // -- KSP
450     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
451     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
452     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
453     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
454                             PETSC_DEFAULT); CHKERRQ(ierr);
455     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
456 
457     // -- Preconditioning
458     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
459     ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr);
460     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
461 
462     if (appCtx->multigridChoice == MULTIGRID_NONE) {
463       // ---- No Multigrid
464       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
465       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
466     } else if (appCtx->degree == 1) {
467       // ---- AMG for degree 1
468       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
469     } else {
470       // ---- PCMG
471       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
472 
473       // ------ PCMG levels
474       ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr);
475       for (PetscInt level = 0; level < numLevels; level++) {
476         // -------- Smoother
477         KSP kspSmoother, kspEst;
478         PC pcSmoother;
479 
480         // ---------- Smoother KSP
481         ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr);
482         ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr);
483         ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr);
484 
485         // ---------- Chebyshev options
486         ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr);
487         ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1);
488         CHKERRQ(ierr);
489         ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr);
490         ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr);
491         ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE);
492         CHKERRQ(ierr);
493         ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]);
494         CHKERRQ(ierr);
495 
496         // ---------- Smoother preconditioner
497         ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr);
498         ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr);
499         ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
500 
501         // -------- Work vector
502         if (level != fineLevel) {
503           ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr);
504         }
505 
506         // -------- Level prolongation/restriction operator
507         if (level > 0) {
508           ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]);
509           CHKERRQ(ierr);
510           ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]);
511           CHKERRQ(ierr);
512         }
513       }
514 
515       // ------ PCMG coarse solve
516       KSP kspCoarse;
517       PC pcCoarse;
518 
519       // -------- Coarse KSP
520       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
521       ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr);
522       ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse);
523       CHKERRQ(ierr);
524       ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr);
525 
526       // -------- Coarse preconditioner
527       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
528       ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr);
529       ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr);
530 
531       ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr);
532       ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr);
533 
534       // ------ PCMG options
535       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
536       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
537       ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr);
538     }
539     ierr = KSPSetFromOptions(ksp);
540     ierr = PCSetFromOptions(pc);
541   }
542   {
543     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
544     SNESLineSearch linesearch;
545 
546     ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr);
547     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr);
548   }
549 
550   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
551 
552   // Performance logging
553   ierr = PetscLogStagePop();
554 
555   // ---------------------------------------------------------------------------
556   // Set initial guess
557   // ---------------------------------------------------------------------------
558   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
559 
560   // View solution
561   if (appCtx->viewSoln) {
562     ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr);
563   }
564 
565   // ---------------------------------------------------------------------------
566   // Solve SNES
567   // ---------------------------------------------------------------------------
568   PetscBool snesMonitor = PETSC_FALSE;
569   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor);
570   CHKERRQ(ierr);
571 
572   // Performance logging
573   ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve);
574   CHKERRQ(ierr);
575   ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr);
576 
577   // Timing
578   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
579   startTime = MPI_Wtime();
580 
581   // Solve for each load increment
582   PetscInt increment;
583   for (increment = 1; increment <= appCtx->numIncrements; increment++) {
584     // -- Log increment count
585     if (snesMonitor) {
586       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
587       CHKERRQ(ierr);
588     }
589 
590     // -- Scale the problem
591     PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements,
592                 scalingFactor = loadIncrement /
593                                 (increment == 1 ? 1 : resCtx->loadIncrement);
594     resCtx->loadIncrement = loadIncrement;
595     if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) {
596       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
597     }
598 
599     // -- Solve
600     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
601 
602     // -- View solution
603     if (appCtx->viewSoln) {
604       ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr);
605     }
606 
607     // -- Update SNES iteration count
608     PetscInt its;
609     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
610     snesIts += its;
611 
612     // -- Check for divergence
613     SNESConvergedReason reason;
614     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
615     if (reason < 0)
616       break;
617   }
618 
619   // Timing
620   elapsedTime = MPI_Wtime() - startTime;
621 
622   // Performance logging
623   ierr = PetscLogStagePop();
624 
625   // ---------------------------------------------------------------------------
626   // Output summary
627   // ---------------------------------------------------------------------------
628   if (!appCtx->testMode) {
629     // -- SNES
630     SNESType snesType;
631     SNESConvergedReason reason;
632     PetscReal rnorm;
633     ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr);
634     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
635     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
636     ierr = PetscPrintf(comm,
637                        "  SNES:\n"
638                        "    SNES Type                          : %s\n"
639                        "    SNES Convergence                   : %s\n"
640                        "    Number of Load Increments          : %d\n"
641                        "    Completed Load Increments          : %d\n"
642                        "    Total SNES Iterations              : %D\n"
643                        "    Final rnorm                        : %e\n",
644                        snesType, SNESConvergedReasons[reason],
645                        appCtx->numIncrements, increment - 1,
646                        snesIts, (double)rnorm); CHKERRQ(ierr);
647 
648     // -- KSP
649     KSP ksp;
650     KSPType kspType;
651     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
652     ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr);
653     ierr = PetscPrintf(comm,
654                        "  Linear Solver:\n"
655                        "    KSP Type                           : %s\n",
656                        kspType); CHKERRQ(ierr);
657 
658     // -- PC
659     PC pc;
660     PCType pcType;
661     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
662     ierr = PCGetType(pc, &pcType); CHKERRQ(ierr);
663     ierr = PetscPrintf(comm,
664                        "    PC Type                            : %s\n",
665                        pcType); CHKERRQ(ierr);
666 
667     if (!strcmp(pcType, PCMG)) {
668       PCMGType pcmgType;
669       ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr);
670       ierr = PetscPrintf(comm,
671                          "  P-Multigrid:\n"
672                          "    PCMG Type                          : %s\n"
673                          "    PCMG Cycle Type                    : %s\n",
674                          PCMGTypes[pcmgType],
675                          PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr);
676 
677       // -- Coarse Solve
678       KSP kspCoarse;
679       PC pcCoarse;
680       PCType pcType;
681 
682       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
683       ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr);
684       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
685       ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr);
686       ierr = PetscPrintf(comm,
687                          "    Coarse Solve:\n"
688                          "      KSP Type                         : %s\n"
689                          "      PC Type                          : %s\n",
690                          kspType, pcType); CHKERRQ(ierr);
691     }
692   }
693 
694   // ---------------------------------------------------------------------------
695   // Compute solve time
696   // ---------------------------------------------------------------------------
697   if (!appCtx->testMode) {
698     ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm);
699     CHKERRQ(ierr);
700     ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm);
701     CHKERRQ(ierr);
702     ierr = PetscPrintf(comm,
703                        "  Performance:\n"
704                        "    SNES Solve Time                    : %g (%g) sec\n",
705                        maxTime, minTime); CHKERRQ(ierr);
706   }
707 
708   // ---------------------------------------------------------------------------
709   // Compute error
710   // ---------------------------------------------------------------------------
711   if (appCtx->forcingChoice == FORCE_MMS) {
712     CeedScalar l2Error = 1., l2Unorm = 1.;
713     const CeedScalar *truearray;
714     Vec errorVec, trueVec;
715 
716     // -- Work vectors
717     ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr);
718     ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr);
719     ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr);
720     ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr);
721 
722     // -- Assemble global true solution vector
723     CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST,
724                            &truearray);
725     ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr);
726     ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec);
727     CHKERRQ(ierr);
728     ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr);
729     CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray);
730 
731     // -- Compute L2 error
732     ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr);
733     ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr);
734     ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr);
735     l2Error /= l2Unorm;
736 
737     // -- Output
738     if (!appCtx->testMode || l2Error > 0.05) {
739       ierr = PetscPrintf(comm,
740                          "    L2 Error                           : %e\n",
741                          l2Error); CHKERRQ(ierr);
742     }
743 
744     // -- Cleanup
745     ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
746     ierr = VecDestroy(&trueVec); CHKERRQ(ierr);
747   }
748 
749   // ---------------------------------------------------------------------------
750   // Compute energy
751   // ---------------------------------------------------------------------------
752   if (!appCtx->testMode) {
753     // -- Compute L2 error
754     CeedScalar energy;
755     ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy,
756                                U, &energy); CHKERRQ(ierr);
757 
758     // -- Output
759     ierr = PetscPrintf(comm,
760                        "    Strain Energy                      : %e\n",
761                        energy); CHKERRQ(ierr);
762   }
763 
764   // ---------------------------------------------------------------------------
765   // Output diagnostic quantities
766   // ---------------------------------------------------------------------------
767   if (appCtx->viewSoln || appCtx->viewFinalSoln) {
768     // -- Setup context
769     UserMult diagnosticCtx;
770     ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr);
771     ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr);
772     diagnosticCtx->dm = dmDiagnostic;
773     diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic;
774 
775     // -- Compute and output
776     ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U,
777                                     ceedData[fineLevel]->ErestrictDiagnostic);
778     CHKERRQ(ierr);
779 
780     // -- Cleanup
781     ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr);
782   }
783 
784   // ---------------------------------------------------------------------------
785   // Free objects
786   // ---------------------------------------------------------------------------
787   // Data in arrays per level
788   for (PetscInt level = 0; level < numLevels; level++) {
789     // Vectors
790     ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr);
791     ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr);
792 
793     // Jacobian matrix and data
794     ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr);
795     ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr);
796     ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr);
797 
798     // Prolongation/Restriction matrix and data
799     if (level > 0) {
800       ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr);
801       ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr);
802       ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr);
803     }
804 
805     // DM
806     ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr);
807 
808     // libCEED objects
809     ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr);
810   }
811 
812   // Arrays
813   ierr = PetscFree(Ug); CHKERRQ(ierr);
814   ierr = PetscFree(Uloc); CHKERRQ(ierr);
815   ierr = PetscFree(Ugsz); CHKERRQ(ierr);
816   ierr = PetscFree(Ulsz); CHKERRQ(ierr);
817   ierr = PetscFree(Ulocsz); CHKERRQ(ierr);
818   ierr = PetscFree(jacobCtx); CHKERRQ(ierr);
819   ierr = PetscFree(jacobMat); CHKERRQ(ierr);
820   ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr);
821   ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr);
822   ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr);
823   ierr = PetscFree(ceedData); CHKERRQ(ierr);
824 
825   // libCEED objects
826   CeedQFunctionDestroy(&qfRestrict);
827   CeedQFunctionDestroy(&qfProlong);
828   CeedDestroy(&ceed);
829   CeedDestroy(&ceedFine);
830 
831   // PETSc objects
832   ierr = VecDestroy(&U); CHKERRQ(ierr);
833   ierr = VecDestroy(&R); CHKERRQ(ierr);
834   ierr = VecDestroy(&Rloc); CHKERRQ(ierr);
835   ierr = VecDestroy(&F); CHKERRQ(ierr);
836   ierr = VecDestroy(&Floc); CHKERRQ(ierr);
837   ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr);
838   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
839   ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr);
840   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
841   ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr);
842   ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr);
843   ierr = PetscFree(levelDMs); CHKERRQ(ierr);
844 
845   // Structs
846   ierr = PetscFree(resCtx); CHKERRQ(ierr);
847   ierr = PetscFree(formJacobCtx); CHKERRQ(ierr);
848   ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr);
849   ierr = PetscFree(appCtx); CHKERRQ(ierr);
850   ierr = PetscFree(phys); CHKERRQ(ierr);
851   ierr = PetscFree(physSmoother); CHKERRQ(ierr);
852   ierr = PetscFree(units); CHKERRQ(ierr);
853 
854   return PetscFinalize();
855 }
856