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