xref: /libCEED/examples/solids/elasticity.c (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ccaff030SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ccaff030SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7ccaff030SJeremy L Thompson 
8ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Elasticity
9ccaff030SJeremy L Thompson //
10ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve elasticity problems.
11ccaff030SJeremy L Thompson //
12ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex.
13ccaff030SJeremy L Thompson //
14ccaff030SJeremy L Thompson // Build with:
15ccaff030SJeremy L Thompson //
16ccaff030SJeremy L Thompson //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
17ccaff030SJeremy L Thompson //
18ccaff030SJeremy L Thompson // Sample runs:
19ccaff030SJeremy L Thompson //
20eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms
212b730f8bSJeremy L Thompson //     ./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
222b730f8bSJeremy L Thompson //     /cpu/self
232b730f8bSJeremy L Thompson //     ./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
242b730f8bSJeremy L Thompson //     -ceed /gpu/cuda
25ccaff030SJeremy L Thompson //
26ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
27ccaff030SJeremy L Thompson //
28fa5adaf5SJeremy L Thompson //TESTARGS(name="linear elasticity, MMS")        -ceed {ceed_resource} -test -degree 3 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3
29c8565611SJeremy L Thompson //TESTARGS(name="Neo-Hookean hyperelasticity")   -ceed {ceed_resource} -test -problem FS-NH -E 2.8 -nu 0.4 -degree 2 -dm_plex_box_faces 2,2,2 -num_steps 1 -bc_clamp 6 -bc_traction 5 -bc_traction_5 0,0,-.5 -expect_final_strain_energy 2.124627916174e-01
30c8565611SJeremy L Thompson //TESTARGS(name="Mooney-Rivlin hyperelasticity") -ceed {ceed_resource} -test -problem FS-MR -mu_1 .5 -mu_2 .5 -nu 0.4 -degree 2 -dm_plex_box_faces 2,2,2 -num_steps 1 -bc_clamp 6 -bc_traction 5 -bc_traction_5 0,0,-.5 -expect_final_strain_energy 2.339138880207e-01
31ccaff030SJeremy L Thompson 
32ccaff030SJeremy L Thompson /// @file
33ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
34ccaff030SJeremy L Thompson 
35ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
36ccaff030SJeremy L Thompson 
37ccaff030SJeremy L Thompson #include "elasticity.h"
38ccaff030SJeremy L Thompson 
3949aac155SJeremy L Thompson #include <ceed.h>
4049aac155SJeremy L Thompson #include <petscdmplex.h>
4149aac155SJeremy L Thompson #include <petscsnes.h>
4249aac155SJeremy L Thompson 
main(int argc,char ** argv)43ccaff030SJeremy L Thompson int main(int argc, char **argv) {
44ccaff030SJeremy L Thompson   MPI_Comm comm;
45ccaff030SJeremy L Thompson   // Context structs
46d1d35e2fSjeremylt   AppCtx           app_ctx;            // Contains problem options
475754ecacSJeremy L Thompson   ProblemFunctions problem_functions;  // Setup functions for each problem
48ccaff030SJeremy L Thompson   Units            units;              // Contains units scaling
49ccaff030SJeremy L Thompson   // PETSc objects
502b730f8bSJeremy L Thompson   PetscLogStage stage_dm_setup, stage_libceed_setup, stage_snes_setup, stage_snes_solve;
51d1d35e2fSjeremylt   DM            dm_orig;                   // Distributed DM to clone
52d1d35e2fSjeremylt   DM            dm_energy, dm_diagnostic;  // DMs for postprocessing
53d1d35e2fSjeremylt   DM           *level_dms;
54d1d35e2fSjeremylt   Vec           U, *U_g, *U_loc;     // U: solution, R: residual, F: forcing
55d1d35e2fSjeremylt   Vec           R, R_loc, F, F_loc;  // g: global, loc: local
56d1d35e2fSjeremylt   Vec           neumann_bcs = NULL, bcs_loc = NULL;
5717f0843fSJeremy L Thompson   SNES          snes;
58d1d35e2fSjeremylt   Mat          *jacob_mat, jacob_mat_coarse, *prolong_restr_mat;
59ccaff030SJeremy L Thompson   // PETSc data
60d1d35e2fSjeremylt   UserMult              res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx;
61d1d35e2fSjeremylt   FormJacobCtx          form_jacob_ctx;
62d1d35e2fSjeremylt   UserMultProlongRestr *prolong_restr_ctx;
63d1d35e2fSjeremylt   PCMGCycleType         pcmg_cycle_type = PC_MG_CYCLE_V;
64ccaff030SJeremy L Thompson   // libCEED objects
65d99fa3c5SJeremy L Thompson   Ceed                 ceed;
66d1d35e2fSjeremylt   CeedData            *ceed_data;
67d1d35e2fSjeremylt   CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL;
68ccaff030SJeremy L Thompson   // Parameters
69d1d35e2fSjeremylt   PetscInt  num_comp_u = 3;                  // 3 DoFs in 3D
70d1d35e2fSjeremylt   PetscInt  num_comp_e = 1, num_comp_d = 5;  // 1 energy output, 5 diagnostic
71d1d35e2fSjeremylt   PetscInt  num_levels = 1, fine_level = 0;
72bf0f51feSjeremylt   PetscInt *U_g_size, *U_l_size, *U_loc_size;
73bf0f51feSjeremylt   PetscInt  snes_its = 0, ksp_its = 0;
74d1d35e2fSjeremylt   double    start_time, elapsed_time, min_time, max_time;
75ccaff030SJeremy L Thompson 
762b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
77ccaff030SJeremy L Thompson 
78ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
79ccaff030SJeremy L Thompson   // Process command line options
80ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
81ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
82ccaff030SJeremy L Thompson 
83ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
842b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &app_ctx));
852b730f8bSJeremy L Thompson   PetscCall(ProcessCommandLineOptions(comm, app_ctx));
862b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &problem_functions));
872b730f8bSJeremy L Thompson   PetscCall(RegisterProblems(problem_functions));
88d1d35e2fSjeremylt   num_levels = app_ctx->num_levels;
89d1d35e2fSjeremylt   fine_level = num_levels - 1;
90ccaff030SJeremy L Thompson 
91ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
92d1d35e2fSjeremylt   // Initialize libCEED
9362e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
94d1d35e2fSjeremylt   // Initialize backend
95d1d35e2fSjeremylt   CeedInit(app_ctx->ceed_resource, &ceed);
9662e9c006SJeremy L Thompson 
9762e9c006SJeremy L Thompson   // Check preferred MemType
98d1d35e2fSjeremylt   CeedMemType mem_type_backend;
99d1d35e2fSjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1005754ecacSJeremy L Thompson   // Setup physics context and wrap in libCEED object
1015754ecacSJeremy L Thompson   {
1025754ecacSJeremy L Thompson     PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *);
1032b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name, &SetupPhysics));
1042b730f8bSJeremy L Thompson     if (!SetupPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found", app_ctx->name);
1052b730f8bSJeremy L Thompson     PetscCall((*SetupPhysics)(comm, ceed, &units, &ctx_phys));
1062b730f8bSJeremy L Thompson     PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext, CeedQFunctionContext *);
1072b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupSmootherPhysics, app_ctx->name, &SetupSmootherPhysics));
1082b730f8bSJeremy L Thompson     if (!SetupSmootherPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found", app_ctx->name);
1092b730f8bSJeremy L Thompson     PetscCall((*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother));
110777ff853SJeremy L Thompson   }
111777ff853SJeremy L Thompson 
11262e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
113ccaff030SJeremy L Thompson   // Setup DM
114ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
115ccaff030SJeremy L Thompson   // Performance logging
1162b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup));
1172b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_dm_setup));
118ccaff030SJeremy L Thompson 
119ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
1202b730f8bSJeremy L Thompson   PetscCall(CreateDistributedDM(comm, app_ctx, &dm_orig));
121b68a8d79SJed Brown   VecType vectype;
122d1d35e2fSjeremylt   switch (mem_type_backend) {
1232b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
1242b730f8bSJeremy L Thompson       vectype = VECSTANDARD;
1252b730f8bSJeremy L Thompson       break;
126b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
127b68a8d79SJed Brown       const char *resolved;
128b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
129b68a8d79SJed Brown       if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
130b68a8d79SJed Brown       else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
131b68a8d79SJed Brown       else vectype = VECSTANDARD;
132b68a8d79SJed Brown     }
133b68a8d79SJed Brown   }
1342b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_orig, vectype));
1352b730f8bSJeremy L Thompson   PetscCall(DMPlexDistributeSetDefault(dm_orig, PETSC_FALSE));
1362b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm_orig));
137ccaff030SJeremy L Thompson 
138ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
1392b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &level_dms));
140d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
1412b730f8bSJeremy L Thompson     PetscCall(DMClone(dm_orig, &level_dms[level]));
1422b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm_orig, &vectype));
1432b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(level_dms[level], vectype));
1442b730f8bSJeremy L Thompson     PetscCall(SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level], PETSC_TRUE, num_comp_u));
145a3c02c40SJeremy L Thompson     // -- Label field components for viewing
146a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
147a3c02c40SJeremy L Thompson     PetscSection section;
1482b730f8bSJeremy L Thompson     PetscCall(DMGetLocalSection(level_dms[level], &section));
1492b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetFieldName(section, 0, "Displacement"));
1502b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX"));
1512b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY"));
1522b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"));
153a3c02c40SJeremy L Thompson   }
154a3c02c40SJeremy L Thompson 
155a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
1562b730f8bSJeremy L Thompson   PetscCall(DMClone(dm_orig, &dm_energy));
1572b730f8bSJeremy L Thompson   PetscCall(SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_e));
1582b730f8bSJeremy L Thompson   PetscCall(DMClone(dm_orig, &dm_diagnostic));
1592b730f8bSJeremy L Thompson   PetscCall(SetupDMByDegree(dm_diagnostic, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_u + num_comp_d));
1602b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_energy, vectype));
1612b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_diagnostic, vectype));
162a3c02c40SJeremy L Thompson   {
163a3c02c40SJeremy L Thompson     // -- Label field components for viewing
164a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
165a3c02c40SJeremy L Thompson     PetscSection section;
1662b730f8bSJeremy L Thompson     PetscCall(DMGetLocalSection(dm_diagnostic, &section));
1672b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetFieldName(section, 0, "Diagnostics"));
1682b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX"));
1692b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY"));
1702b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"));
1712b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 3, "Pressure"));
1722b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"));
1732b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 5, "TraceE2"));
1742b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 6, "detJ"));
1752b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"));
176ccaff030SJeremy L Thompson   }
177ccaff030SJeremy L Thompson 
178ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
179ccaff030SJeremy L Thompson   // Setup solution and work vectors
180ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
181ccaff030SJeremy L Thompson   // Allocate arrays
1822b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_g));
1832b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_loc));
1842b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_g_size));
1852b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_l_size));
1862b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_loc_size));
187ccaff030SJeremy L Thompson 
188ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
189d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
190ccaff030SJeremy L Thompson     // -- Create global unknown vector U
1912b730f8bSJeremy L Thompson     PetscCall(DMCreateGlobalVector(level_dms[level], &U_g[level]));
1922b730f8bSJeremy L Thompson     PetscCall(VecGetSize(U_g[level], &U_g_size[level]));
193ccaff030SJeremy L Thompson     // Note: Local size for matShell
1942b730f8bSJeremy L Thompson     PetscCall(VecGetLocalSize(U_g[level], &U_l_size[level]));
195ccaff030SJeremy L Thompson 
196d1d35e2fSjeremylt     // -- Create local unknown vector U_loc
1972b730f8bSJeremy L Thompson     PetscCall(DMCreateLocalVector(level_dms[level], &U_loc[level]));
198ccaff030SJeremy L Thompson     // Note: local size for libCEED
1992b730f8bSJeremy L Thompson     PetscCall(VecGetSize(U_loc[level], &U_loc_size[level]));
200ccaff030SJeremy L Thompson   }
201ccaff030SJeremy L Thompson 
202ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
2032b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &U));
2042b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &R));
2052b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &F));
2062b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_loc[fine_level], &R_loc));
2072b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_loc[fine_level], &F_loc));
208ccaff030SJeremy L Thompson 
209ccaff030SJeremy L Thompson   // Performance logging
2102b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
211ccaff030SJeremy L Thompson 
212ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
213ccaff030SJeremy L Thompson   // Set up libCEED
214ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
215ccaff030SJeremy L Thompson   // Performance logging
2162b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup));
2172b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_libceed_setup));
218ccaff030SJeremy L Thompson 
219ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
220d1d35e2fSjeremylt   CeedVector   force_ceed;
221ccaff030SJeremy L Thompson   CeedScalar  *f;
222d1d35e2fSjeremylt   PetscMemType force_mem_type;
223d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
2242b730f8bSJeremy L Thompson     PetscCall(VecGetArrayAndMemType(F_loc, &f, &force_mem_type));
225d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed);
226d1d35e2fSjeremylt     CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f);
227ccaff030SJeremy L Thompson   }
228ccaff030SJeremy L Thompson 
229fe394131Sjeremylt   // -- Create libCEED local Neumann BCs vector
230d1d35e2fSjeremylt   CeedVector   neumann_ceed;
231fe394131Sjeremylt   CeedScalar  *n;
232d1d35e2fSjeremylt   PetscMemType nummann_mem_type;
233d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
2342b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &neumann_bcs));
2352b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U_loc[fine_level], &bcs_loc));
2362b730f8bSJeremy L Thompson     PetscCall(VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type));
237d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed);
2382b730f8bSJeremy L Thompson     CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type), CEED_USE_POINTER, n);
239fe394131Sjeremylt   }
240fe394131Sjeremylt 
241ccaff030SJeremy L Thompson   // -- Setup libCEED objects
2422b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &ceed_data));
243d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
2442b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &ceed_data[fine_level]));
245ccaff030SJeremy L Thompson   {
2462b730f8bSJeremy L Thompson     PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx, CeedQFunctionContext, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector,
2472b730f8bSJeremy L Thompson                                             CeedVector, CeedData *);
2482b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupLibceedFineLevel, app_ctx->name, &SetupLibceedFineLevel));
2492b730f8bSJeremy L Thompson     if (!SetupLibceedFineLevel) SETERRQ(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found", app_ctx->name);
2502b730f8bSJeremy L Thompson     PetscCall((*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic, ceed, app_ctx, ctx_phys, fine_level, num_comp_u,
2512b730f8bSJeremy L Thompson                                        U_g_size[fine_level], U_loc_size[fine_level], force_ceed, neumann_ceed, ceed_data));
252ccaff030SJeremy L Thompson   }
253d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
254d1d35e2fSjeremylt   for (PetscInt level = num_levels - 2; level >= 0; level--) {
2552b730f8bSJeremy L Thompson     PetscCall(PetscCalloc1(1, &ceed_data[level]));
256d99fa3c5SJeremy L Thompson 
257d99fa3c5SJeremy L Thompson     // Get global communication restriction
2582b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_g[level + 1]));
2592b730f8bSJeremy L Thompson     PetscCall(VecSet(U_loc[level + 1], 1.0));
2602b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[level + 1], U_loc[level + 1], ADD_VALUES, U_g[level + 1]));
2612b730f8bSJeremy L Thompson     PetscCall(DMGlobalToLocal(level_dms[level + 1], U_g[level + 1], INSERT_VALUES, U_loc[level + 1]));
262d99fa3c5SJeremy L Thompson 
263d99fa3c5SJeremy L Thompson     // Place in libCEED array
264d99fa3c5SJeremy L Thompson     const PetscScalar *m;
265d1d35e2fSjeremylt     PetscMemType       m_mem_type;
2662b730f8bSJeremy L Thompson     PetscCall(VecGetArrayReadAndMemType(U_loc[level + 1], &m, &m_mem_type));
2672b730f8bSJeremy L Thompson     CeedVectorSetArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, (CeedScalar *)m);
268ccaff030SJeremy L Thompson 
2691c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
2702b730f8bSJeremy L Thompson     PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector, CeedData *);
2712b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupLibceedLevel, app_ctx->name, &SetupLibceedLevel));
2722b730f8bSJeremy L Thompson     if (!SetupLibceedLevel) SETERRQ(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found", app_ctx->name);
2732b730f8bSJeremy L Thompson     PetscCall((*SetupLibceedLevel)(level_dms[level], ceed, app_ctx, level, num_comp_u, U_g_size[level], U_loc_size[level],
2742b730f8bSJeremy L Thompson                                    ceed_data[level + 1]->x_ceed, ceed_data));
275d99fa3c5SJeremy L Thompson 
276d99fa3c5SJeremy L Thompson     // Restore PETSc vector
2772b730f8bSJeremy L Thompson     CeedVectorTakeArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), (CeedScalar **)&m);
2782b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayReadAndMemType(U_loc[level + 1], &m));
2792b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_g[level + 1]));
2802b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_loc[level + 1]));
281ccaff030SJeremy L Thompson   }
282ccaff030SJeremy L Thompson 
283ccaff030SJeremy L Thompson   // Performance logging
2842b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
285ccaff030SJeremy L Thompson 
286ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
287fe394131Sjeremylt   // Setup global forcing and Neumann BC vectors
288ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
2892b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(F));
290ccaff030SJeremy L Thompson 
291d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
292d1d35e2fSjeremylt     CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL);
2932b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayAndMemType(F_loc, &f));
2942b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F));
295d1d35e2fSjeremylt     CeedVectorDestroy(&force_ceed);
296ccaff030SJeremy L Thompson   }
297ccaff030SJeremy L Thompson 
298d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
2992b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(neumann_bcs));
300d1d35e2fSjeremylt     CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL);
3012b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayAndMemType(bcs_loc, &n));
3022b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs));
303d1d35e2fSjeremylt     CeedVectorDestroy(&neumann_ceed);
304fe394131Sjeremylt   }
305fe394131Sjeremylt 
306ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
307ccaff030SJeremy L Thompson   // Print problem summary
308ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
309d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
31077841947SLeila Ghaffari     char        hostname[PETSC_MAX_PATH_LEN];
3114d00b080SJeremy L Thompson     const char *used_resource;
3124d00b080SJeremy L Thompson     PetscMPIInt comm_size;
3134d00b080SJeremy L Thompson 
3144d00b080SJeremy L Thompson     CeedGetResource(ceed, &used_resource);
3152b730f8bSJeremy L Thompson     PetscCall(PetscGetHostName(hostname, sizeof hostname));
3162b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_size(comm, &comm_size));
317ccaff030SJeremy L Thompson 
3182b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
31917fd649aSJeremy L Thompson                           "\n-- Elasticity Example - libCEED + PETSc --\n"
32077841947SLeila Ghaffari                           "  MPI:\n"
32177841947SLeila Ghaffari                           "    Hostname                           : %s\n"
32277841947SLeila Ghaffari                           "    Total ranks                        : %d\n"
323ccaff030SJeremy L Thompson                           "  libCEED:\n"
32462e9c006SJeremy L Thompson                           "    libCEED Backend                    : %s\n"
325b68a8d79SJed Brown                           "    libCEED Backend MemType            : %s\n",
3264d00b080SJeremy L Thompson                           hostname, comm_size, used_resource, CeedMemTypes[mem_type_backend]));
327ccaff030SJeremy L Thompson 
32862e9c006SJeremy L Thompson     VecType vecType;
3292b730f8bSJeremy L Thompson     PetscCall(VecGetType(U, &vecType));
3302b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
33162e9c006SJeremy L Thompson                           "  PETSc:\n"
33262e9c006SJeremy L Thompson                           "    PETSc Vec Type                     : %s\n",
3332b730f8bSJeremy L Thompson                           vecType));
33462e9c006SJeremy L Thompson 
3352b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
336ccaff030SJeremy L Thompson                           "  Problem:\n"
337ccaff030SJeremy L Thompson                           "    Problem Name                       : %s\n"
338ccaff030SJeremy L Thompson                           "    Forcing Function                   : %s\n"
339ccaff030SJeremy L Thompson                           "  Mesh:\n"
340ccaff030SJeremy L Thompson                           "    File                               : %s\n"
3414d00b080SJeremy L Thompson                           "    Number of 1D Basis Nodes (p)       : %" PetscInt_FMT "\n"
3424d00b080SJeremy L Thompson                           "    Number of 1D Quadrature Points (q) : %" PetscInt_FMT "\n"
34308140895SJed Brown                           "    Global nodes                       : %" PetscInt_FMT "\n"
34408140895SJed Brown                           "    Owned nodes                        : %" PetscInt_FMT "\n"
34508140895SJed Brown                           "    DoF per node                       : %" PetscInt_FMT "\n"
346ccaff030SJeremy L Thompson                           "  Multigrid:\n"
347ccaff030SJeremy L Thompson                           "    Type                               : %s\n"
3484d00b080SJeremy L Thompson                           "    Number of Levels                   : %" PetscInt_FMT "\n",
3492b730f8bSJeremy L Thompson                           app_ctx->name_for_disp, forcing_types_for_disp[app_ctx->forcing_choice],
3502b730f8bSJeremy L Thompson                           app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh", app_ctx->degree + 1, app_ctx->degree + 1,
3512b730f8bSJeremy L Thompson                           U_g_size[fine_level] / num_comp_u, U_l_size[fine_level] / num_comp_u, num_comp_u,
3522b730f8bSJeremy L Thompson                           (app_ctx->degree == 1 && app_ctx->multigrid_choice != MULTIGRID_NONE) ? "Algebraic multigrid"
3532b730f8bSJeremy L Thompson                                                                                                 : multigrid_types_for_disp[app_ctx->multigrid_choice],
3542b730f8bSJeremy L Thompson                           (app_ctx->degree == 1 || app_ctx->multigrid_choice == MULTIGRID_NONE) ? 0 : num_levels));
355ccaff030SJeremy L Thompson 
356d1d35e2fSjeremylt     if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
357e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
3584d00b080SJeremy L Thompson         PetscInt level = i ? fine_level : 0;
3594d00b080SJeremy L Thompson 
3602b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
36108140895SJed Brown                               "    Level %" PetscInt_FMT " (%s):\n"
3624d00b080SJeremy L Thompson                               "      Number of 1D Basis Nodes (p)     : %" PetscInt_FMT "\n"
36308140895SJed Brown                               "      Global Nodes                     : %" PetscInt_FMT "\n"
36408140895SJed Brown                               "      Owned Nodes                      : %" PetscInt_FMT "\n",
3652b730f8bSJeremy L Thompson                               level, i ? "fine" : "coarse", app_ctx->level_degrees[level] + 1, U_g_size[level] / num_comp_u,
3662b730f8bSJeremy L Thompson                               U_l_size[level] / num_comp_u));
367ccaff030SJeremy L Thompson       }
368ccaff030SJeremy L Thompson     }
369ccaff030SJeremy L Thompson   }
370ccaff030SJeremy L Thompson 
371ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
372ccaff030SJeremy L Thompson   // Setup SNES
373ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
374ccaff030SJeremy L Thompson   // Performance logging
3752b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup));
3762b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_snes_setup));
377ccaff030SJeremy L Thompson 
378ccaff030SJeremy L Thompson   // Create SNES
3792b730f8bSJeremy L Thompson   PetscCall(SNESCreate(comm, &snes));
3802b730f8bSJeremy L Thompson   PetscCall(SNESSetDM(snes, level_dms[fine_level]));
381ccaff030SJeremy L Thompson 
382ccaff030SJeremy L Thompson   // -- Jacobian evaluators
3832b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &jacob_ctx));
3842b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &jacob_mat));
385d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
386ccaff030SJeremy L Thompson     // -- Jacobian context for level
3872b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &jacob_ctx[level]));
3882b730f8bSJeremy L Thompson     PetscCall(SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level], U_loc[level], ceed_data[level], ceed, ctx_phys, ctx_phys_smoother,
3892b730f8bSJeremy L Thompson                                jacob_ctx[level]));
390ccaff030SJeremy L Thompson 
391ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
3922b730f8bSJeremy L Thompson     PetscCall(MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level], U_g_size[level], jacob_ctx[level], &jacob_mat[level]));
3932b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_MULT, (void (*)(void))ApplyJacobian_Ceed));
3942b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL, (void (*)(void))GetDiag_Ceed));
3952b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(jacob_mat[level], vectype));
396ccaff030SJeremy L Thompson   }
397ea61e9acSJeremy L Thompson   // Note: FormJacobian updates Jacobian matrices on each level and assembles the Jpre matrix, if needed
3982b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &form_jacob_ctx));
399d1d35e2fSjeremylt   form_jacob_ctx->jacob_ctx  = jacob_ctx;
400d1d35e2fSjeremylt   form_jacob_ctx->num_levels = num_levels;
401d1d35e2fSjeremylt   form_jacob_ctx->jacob_mat  = jacob_mat;
402ccaff030SJeremy L Thompson 
403ccaff030SJeremy L Thompson   // -- Residual evaluation function
4042b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &res_ctx));
4052b730f8bSJeremy L Thompson   PetscCall(PetscMemcpy(res_ctx, jacob_ctx[fine_level], sizeof(*jacob_ctx[fine_level])));
4065754ecacSJeremy L Thompson   res_ctx->op = ceed_data[fine_level]->op_residual;
4075754ecacSJeremy L Thompson   res_ctx->qf = ceed_data[fine_level]->qf_residual;
4082b730f8bSJeremy L Thompson   if (app_ctx->bc_traction_count > 0) res_ctx->neumann_bcs = neumann_bcs;
4092b730f8bSJeremy L Thompson   else res_ctx->neumann_bcs = NULL;
4102b730f8bSJeremy L Thompson   PetscCall(SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx));
411ccaff030SJeremy L Thompson 
412ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
4132b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &prolong_restr_ctx));
4142b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &prolong_restr_mat));
415d1d35e2fSjeremylt   for (PetscInt level = 1; level < num_levels; level++) {
416ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
4172b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &prolong_restr_ctx[level]));
4182b730f8bSJeremy L Thompson     PetscCall(SetupProlongRestrictCtx(comm, app_ctx, level_dms[level - 1], level_dms[level], U_g[level], U_loc[level - 1], U_loc[level],
4192b730f8bSJeremy L Thompson                                       ceed_data[level - 1], ceed_data[level], ceed, prolong_restr_ctx[level]));
420ccaff030SJeremy L Thompson 
421ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
4222b730f8bSJeremy L Thompson     PetscCall(MatCreateShell(comm, U_l_size[level], U_l_size[level - 1], U_g_size[level], U_g_size[level - 1], prolong_restr_ctx[level],
4232b730f8bSJeremy L Thompson                              &prolong_restr_mat[level]));
424ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
4252b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT, (void (*)(void))Prolong_Ceed));
4262b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE, (void (*)(void))Restrict_Ceed));
4272b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(prolong_restr_mat[level], vectype));
428ccaff030SJeremy L Thompson   }
429ccaff030SJeremy L Thompson 
430ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
43117f0843fSJeremy L Thompson   // Setup for AMG coarse solve
432ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
433d1d35e2fSjeremylt   if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
434e3e3df41Sjeremylt     // -- Jacobian Matrix
4352b730f8bSJeremy L Thompson     PetscCall(DMCreateMatrix(level_dms[0], &jacob_mat_coarse));
436e3e3df41Sjeremylt 
437d1d35e2fSjeremylt     if (app_ctx->degree > 1) {
43817f0843fSJeremy L Thompson       // -- Assemble sparsity pattern
4393047f789SJeremy L Thompson       PetscCount             num_entries;
4404d00b080SJeremy L Thompson       PetscInt              *rows_petsc, *cols_petsc;
4414d00b080SJeremy L Thompson       CeedInt               *rows_ceed, *cols_ceed;
44217f0843fSJeremy L Thompson       CeedVector             coo_values;
44317f0843fSJeremy L Thompson       ISLocalToGlobalMapping ltog_row, ltog_col;
4444d00b080SJeremy L Thompson 
4454d00b080SJeremy L Thompson       CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries, &rows_ceed, &cols_ceed);
4464d00b080SJeremy L Thompson       PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
4474d00b080SJeremy L Thompson       PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
4482b730f8bSJeremy L Thompson       PetscCall(MatGetLocalToGlobalMapping(jacob_mat_coarse, &ltog_row, &ltog_col));
4494d00b080SJeremy L Thompson       PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows_petsc, rows_petsc));
4504d00b080SJeremy L Thompson       PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols_petsc, cols_petsc));
4514d00b080SJeremy L Thompson       PetscCall(MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows_petsc, cols_petsc));
4524d00b080SJeremy L Thompson       free(rows_petsc);
4534d00b080SJeremy L Thompson       free(cols_petsc);
45417f0843fSJeremy L Thompson       CeedVectorCreate(ceed, num_entries, &coo_values);
455ccaff030SJeremy L Thompson 
456d1d35e2fSjeremylt       // -- Update form_jacob_ctx
45717f0843fSJeremy L Thompson       form_jacob_ctx->coo_values       = coo_values;
45817f0843fSJeremy L Thompson       form_jacob_ctx->op_coarse        = ceed_data[0]->op_jacobian;
459d1d35e2fSjeremylt       form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse;
460e3e3df41Sjeremylt     }
461e3e3df41Sjeremylt   }
462e3e3df41Sjeremylt 
463e3e3df41Sjeremylt   // Set Jacobian function
464d1d35e2fSjeremylt   if (app_ctx->degree > 1) {
4652b730f8bSJeremy L Thompson     PetscCall(SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level], FormJacobian, form_jacob_ctx));
466e3e3df41Sjeremylt   } else {
4672b730f8bSJeremy L Thompson     PetscCall(SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse, SNESComputeJacobianDefaultColor, NULL));
468e3e3df41Sjeremylt   }
469ccaff030SJeremy L Thompson 
470ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
471ccaff030SJeremy L Thompson   // Setup KSP
472ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
473ccaff030SJeremy L Thompson   {
474ccaff030SJeremy L Thompson     PC  pc;
475ccaff030SJeremy L Thompson     KSP ksp;
476ccaff030SJeremy L Thompson 
477ccaff030SJeremy L Thompson     // -- KSP
4782b730f8bSJeremy L Thompson     PetscCall(SNESGetKSP(snes, &ksp));
4792b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
4802b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
4812b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
4822b730f8bSJeremy L Thompson     PetscCall(KSPSetOptionsPrefix(ksp, "outer_"));
483ccaff030SJeremy L Thompson 
484ccaff030SJeremy L Thompson     // -- Preconditioning
4852b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
4862b730f8bSJeremy L Thompson     PetscCall(PCSetDM(pc, level_dms[fine_level]));
4872b730f8bSJeremy L Thompson     PetscCall(PCSetOptionsPrefix(pc, "outer_"));
488ccaff030SJeremy L Thompson 
489d1d35e2fSjeremylt     if (app_ctx->multigrid_choice == MULTIGRID_NONE) {
490ccaff030SJeremy L Thompson       // ---- No Multigrid
4912b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
4922b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
493d1d35e2fSjeremylt     } else if (app_ctx->degree == 1) {
494f892e6d1Sjeremylt       // ---- AMG for degree 1
4952b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCGAMG));
496ccaff030SJeremy L Thompson     } else {
497ccaff030SJeremy L Thompson       // ---- PCMG
4982b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCMG));
499ccaff030SJeremy L Thompson 
500ccaff030SJeremy L Thompson       // ------ PCMG levels
5012b730f8bSJeremy L Thompson       PetscCall(PCMGSetLevels(pc, num_levels, NULL));
502d1d35e2fSjeremylt       for (PetscInt level = 0; level < num_levels; level++) {
503ccaff030SJeremy L Thompson         // -------- Smoother
504d1d35e2fSjeremylt         KSP ksp_smoother, ksp_est;
505d1d35e2fSjeremylt         PC  pc_smoother;
506ccaff030SJeremy L Thompson 
507ccaff030SJeremy L Thompson         // ---------- Smoother KSP
5082b730f8bSJeremy L Thompson         PetscCall(PCMGGetSmoother(pc, level, &ksp_smoother));
5092b730f8bSJeremy L Thompson         PetscCall(KSPSetDM(ksp_smoother, level_dms[level]));
5102b730f8bSJeremy L Thompson         PetscCall(KSPSetDMActive(ksp_smoother, PETSC_FALSE));
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson         // ---------- Chebyshev options
5132b730f8bSJeremy L Thompson         PetscCall(KSPSetType(ksp_smoother, KSPCHEBYSHEV));
5142b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1));
5152b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est));
5162b730f8bSJeremy L Thompson         PetscCall(KSPSetType(ksp_est, KSPCG));
5172b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE));
5182b730f8bSJeremy L Thompson         PetscCall(KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]));
519ccaff030SJeremy L Thompson 
520ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
5212b730f8bSJeremy L Thompson         PetscCall(KSPGetPC(ksp_smoother, &pc_smoother));
5222b730f8bSJeremy L Thompson         PetscCall(PCSetType(pc_smoother, PCJACOBI));
5232b730f8bSJeremy L Thompson         PetscCall(PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL));
524ccaff030SJeremy L Thompson 
525ccaff030SJeremy L Thompson         // -------- Work vector
526d1d35e2fSjeremylt         if (level != fine_level) {
5272b730f8bSJeremy L Thompson           PetscCall(PCMGSetX(pc, level, U_g[level]));
528ccaff030SJeremy L Thompson         }
529ccaff030SJeremy L Thompson 
530ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
531ccaff030SJeremy L Thompson         if (level > 0) {
5322b730f8bSJeremy L Thompson           PetscCall(PCMGSetInterpolation(pc, level, prolong_restr_mat[level]));
5332b730f8bSJeremy L Thompson           PetscCall(PCMGSetRestriction(pc, level, prolong_restr_mat[level]));
534ccaff030SJeremy L Thompson         }
535ccaff030SJeremy L Thompson       }
536ccaff030SJeremy L Thompson 
537ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
538d1d35e2fSjeremylt       KSP ksp_coarse;
539d1d35e2fSjeremylt       PC  pc_coarse;
540ccaff030SJeremy L Thompson 
541ccaff030SJeremy L Thompson       // -------- Coarse KSP
5422b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse));
5432b730f8bSJeremy L Thompson       PetscCall(KSPSetType(ksp_coarse, KSPPREONLY));
5442b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse));
5452b730f8bSJeremy L Thompson       PetscCall(KSPSetOptionsPrefix(ksp_coarse, "coarse_"));
546ccaff030SJeremy L Thompson 
547ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
5482b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(ksp_coarse, &pc_coarse));
5492b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc_coarse, PCGAMG));
5502b730f8bSJeremy L Thompson       PetscCall(PCSetOptionsPrefix(pc_coarse, "coarse_"));
551ccaff030SJeremy L Thompson 
5522b730f8bSJeremy L Thompson       PetscCall(KSPSetFromOptions(ksp_coarse));
5532b730f8bSJeremy L Thompson       PetscCall(PCSetFromOptions(pc_coarse));
554ccaff030SJeremy L Thompson 
555ccaff030SJeremy L Thompson       // ------ PCMG options
5562b730f8bSJeremy L Thompson       PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE));
5572b730f8bSJeremy L Thompson       PetscCall(PCMGSetNumberSmooth(pc, 3));
5582b730f8bSJeremy L Thompson       PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type));
559ccaff030SJeremy L Thompson     }
5602b730f8bSJeremy L Thompson     PetscCall(KSPSetFromOptions(ksp));
5612b730f8bSJeremy L Thompson     PetscCall(PCSetFromOptions(pc));
562ccaff030SJeremy L Thompson   }
563038d0cd7Sjeremylt   {
564038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
565d1d35e2fSjeremylt     SNESLineSearch line_search;
566481a4840SJed Brown 
5672b730f8bSJeremy L Thompson     PetscCall(SNESGetLineSearch(snes, &line_search));
5682b730f8bSJeremy L Thompson     PetscCall(SNESLineSearchSetType(line_search, SNESLINESEARCHCP));
569481a4840SJed Brown   }
570481a4840SJed Brown 
5712b730f8bSJeremy L Thompson   PetscCall(SNESSetFromOptions(snes));
572ccaff030SJeremy L Thompson 
573ccaff030SJeremy L Thompson   // Performance logging
5742b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
575ccaff030SJeremy L Thompson 
576ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
577ccaff030SJeremy L Thompson   // Set initial guess
578ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
5792b730f8bSJeremy L Thompson   PetscCall(PetscObjectSetName((PetscObject)U, ""));
5802b730f8bSJeremy L Thompson   PetscCall(VecSet(U, 0.0));
581ccaff030SJeremy L Thompson 
582ccaff030SJeremy L Thompson   // View solution
583d1d35e2fSjeremylt   if (app_ctx->view_soln) {
5842b730f8bSJeremy L Thompson     PetscCall(ViewSolution(comm, app_ctx, U, 0, 0.0));
585ccaff030SJeremy L Thompson   }
586ccaff030SJeremy L Thompson 
587ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
588ccaff030SJeremy L Thompson   // Solve SNES
589ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
590d1d35e2fSjeremylt   PetscBool snes_monitor = PETSC_FALSE;
5912b730f8bSJeremy L Thompson   PetscCall(PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor));
5925f24f028Sjeremylt 
593ccaff030SJeremy L Thompson   // Performance logging
5942b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve));
5952b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_snes_solve));
596ccaff030SJeremy L Thompson 
597ccaff030SJeremy L Thompson   // Timing
5982b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)snes));
599d1d35e2fSjeremylt   start_time = MPI_Wtime();
600ccaff030SJeremy L Thompson 
601ccaff030SJeremy L Thompson   // Solve for each load increment
6025f24f028Sjeremylt   PetscInt increment;
603d1d35e2fSjeremylt   for (increment = 1; increment <= app_ctx->num_increments; increment++) {
6045f24f028Sjeremylt     // -- Log increment count
605d1d35e2fSjeremylt     if (snes_monitor) {
6062b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "%" PetscInt_FMT " Load Increment\n", increment - 1));
6075f24f028Sjeremylt     }
6085f24f028Sjeremylt 
609ccaff030SJeremy L Thompson     // -- Scale the problem
610d1d35e2fSjeremylt     PetscScalar load_increment = 1.0 * increment / app_ctx->num_increments,
6112b730f8bSJeremy L Thompson                 scalingFactor  = load_increment / (increment == 1 ? 1 : res_ctx->load_increment);
612d1d35e2fSjeremylt     res_ctx->load_increment    = load_increment;
613d1d35e2fSjeremylt     if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) {
6142b730f8bSJeremy L Thompson       PetscCall(VecScale(F, scalingFactor));
615ccaff030SJeremy L Thompson     }
616ccaff030SJeremy L Thompson 
617ccaff030SJeremy L Thompson     // -- Solve
6182b730f8bSJeremy L Thompson     PetscCall(SNESSolve(snes, F, U));
619ccaff030SJeremy L Thompson 
620ccaff030SJeremy L Thompson     // -- View solution
621d1d35e2fSjeremylt     if (app_ctx->view_soln) {
6222b730f8bSJeremy L Thompson       PetscCall(ViewSolution(comm, app_ctx, U, increment, load_increment));
623ccaff030SJeremy L Thompson     }
624ccaff030SJeremy L Thompson 
625ccaff030SJeremy L Thompson     // -- Update SNES iteration count
626ccaff030SJeremy L Thompson     PetscInt its;
6272b730f8bSJeremy L Thompson     PetscCall(SNESGetIterationNumber(snes, &its));
628bf0f51feSjeremylt     snes_its += its;
6292b730f8bSJeremy L Thompson     PetscCall(SNESGetLinearSolveIterations(snes, &its));
630bf0f51feSjeremylt     ksp_its += its;
6313951731eSjeremylt 
6323951731eSjeremylt     // -- Check for divergence
6333951731eSjeremylt     SNESConvergedReason reason;
6342b730f8bSJeremy L Thompson     PetscCall(SNESGetConvergedReason(snes, &reason));
6352b730f8bSJeremy L Thompson     if (reason < 0) break;
6368355605fSKaren (Ren) Stengel     if (app_ctx->energy_viewer) {
6378355605fSKaren (Ren) Stengel       // -- Log strain energy for current load increment
6388355605fSKaren (Ren) Stengel       CeedScalar energy;
6392b730f8bSJeremy L Thompson       PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy));
6408355605fSKaren (Ren) Stengel 
6418355605fSKaren (Ren) Stengel       if (!app_ctx->test_mode) {
6428355605fSKaren (Ren) Stengel         // -- Output
6432b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm, "    Strain Energy                      : %.12e\n", energy));
6448355605fSKaren (Ren) Stengel       }
6452b730f8bSJeremy L Thompson       PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment, energy));
6468355605fSKaren (Ren) Stengel     }
647ccaff030SJeremy L Thompson   }
648ccaff030SJeremy L Thompson 
649ccaff030SJeremy L Thompson   // Timing
650d1d35e2fSjeremylt   elapsed_time = MPI_Wtime() - start_time;
651ccaff030SJeremy L Thompson 
652ccaff030SJeremy L Thompson   // Performance logging
6532b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
654ccaff030SJeremy L Thompson 
655ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
656ccaff030SJeremy L Thompson   // Output summary
657ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
658d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
659ccaff030SJeremy L Thompson     // -- SNES
660d1d35e2fSjeremylt     SNESType            snes_type;
661ccaff030SJeremy L Thompson     SNESConvergedReason reason;
662ccaff030SJeremy L Thompson     PetscReal           rnorm;
6632b730f8bSJeremy L Thompson     PetscCall(SNESGetType(snes, &snes_type));
6642b730f8bSJeremy L Thompson     PetscCall(SNESGetConvergedReason(snes, &reason));
6652b730f8bSJeremy L Thompson     PetscCall(SNESGetFunctionNorm(snes, &rnorm));
6662b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
667ccaff030SJeremy L Thompson                           "  SNES:\n"
668ccaff030SJeremy L Thompson                           "    SNES Type                          : %s\n"
669ccaff030SJeremy L Thompson                           "    SNES Convergence                   : %s\n"
670990fdeb6SJeremy L Thompson                           "    Number of Load Increments          : %" PetscInt_FMT "\n"
671990fdeb6SJeremy L Thompson                           "    Completed Load Increments          : %" PetscInt_FMT "\n"
67208140895SJed Brown                           "    Total SNES Iterations              : %" PetscInt_FMT "\n"
673ccaff030SJeremy L Thompson                           "    Final rnorm                        : %e\n",
6742b730f8bSJeremy L Thompson                           snes_type, SNESConvergedReasons[reason], app_ctx->num_increments, increment - 1, snes_its, (double)rnorm));
675ccaff030SJeremy L Thompson 
676ccaff030SJeremy L Thompson     // -- KSP
677ccaff030SJeremy L Thompson     KSP     ksp;
678d1d35e2fSjeremylt     KSPType ksp_type;
6792b730f8bSJeremy L Thompson     PetscCall(SNESGetKSP(snes, &ksp));
6802b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
6812b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
682ccaff030SJeremy L Thompson                           "  Linear Solver:\n"
6837418ca88SJeremy L Thompson                           "    KSP Type                           : %s\n"
68408140895SJed Brown                           "    Total KSP Iterations               : %" PetscInt_FMT "\n",
6852b730f8bSJeremy L Thompson                           ksp_type, ksp_its));
686ccaff030SJeremy L Thompson 
687ccaff030SJeremy L Thompson     // -- PC
688ccaff030SJeremy L Thompson     PC     pc;
689d1d35e2fSjeremylt     PCType pc_type;
6902b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
6912b730f8bSJeremy L Thompson     PetscCall(PCGetType(pc, &pc_type));
6922b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "    PC Type                            : %s\n", pc_type));
693e3e3df41Sjeremylt 
694d1d35e2fSjeremylt     if (!strcmp(pc_type, PCMG)) {
695d1d35e2fSjeremylt       PCMGType pcmg_type;
6962b730f8bSJeremy L Thompson       PetscCall(PCMGGetType(pc, &pcmg_type));
6972b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
698ccaff030SJeremy L Thompson                             "  P-Multigrid:\n"
699ccaff030SJeremy L Thompson                             "    PCMG Type                          : %s\n"
700ccaff030SJeremy L Thompson                             "    PCMG Cycle Type                    : %s\n",
7012b730f8bSJeremy L Thompson                             PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type]));
702ccaff030SJeremy L Thompson 
703ccaff030SJeremy L Thompson       // -- Coarse Solve
704d1d35e2fSjeremylt       KSP    ksp_coarse;
705d1d35e2fSjeremylt       PC     pc_coarse;
706d1d35e2fSjeremylt       PCType pc_type;
707ccaff030SJeremy L Thompson 
7082b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse));
7092b730f8bSJeremy L Thompson       PetscCall(KSPGetType(ksp_coarse, &ksp_type));
7102b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(ksp_coarse, &pc_coarse));
7112b730f8bSJeremy L Thompson       PetscCall(PCGetType(pc_coarse, &pc_type));
7122b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
713ccaff030SJeremy L Thompson                             "    Coarse Solve:\n"
714ccaff030SJeremy L Thompson                             "      KSP Type                         : %s\n"
715ccaff030SJeremy L Thompson                             "      PC Type                          : %s\n",
7162b730f8bSJeremy L Thompson                             ksp_type, pc_type));
717ccaff030SJeremy L Thompson     }
718ccaff030SJeremy L Thompson   }
719ccaff030SJeremy L Thompson 
720ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
721ccaff030SJeremy L Thompson   // Compute solve time
722ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
723d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
7242b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm));
7252b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm));
7262b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
727ccaff030SJeremy L Thompson                           "  Performance:\n"
7287418ca88SJeremy L Thompson                           "    SNES Solve Time                    : %g (%g) sec\n"
7297418ca88SJeremy L Thompson                           "    DoFs/Sec in SNES                   : %g (%g) million\n",
7302b730f8bSJeremy L Thompson                           max_time, min_time, 1e-6 * U_g_size[fine_level] * ksp_its / max_time, 1e-6 * U_g_size[fine_level] * ksp_its / min_time));
731ccaff030SJeremy L Thompson   }
732ccaff030SJeremy L Thompson 
733ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
734ccaff030SJeremy L Thompson   // Compute error
735ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
736d1d35e2fSjeremylt   if (app_ctx->forcing_choice == FORCE_MMS) {
737d1d35e2fSjeremylt     CeedScalar        l2_error = 1., l2_U_norm = 1.;
738d1d35e2fSjeremylt     const CeedScalar *true_array;
739d1d35e2fSjeremylt     Vec               error_vec, true_vec;
740ccaff030SJeremy L Thompson 
741ccaff030SJeremy L Thompson     // -- Work vectors
7422b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &error_vec));
7432b730f8bSJeremy L Thompson     PetscCall(VecSet(error_vec, 0.0));
7442b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &true_vec));
7452b730f8bSJeremy L Thompson     PetscCall(VecSet(true_vec, 0.0));
746ccaff030SJeremy L Thompson 
747ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
7482b730f8bSJeremy L Thompson     CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln, CEED_MEM_HOST, &true_array);
7492b730f8bSJeremy L Thompson     PetscCall(VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array));
7502b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec));
7512b730f8bSJeremy L Thompson     PetscCall(VecResetArray(res_ctx->Y_loc));
752d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array);
753ccaff030SJeremy L Thompson 
754ccaff030SJeremy L Thompson     // -- Compute L2 error
7552b730f8bSJeremy L Thompson     PetscCall(VecWAXPY(error_vec, -1.0, U, true_vec));
7562b730f8bSJeremy L Thompson     PetscCall(VecNorm(error_vec, NORM_2, &l2_error));
7572b730f8bSJeremy L Thompson     PetscCall(VecNorm(U, NORM_2, &l2_U_norm));
758d1d35e2fSjeremylt     l2_error /= l2_U_norm;
759ccaff030SJeremy L Thompson 
760ccaff030SJeremy L Thompson     // -- Output
761d1d35e2fSjeremylt     if (!app_ctx->test_mode || l2_error > 0.05) {
7622b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "    L2 Error                           : %e\n", l2_error));
763ccaff030SJeremy L Thompson     }
764ccaff030SJeremy L Thompson 
765ccaff030SJeremy L Thompson     // -- Cleanup
7662b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&error_vec));
7672b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&true_vec));
7682d93065eSjeremylt   }
7692d93065eSjeremylt 
7702d93065eSjeremylt   // ---------------------------------------------------------------------------
7712d93065eSjeremylt   // Compute energy
7722d93065eSjeremylt   // ---------------------------------------------------------------------------
7738355605fSKaren (Ren) Stengel   PetscReal energy;
7742b730f8bSJeremy L Thompson   PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy));
7758355605fSKaren (Ren) Stengel   if (!app_ctx->test_mode) {
7762d93065eSjeremylt     // -- Output
7772b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "    Strain Energy                      : %.12e\n", energy));
778ccaff030SJeremy L Thompson   }
7792b730f8bSJeremy L Thompson   PetscCall(RegressionTests_solids(app_ctx, energy));
780ccaff030SJeremy L Thompson 
781ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
7825c25879aSJeremy L Thompson   // Output diagnostic quantities
7835c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
784d1d35e2fSjeremylt   if (app_ctx->view_soln || app_ctx->view_final_soln) {
7855c25879aSJeremy L Thompson     // -- Setup context
786d1d35e2fSjeremylt     UserMult diagnostic_ctx;
7872b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &diagnostic_ctx));
7882b730f8bSJeremy L Thompson     PetscCall(PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)));
789d1d35e2fSjeremylt     diagnostic_ctx->dm = dm_diagnostic;
790d1d35e2fSjeremylt     diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic;
7915c25879aSJeremy L Thompson 
7925c25879aSJeremy L Thompson     // -- Compute and output
7932b730f8bSJeremy L Thompson     PetscCall(ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx, app_ctx, U, ceed_data[fine_level]->elem_restr_diagnostic));
7945c25879aSJeremy L Thompson 
7955c25879aSJeremy L Thompson     // -- Cleanup
7962b730f8bSJeremy L Thompson     PetscCall(PetscFree(diagnostic_ctx));
7975c25879aSJeremy L Thompson   }
7985c25879aSJeremy L Thompson 
7995c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
800ccaff030SJeremy L Thompson   // Free objects
801ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
802ccaff030SJeremy L Thompson   // Data in arrays per level
803d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
804ccaff030SJeremy L Thompson     // Vectors
8052b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&U_g[level]));
8062b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&U_loc[level]));
807ccaff030SJeremy L Thompson 
808ccaff030SJeremy L Thompson     // Jacobian matrix and data
8092b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&jacob_ctx[level]->Y_loc));
8102b730f8bSJeremy L Thompson     PetscCall(MatDestroy(&jacob_mat[level]));
8112b730f8bSJeremy L Thompson     PetscCall(PetscFree(jacob_ctx[level]));
812ccaff030SJeremy L Thompson 
813ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
814ccaff030SJeremy L Thompson     if (level > 0) {
8152b730f8bSJeremy L Thompson       PetscCall(PetscFree(prolong_restr_ctx[level]));
8162b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&prolong_restr_mat[level]));
817ccaff030SJeremy L Thompson     }
818ccaff030SJeremy L Thompson 
819ccaff030SJeremy L Thompson     // DM
8202b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&level_dms[level]));
821ccaff030SJeremy L Thompson 
822ccaff030SJeremy L Thompson     // libCEED objects
8232b730f8bSJeremy L Thompson     PetscCall(CeedDataDestroy(level, ceed_data[level]));
824ccaff030SJeremy L Thompson   }
825ccaff030SJeremy L Thompson 
8262b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&app_ctx->energy_viewer));
8278355605fSKaren (Ren) Stengel 
828ccaff030SJeremy L Thompson   // Arrays
8292b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_g));
8302b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_loc));
8312b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_g_size));
8322b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_l_size));
8332b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_loc_size));
8342b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_ctx));
8352b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_mat));
8362b730f8bSJeremy L Thompson   PetscCall(PetscFree(prolong_restr_ctx));
8372b730f8bSJeremy L Thompson   PetscCall(PetscFree(prolong_restr_mat));
8382b730f8bSJeremy L Thompson   PetscCall(PetscFree(app_ctx->level_degrees));
8392b730f8bSJeremy L Thompson   PetscCall(PetscFree(ceed_data));
840ccaff030SJeremy L Thompson 
841ccaff030SJeremy L Thompson   // libCEED objects
84217f0843fSJeremy L Thompson   CeedVectorDestroy(&form_jacob_ctx->coo_values);
843d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys);
844d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys_smoother);
845ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
846ccaff030SJeremy L Thompson 
847ccaff030SJeremy L Thompson   // PETSc objects
8482b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&U));
8492b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&R));
8502b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&R_loc));
8512b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&F));
8522b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&F_loc));
8532b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&neumann_bcs));
8542b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&bcs_loc));
8552b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&jacob_mat_coarse));
8562b730f8bSJeremy L Thompson   PetscCall(SNESDestroy(&snes));
8572b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_orig));
8582b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_energy));
8592b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_diagnostic));
8602b730f8bSJeremy L Thompson   PetscCall(PetscFree(level_dms));
861ccaff030SJeremy L Thompson 
8625754ecacSJeremy L Thompson   // -- Function list
8632b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupPhysics));
8642b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupSmootherPhysics));
8652b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel));
8662b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedLevel));
8675754ecacSJeremy L Thompson 
868ccaff030SJeremy L Thompson   // Structs
8692b730f8bSJeremy L Thompson   PetscCall(PetscFree(res_ctx));
8702b730f8bSJeremy L Thompson   PetscCall(PetscFree(form_jacob_ctx));
8712b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_coarse_ctx));
8722b730f8bSJeremy L Thompson   PetscCall(PetscFree(app_ctx));
8732b730f8bSJeremy L Thompson   PetscCall(PetscFree(problem_functions));
8742b730f8bSJeremy L Thompson   PetscCall(PetscFree(units));
875ccaff030SJeremy L Thompson 
876ccaff030SJeremy L Thompson   return PetscFinalize();
877ccaff030SJeremy L Thompson }
878