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