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