xref: /libCEED/examples/solids/elasticity.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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 //
10ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve
11ccaff030SJeremy L Thompson //   elasticity problems.
12ccaff030SJeremy L Thompson //
13ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex.
14ccaff030SJeremy L Thompson //
15ccaff030SJeremy L Thompson // Build with:
16ccaff030SJeremy L Thompson //
17ccaff030SJeremy L Thompson //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18ccaff030SJeremy L Thompson //
19ccaff030SJeremy L Thompson // Sample runs:
20ccaff030SJeremy L Thompson //
21eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms
22*2b730f8bSJeremy 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
23*2b730f8bSJeremy L Thompson //     /cpu/self
24*2b730f8bSJeremy 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
25*2b730f8bSJeremy L Thompson //     -ceed /gpu/cuda
26ccaff030SJeremy L Thompson //
27ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
28ccaff030SJeremy L Thompson //
298355605fSKaren (Ren) Stengel //TESTARGS(name="solids-Linear-MMS") -ceed {ceed_resource} -test -degree 3 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3
308355605fSKaren (Ren) Stengel //TESTARGS(name="solids-NH1-1") -ceed {ceed_resource} -test -problem FSInitial-NH1 -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
318355605fSKaren (Ren) Stengel //TESTARGS(name="solids-MR1-1") -ceed {ceed_resource} -test -problem FSInitial-MR1 -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
32ccaff030SJeremy L Thompson 
33ccaff030SJeremy L Thompson /// @file
34ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson #include "elasticity.h"
39ccaff030SJeremy L Thompson 
40ccaff030SJeremy L Thompson int main(int argc, char **argv) {
41ccaff030SJeremy L Thompson   MPI_Comm comm;
42ccaff030SJeremy L Thompson   // Context structs
43d1d35e2fSjeremylt   AppCtx           app_ctx;            // Contains problem options
445754ecacSJeremy L Thompson   ProblemFunctions problem_functions;  // Setup functions for each problem
45ccaff030SJeremy L Thompson   Units            units;              // Contains units scaling
46ccaff030SJeremy L Thompson   // PETSc objects
47*2b730f8bSJeremy L Thompson   PetscLogStage stage_dm_setup, stage_libceed_setup, stage_snes_setup, stage_snes_solve;
48d1d35e2fSjeremylt   DM            dm_orig;                   // Distributed DM to clone
49d1d35e2fSjeremylt   DM            dm_energy, dm_diagnostic;  // DMs for postprocessing
50d1d35e2fSjeremylt   DM           *level_dms;
51d1d35e2fSjeremylt   Vec           U, *U_g, *U_loc;     // U: solution, R: residual, F: forcing
52d1d35e2fSjeremylt   Vec           R, R_loc, F, F_loc;  // g: global, loc: local
53d1d35e2fSjeremylt   Vec           neumann_bcs = NULL, bcs_loc = NULL;
5417f0843fSJeremy L Thompson   SNES          snes;
55d1d35e2fSjeremylt   Mat          *jacob_mat, jacob_mat_coarse, *prolong_restr_mat;
56ccaff030SJeremy L Thompson   // PETSc data
57d1d35e2fSjeremylt   UserMult              res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx;
58d1d35e2fSjeremylt   FormJacobCtx          form_jacob_ctx;
59d1d35e2fSjeremylt   UserMultProlongRestr *prolong_restr_ctx;
60d1d35e2fSjeremylt   PCMGCycleType         pcmg_cycle_type = PC_MG_CYCLE_V;
61ccaff030SJeremy L Thompson   // libCEED objects
62d99fa3c5SJeremy L Thompson   Ceed                 ceed;
63d1d35e2fSjeremylt   CeedData            *ceed_data;
64d1d35e2fSjeremylt   CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL;
65ccaff030SJeremy L Thompson   // Parameters
66d1d35e2fSjeremylt   PetscInt  num_comp_u = 3;                  // 3 DoFs in 3D
67d1d35e2fSjeremylt   PetscInt  num_comp_e = 1, num_comp_d = 5;  // 1 energy output, 5 diagnostic
68d1d35e2fSjeremylt   PetscInt  num_levels = 1, fine_level = 0;
69bf0f51feSjeremylt   PetscInt *U_g_size, *U_l_size, *U_loc_size;
70bf0f51feSjeremylt   PetscInt  snes_its = 0, ksp_its = 0;
71d1d35e2fSjeremylt   double    start_time, elapsed_time, min_time, max_time;
72ccaff030SJeremy L Thompson 
73*2b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
74ccaff030SJeremy L Thompson 
75ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
76ccaff030SJeremy L Thompson   // Process command line options
77ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
78ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
79ccaff030SJeremy L Thompson 
80ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
81*2b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &app_ctx));
82*2b730f8bSJeremy L Thompson   PetscCall(ProcessCommandLineOptions(comm, app_ctx));
83*2b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &problem_functions));
84*2b730f8bSJeremy L Thompson   PetscCall(RegisterProblems(problem_functions));
85d1d35e2fSjeremylt   num_levels = app_ctx->num_levels;
86d1d35e2fSjeremylt   fine_level = num_levels - 1;
87ccaff030SJeremy L Thompson 
88ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
89d1d35e2fSjeremylt   // Initialize libCEED
9062e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
91d1d35e2fSjeremylt   // Initialize backend
92d1d35e2fSjeremylt   CeedInit(app_ctx->ceed_resource, &ceed);
9362e9c006SJeremy L Thompson 
9462e9c006SJeremy L Thompson   // Check preferred MemType
95d1d35e2fSjeremylt   CeedMemType mem_type_backend;
96d1d35e2fSjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
975754ecacSJeremy L Thompson   // Setup physics context and wrap in libCEED object
985754ecacSJeremy L Thompson   {
995754ecacSJeremy L Thompson     PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *);
100*2b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name, &SetupPhysics));
101*2b730f8bSJeremy L Thompson     if (!SetupPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found", app_ctx->name);
102*2b730f8bSJeremy L Thompson     PetscCall((*SetupPhysics)(comm, ceed, &units, &ctx_phys));
103*2b730f8bSJeremy L Thompson     PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext, CeedQFunctionContext *);
104*2b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupSmootherPhysics, app_ctx->name, &SetupSmootherPhysics));
105*2b730f8bSJeremy L Thompson     if (!SetupSmootherPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found", app_ctx->name);
106*2b730f8bSJeremy L Thompson     PetscCall((*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother));
107777ff853SJeremy L Thompson   }
108777ff853SJeremy L Thompson 
10962e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
110ccaff030SJeremy L Thompson   // Setup DM
111ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
112ccaff030SJeremy L Thompson   // Performance logging
113*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup));
114*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_dm_setup));
115ccaff030SJeremy L Thompson 
116ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
117*2b730f8bSJeremy L Thompson   PetscCall(CreateDistributedDM(comm, app_ctx, &dm_orig));
118b68a8d79SJed Brown   VecType vectype;
119d1d35e2fSjeremylt   switch (mem_type_backend) {
120*2b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
121*2b730f8bSJeremy L Thompson       vectype = VECSTANDARD;
122*2b730f8bSJeremy L Thompson       break;
123b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
124b68a8d79SJed Brown       const char *resolved;
125b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
126b68a8d79SJed Brown       if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
127b68a8d79SJed Brown       else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
128b68a8d79SJed Brown       else vectype = VECSTANDARD;
129b68a8d79SJed Brown     }
130b68a8d79SJed Brown   }
131*2b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_orig, vectype));
132*2b730f8bSJeremy L Thompson   PetscCall(DMPlexDistributeSetDefault(dm_orig, PETSC_FALSE));
133*2b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm_orig));
134ccaff030SJeremy L Thompson 
135ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
136*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &level_dms));
137d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
138*2b730f8bSJeremy L Thompson     PetscCall(DMClone(dm_orig, &level_dms[level]));
139*2b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm_orig, &vectype));
140*2b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(level_dms[level], vectype));
141*2b730f8bSJeremy L Thompson     PetscCall(SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level], PETSC_TRUE, num_comp_u));
142a3c02c40SJeremy L Thompson     // -- Label field components for viewing
143a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
144a3c02c40SJeremy L Thompson     PetscSection section;
145*2b730f8bSJeremy L Thompson     PetscCall(DMGetLocalSection(level_dms[level], &section));
146*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetFieldName(section, 0, "Displacement"));
147*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX"));
148*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY"));
149*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"));
150a3c02c40SJeremy L Thompson   }
151a3c02c40SJeremy L Thompson 
152a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
153*2b730f8bSJeremy L Thompson   PetscCall(DMClone(dm_orig, &dm_energy));
154*2b730f8bSJeremy L Thompson   PetscCall(SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_e));
155*2b730f8bSJeremy L Thompson   PetscCall(DMClone(dm_orig, &dm_diagnostic));
156*2b730f8bSJeremy L Thompson   PetscCall(SetupDMByDegree(dm_diagnostic, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_u + num_comp_d));
157*2b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_energy, vectype));
158*2b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm_diagnostic, vectype));
159a3c02c40SJeremy L Thompson   {
160a3c02c40SJeremy L Thompson     // -- Label field components for viewing
161a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
162a3c02c40SJeremy L Thompson     PetscSection section;
163*2b730f8bSJeremy L Thompson     PetscCall(DMGetLocalSection(dm_diagnostic, &section));
164*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetFieldName(section, 0, "Diagnostics"));
165*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX"));
166*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY"));
167*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"));
168*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 3, "Pressure"));
169*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"));
170*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 5, "TraceE2"));
171*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 6, "detJ"));
172*2b730f8bSJeremy L Thompson     PetscCall(PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"));
173ccaff030SJeremy L Thompson   }
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
176ccaff030SJeremy L Thompson   // Setup solution and work vectors
177ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
178ccaff030SJeremy L Thompson   // Allocate arrays
179*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_g));
180*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_loc));
181*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_g_size));
182*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_l_size));
183*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &U_loc_size));
184ccaff030SJeremy L Thompson 
185ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
186d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
187ccaff030SJeremy L Thompson     // -- Create global unknown vector U
188*2b730f8bSJeremy L Thompson     PetscCall(DMCreateGlobalVector(level_dms[level], &U_g[level]));
189*2b730f8bSJeremy L Thompson     PetscCall(VecGetSize(U_g[level], &U_g_size[level]));
190ccaff030SJeremy L Thompson     // Note: Local size for matShell
191*2b730f8bSJeremy L Thompson     PetscCall(VecGetLocalSize(U_g[level], &U_l_size[level]));
192ccaff030SJeremy L Thompson 
193d1d35e2fSjeremylt     // -- Create local unknown vector U_loc
194*2b730f8bSJeremy L Thompson     PetscCall(DMCreateLocalVector(level_dms[level], &U_loc[level]));
195ccaff030SJeremy L Thompson     // Note: local size for libCEED
196*2b730f8bSJeremy L Thompson     PetscCall(VecGetSize(U_loc[level], &U_loc_size[level]));
197ccaff030SJeremy L Thompson   }
198ccaff030SJeremy L Thompson 
199ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
200*2b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &U));
201*2b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &R));
202*2b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_g[fine_level], &F));
203*2b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_loc[fine_level], &R_loc));
204*2b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(U_loc[fine_level], &F_loc));
205ccaff030SJeremy L Thompson 
206ccaff030SJeremy L Thompson   // Performance logging
207*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
208ccaff030SJeremy L Thompson 
209ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
210ccaff030SJeremy L Thompson   // Set up libCEED
211ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
212ccaff030SJeremy L Thompson   // Performance logging
213*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup));
214*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_libceed_setup));
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
217d1d35e2fSjeremylt   CeedVector   force_ceed;
218ccaff030SJeremy L Thompson   CeedScalar  *f;
219d1d35e2fSjeremylt   PetscMemType force_mem_type;
220d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
221*2b730f8bSJeremy L Thompson     PetscCall(VecGetArrayAndMemType(F_loc, &f, &force_mem_type));
222d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed);
223d1d35e2fSjeremylt     CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f);
224ccaff030SJeremy L Thompson   }
225ccaff030SJeremy L Thompson 
226fe394131Sjeremylt   // -- Create libCEED local Neumann BCs vector
227d1d35e2fSjeremylt   CeedVector   neumann_ceed;
228fe394131Sjeremylt   CeedScalar  *n;
229d1d35e2fSjeremylt   PetscMemType nummann_mem_type;
230d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
231*2b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &neumann_bcs));
232*2b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U_loc[fine_level], &bcs_loc));
233*2b730f8bSJeremy L Thompson     PetscCall(VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type));
234d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed);
235*2b730f8bSJeremy L Thompson     CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type), CEED_USE_POINTER, n);
236fe394131Sjeremylt   }
237fe394131Sjeremylt 
238ccaff030SJeremy L Thompson   // -- Setup libCEED objects
239*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &ceed_data));
240d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
241*2b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &ceed_data[fine_level]));
242ccaff030SJeremy L Thompson   {
243*2b730f8bSJeremy L Thompson     PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx, CeedQFunctionContext, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector,
244*2b730f8bSJeremy L Thompson                                             CeedVector, CeedData *);
245*2b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupLibceedFineLevel, app_ctx->name, &SetupLibceedFineLevel));
246*2b730f8bSJeremy L Thompson     if (!SetupLibceedFineLevel) SETERRQ(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found", app_ctx->name);
247*2b730f8bSJeremy L Thompson     PetscCall((*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic, ceed, app_ctx, ctx_phys, fine_level, num_comp_u,
248*2b730f8bSJeremy L Thompson                                        U_g_size[fine_level], U_loc_size[fine_level], force_ceed, neumann_ceed, ceed_data));
249ccaff030SJeremy L Thompson   }
250d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
251d1d35e2fSjeremylt   for (PetscInt level = num_levels - 2; level >= 0; level--) {
252*2b730f8bSJeremy L Thompson     PetscCall(PetscCalloc1(1, &ceed_data[level]));
253d99fa3c5SJeremy L Thompson 
254d99fa3c5SJeremy L Thompson     // Get global communication restriction
255*2b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_g[level + 1]));
256*2b730f8bSJeremy L Thompson     PetscCall(VecSet(U_loc[level + 1], 1.0));
257*2b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[level + 1], U_loc[level + 1], ADD_VALUES, U_g[level + 1]));
258*2b730f8bSJeremy L Thompson     PetscCall(DMGlobalToLocal(level_dms[level + 1], U_g[level + 1], INSERT_VALUES, U_loc[level + 1]));
259d99fa3c5SJeremy L Thompson 
260d99fa3c5SJeremy L Thompson     // Place in libCEED array
261d99fa3c5SJeremy L Thompson     const PetscScalar *m;
262d1d35e2fSjeremylt     PetscMemType       m_mem_type;
263*2b730f8bSJeremy L Thompson     PetscCall(VecGetArrayReadAndMemType(U_loc[level + 1], &m, &m_mem_type));
264*2b730f8bSJeremy L Thompson     CeedVectorSetArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, (CeedScalar *)m);
265ccaff030SJeremy L Thompson 
2661c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
267*2b730f8bSJeremy L Thompson     PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector, CeedData *);
268*2b730f8bSJeremy L Thompson     PetscCall(PetscFunctionListFind(problem_functions->setupLibceedLevel, app_ctx->name, &SetupLibceedLevel));
269*2b730f8bSJeremy L Thompson     if (!SetupLibceedLevel) SETERRQ(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found", app_ctx->name);
270*2b730f8bSJeremy L Thompson     PetscCall((*SetupLibceedLevel)(level_dms[level], ceed, app_ctx, level, num_comp_u, U_g_size[level], U_loc_size[level],
271*2b730f8bSJeremy L Thompson                                    ceed_data[level + 1]->x_ceed, ceed_data));
272d99fa3c5SJeremy L Thompson 
273d99fa3c5SJeremy L Thompson     // Restore PETSc vector
274*2b730f8bSJeremy L Thompson     CeedVectorTakeArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), (CeedScalar **)&m);
275*2b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayReadAndMemType(U_loc[level + 1], &m));
276*2b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_g[level + 1]));
277*2b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(U_loc[level + 1]));
278ccaff030SJeremy L Thompson   }
279ccaff030SJeremy L Thompson 
280ccaff030SJeremy L Thompson   // Performance logging
281*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
282ccaff030SJeremy L Thompson 
283ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
284fe394131Sjeremylt   // Setup global forcing and Neumann BC vectors
285ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
286*2b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(F));
287ccaff030SJeremy L Thompson 
288d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
289d1d35e2fSjeremylt     CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL);
290*2b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayAndMemType(F_loc, &f));
291*2b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F));
292d1d35e2fSjeremylt     CeedVectorDestroy(&force_ceed);
293ccaff030SJeremy L Thompson   }
294ccaff030SJeremy L Thompson 
295d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
296*2b730f8bSJeremy L Thompson     PetscCall(VecZeroEntries(neumann_bcs));
297d1d35e2fSjeremylt     CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL);
298*2b730f8bSJeremy L Thompson     PetscCall(VecRestoreArrayAndMemType(bcs_loc, &n));
299*2b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs));
300d1d35e2fSjeremylt     CeedVectorDestroy(&neumann_ceed);
301fe394131Sjeremylt   }
302fe394131Sjeremylt 
303ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
304ccaff030SJeremy L Thompson   // Print problem summary
305ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
306d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
307ccaff030SJeremy L Thompson     const char *usedresource;
308ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
30977841947SLeila Ghaffari     char hostname[PETSC_MAX_PATH_LEN];
310*2b730f8bSJeremy L Thompson     PetscCall(PetscGetHostName(hostname, sizeof hostname));
31177841947SLeila Ghaffari     PetscInt comm_size;
312*2b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_size(comm, &comm_size));
313ccaff030SJeremy L Thompson 
314*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
31517fd649aSJeremy L Thompson                           "\n-- Elasticity Example - libCEED + PETSc --\n"
31677841947SLeila Ghaffari                           "  MPI:\n"
31777841947SLeila Ghaffari                           "    Hostname                           : %s\n"
31877841947SLeila Ghaffari                           "    Total ranks                        : %d\n"
319ccaff030SJeremy L Thompson                           "  libCEED:\n"
32062e9c006SJeremy L Thompson                           "    libCEED Backend                    : %s\n"
321b68a8d79SJed Brown                           "    libCEED Backend MemType            : %s\n",
322*2b730f8bSJeremy L Thompson                           hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend]));
323ccaff030SJeremy L Thompson 
32462e9c006SJeremy L Thompson     VecType vecType;
325*2b730f8bSJeremy L Thompson     PetscCall(VecGetType(U, &vecType));
326*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
32762e9c006SJeremy L Thompson                           "  PETSc:\n"
32862e9c006SJeremy L Thompson                           "    PETSc Vec Type                     : %s\n",
329*2b730f8bSJeremy L Thompson                           vecType));
33062e9c006SJeremy L Thompson 
331*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
332ccaff030SJeremy L Thompson                           "  Problem:\n"
333ccaff030SJeremy L Thompson                           "    Problem Name                       : %s\n"
334ccaff030SJeremy L Thompson                           "    Forcing Function                   : %s\n"
335ccaff030SJeremy L Thompson                           "  Mesh:\n"
336ccaff030SJeremy L Thompson                           "    File                               : %s\n"
337990fdeb6SJeremy L Thompson                           "    Number of 1D Basis Nodes (p)       : %" CeedInt_FMT "\n"
338990fdeb6SJeremy L Thompson                           "    Number of 1D Quadrature Points (q) : %" CeedInt_FMT "\n"
33908140895SJed Brown                           "    Global nodes                       : %" PetscInt_FMT "\n"
34008140895SJed Brown                           "    Owned nodes                        : %" PetscInt_FMT "\n"
34108140895SJed Brown                           "    DoF per node                       : %" PetscInt_FMT "\n"
342ccaff030SJeremy L Thompson                           "  Multigrid:\n"
343ccaff030SJeremy L Thompson                           "    Type                               : %s\n"
344990fdeb6SJeremy L Thompson                           "    Number of Levels                   : %" CeedInt_FMT "\n",
345*2b730f8bSJeremy L Thompson                           app_ctx->name_for_disp, forcing_types_for_disp[app_ctx->forcing_choice],
346*2b730f8bSJeremy L Thompson                           app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh", app_ctx->degree + 1, app_ctx->degree + 1,
347*2b730f8bSJeremy L Thompson                           U_g_size[fine_level] / num_comp_u, U_l_size[fine_level] / num_comp_u, num_comp_u,
348*2b730f8bSJeremy L Thompson                           (app_ctx->degree == 1 && app_ctx->multigrid_choice != MULTIGRID_NONE) ? "Algebraic multigrid"
349*2b730f8bSJeremy L Thompson                                                                                                 : multigrid_types_for_disp[app_ctx->multigrid_choice],
350*2b730f8bSJeremy L Thompson                           (app_ctx->degree == 1 || app_ctx->multigrid_choice == MULTIGRID_NONE) ? 0 : num_levels));
351ccaff030SJeremy L Thompson 
352d1d35e2fSjeremylt     if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
353e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
354d1d35e2fSjeremylt         CeedInt level = i ? fine_level : 0;
355*2b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
35608140895SJed Brown                               "    Level %" PetscInt_FMT " (%s):\n"
357990fdeb6SJeremy L Thompson                               "      Number of 1D Basis Nodes (p)     : %" CeedInt_FMT "\n"
35808140895SJed Brown                               "      Global Nodes                     : %" PetscInt_FMT "\n"
35908140895SJed Brown                               "      Owned Nodes                      : %" PetscInt_FMT "\n",
360*2b730f8bSJeremy L Thompson                               level, i ? "fine" : "coarse", app_ctx->level_degrees[level] + 1, U_g_size[level] / num_comp_u,
361*2b730f8bSJeremy L Thompson                               U_l_size[level] / num_comp_u));
362ccaff030SJeremy L Thompson       }
363ccaff030SJeremy L Thompson     }
364ccaff030SJeremy L Thompson   }
365ccaff030SJeremy L Thompson 
366ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
367ccaff030SJeremy L Thompson   // Setup SNES
368ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
369ccaff030SJeremy L Thompson   // Performance logging
370*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup));
371*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_snes_setup));
372ccaff030SJeremy L Thompson 
373ccaff030SJeremy L Thompson   // Create SNES
374*2b730f8bSJeremy L Thompson   PetscCall(SNESCreate(comm, &snes));
375*2b730f8bSJeremy L Thompson   PetscCall(SNESSetDM(snes, level_dms[fine_level]));
376ccaff030SJeremy L Thompson 
377ccaff030SJeremy L Thompson   // -- Jacobian evaluators
378*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &jacob_ctx));
379*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &jacob_mat));
380d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
381ccaff030SJeremy L Thompson     // -- Jacobian context for level
382*2b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &jacob_ctx[level]));
383*2b730f8bSJeremy 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,
384*2b730f8bSJeremy L Thompson                                jacob_ctx[level]));
385ccaff030SJeremy L Thompson 
386ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
387*2b730f8bSJeremy 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]));
388*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_MULT, (void (*)(void))ApplyJacobian_Ceed));
389*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL, (void (*)(void))GetDiag_Ceed));
390*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(jacob_mat[level], vectype));
391ccaff030SJeremy L Thompson   }
392e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
393ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
394*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &form_jacob_ctx));
395d1d35e2fSjeremylt   form_jacob_ctx->jacob_ctx  = jacob_ctx;
396d1d35e2fSjeremylt   form_jacob_ctx->num_levels = num_levels;
397d1d35e2fSjeremylt   form_jacob_ctx->jacob_mat  = jacob_mat;
398ccaff030SJeremy L Thompson 
399ccaff030SJeremy L Thompson   // -- Residual evaluation function
400*2b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &res_ctx));
401*2b730f8bSJeremy L Thompson   PetscCall(PetscMemcpy(res_ctx, jacob_ctx[fine_level], sizeof(*jacob_ctx[fine_level])));
4025754ecacSJeremy L Thompson   res_ctx->op = ceed_data[fine_level]->op_residual;
4035754ecacSJeremy L Thompson   res_ctx->qf = ceed_data[fine_level]->qf_residual;
404*2b730f8bSJeremy L Thompson   if (app_ctx->bc_traction_count > 0) res_ctx->neumann_bcs = neumann_bcs;
405*2b730f8bSJeremy L Thompson   else res_ctx->neumann_bcs = NULL;
406*2b730f8bSJeremy L Thompson   PetscCall(SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx));
407ccaff030SJeremy L Thompson 
408ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
409*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &prolong_restr_ctx));
410*2b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(num_levels, &prolong_restr_mat));
411d1d35e2fSjeremylt   for (PetscInt level = 1; level < num_levels; level++) {
412ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
413*2b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &prolong_restr_ctx[level]));
414*2b730f8bSJeremy L Thompson     PetscCall(SetupProlongRestrictCtx(comm, app_ctx, level_dms[level - 1], level_dms[level], U_g[level], U_loc[level - 1], U_loc[level],
415*2b730f8bSJeremy L Thompson                                       ceed_data[level - 1], ceed_data[level], ceed, prolong_restr_ctx[level]));
416ccaff030SJeremy L Thompson 
417ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
418*2b730f8bSJeremy 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],
419*2b730f8bSJeremy L Thompson                              &prolong_restr_mat[level]));
420ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
421*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT, (void (*)(void))Prolong_Ceed));
422*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE, (void (*)(void))Restrict_Ceed));
423*2b730f8bSJeremy L Thompson     PetscCall(MatShellSetVecType(prolong_restr_mat[level], vectype));
424ccaff030SJeremy L Thompson   }
425ccaff030SJeremy L Thompson 
426ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
42717f0843fSJeremy L Thompson   // Setup for AMG coarse solve
428ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
429d1d35e2fSjeremylt   if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
430e3e3df41Sjeremylt     // -- Jacobian Matrix
431*2b730f8bSJeremy L Thompson     PetscCall(DMCreateMatrix(level_dms[0], &jacob_mat_coarse));
432e3e3df41Sjeremylt 
433d1d35e2fSjeremylt     if (app_ctx->degree > 1) {
43417f0843fSJeremy L Thompson       // -- Assemble sparsity pattern
4353047f789SJeremy L Thompson       PetscCount num_entries;
4363047f789SJeremy L Thompson       CeedInt   *rows, *cols;
43717f0843fSJeremy L Thompson       CeedVector coo_values;
438*2b730f8bSJeremy L Thompson       CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries, &rows, &cols);
43917f0843fSJeremy L Thompson       ISLocalToGlobalMapping ltog_row, ltog_col;
440*2b730f8bSJeremy L Thompson       PetscCall(MatGetLocalToGlobalMapping(jacob_mat_coarse, &ltog_row, &ltog_col));
441*2b730f8bSJeremy L Thompson       PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows));
442*2b730f8bSJeremy L Thompson       PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols));
443*2b730f8bSJeremy L Thompson       PetscCall(MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows, cols));
44417f0843fSJeremy L Thompson       free(rows);
44517f0843fSJeremy L Thompson       free(cols);
44617f0843fSJeremy L Thompson       CeedVectorCreate(ceed, num_entries, &coo_values);
447ccaff030SJeremy L Thompson 
448d1d35e2fSjeremylt       // -- Update form_jacob_ctx
44917f0843fSJeremy L Thompson       form_jacob_ctx->coo_values       = coo_values;
45017f0843fSJeremy L Thompson       form_jacob_ctx->op_coarse        = ceed_data[0]->op_jacobian;
451d1d35e2fSjeremylt       form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse;
452e3e3df41Sjeremylt     }
453e3e3df41Sjeremylt   }
454e3e3df41Sjeremylt 
455e3e3df41Sjeremylt   // Set Jacobian function
456d1d35e2fSjeremylt   if (app_ctx->degree > 1) {
457*2b730f8bSJeremy L Thompson     PetscCall(SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level], FormJacobian, form_jacob_ctx));
458e3e3df41Sjeremylt   } else {
459*2b730f8bSJeremy L Thompson     PetscCall(SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse, SNESComputeJacobianDefaultColor, NULL));
460e3e3df41Sjeremylt   }
461ccaff030SJeremy L Thompson 
462ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
463ccaff030SJeremy L Thompson   // Setup KSP
464ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
465ccaff030SJeremy L Thompson   {
466ccaff030SJeremy L Thompson     PC  pc;
467ccaff030SJeremy L Thompson     KSP ksp;
468ccaff030SJeremy L Thompson 
469ccaff030SJeremy L Thompson     // -- KSP
470*2b730f8bSJeremy L Thompson     PetscCall(SNESGetKSP(snes, &ksp));
471*2b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
472*2b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
473*2b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
474*2b730f8bSJeremy L Thompson     PetscCall(KSPSetOptionsPrefix(ksp, "outer_"));
475ccaff030SJeremy L Thompson 
476ccaff030SJeremy L Thompson     // -- Preconditioning
477*2b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
478*2b730f8bSJeremy L Thompson     PetscCall(PCSetDM(pc, level_dms[fine_level]));
479*2b730f8bSJeremy L Thompson     PetscCall(PCSetOptionsPrefix(pc, "outer_"));
480ccaff030SJeremy L Thompson 
481d1d35e2fSjeremylt     if (app_ctx->multigrid_choice == MULTIGRID_NONE) {
482ccaff030SJeremy L Thompson       // ---- No Multigrid
483*2b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
484*2b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
485d1d35e2fSjeremylt     } else if (app_ctx->degree == 1) {
486f892e6d1Sjeremylt       // ---- AMG for degree 1
487*2b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCGAMG));
488ccaff030SJeremy L Thompson     } else {
489ccaff030SJeremy L Thompson       // ---- PCMG
490*2b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCMG));
491ccaff030SJeremy L Thompson 
492ccaff030SJeremy L Thompson       // ------ PCMG levels
493*2b730f8bSJeremy L Thompson       PetscCall(PCMGSetLevels(pc, num_levels, NULL));
494d1d35e2fSjeremylt       for (PetscInt level = 0; level < num_levels; level++) {
495ccaff030SJeremy L Thompson         // -------- Smoother
496d1d35e2fSjeremylt         KSP ksp_smoother, ksp_est;
497d1d35e2fSjeremylt         PC  pc_smoother;
498ccaff030SJeremy L Thompson 
499ccaff030SJeremy L Thompson         // ---------- Smoother KSP
500*2b730f8bSJeremy L Thompson         PetscCall(PCMGGetSmoother(pc, level, &ksp_smoother));
501*2b730f8bSJeremy L Thompson         PetscCall(KSPSetDM(ksp_smoother, level_dms[level]));
502*2b730f8bSJeremy L Thompson         PetscCall(KSPSetDMActive(ksp_smoother, PETSC_FALSE));
503ccaff030SJeremy L Thompson 
504ccaff030SJeremy L Thompson         // ---------- Chebyshev options
505*2b730f8bSJeremy L Thompson         PetscCall(KSPSetType(ksp_smoother, KSPCHEBYSHEV));
506*2b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1));
507*2b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est));
508*2b730f8bSJeremy L Thompson         PetscCall(KSPSetType(ksp_est, KSPCG));
509*2b730f8bSJeremy L Thompson         PetscCall(KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE));
510*2b730f8bSJeremy L Thompson         PetscCall(KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]));
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
513*2b730f8bSJeremy L Thompson         PetscCall(KSPGetPC(ksp_smoother, &pc_smoother));
514*2b730f8bSJeremy L Thompson         PetscCall(PCSetType(pc_smoother, PCJACOBI));
515*2b730f8bSJeremy L Thompson         PetscCall(PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL));
516ccaff030SJeremy L Thompson 
517ccaff030SJeremy L Thompson         // -------- Work vector
518d1d35e2fSjeremylt         if (level != fine_level) {
519*2b730f8bSJeremy L Thompson           PetscCall(PCMGSetX(pc, level, U_g[level]));
520ccaff030SJeremy L Thompson         }
521ccaff030SJeremy L Thompson 
522ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
523ccaff030SJeremy L Thompson         if (level > 0) {
524*2b730f8bSJeremy L Thompson           PetscCall(PCMGSetInterpolation(pc, level, prolong_restr_mat[level]));
525*2b730f8bSJeremy L Thompson           PetscCall(PCMGSetRestriction(pc, level, prolong_restr_mat[level]));
526ccaff030SJeremy L Thompson         }
527ccaff030SJeremy L Thompson       }
528ccaff030SJeremy L Thompson 
529ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
530d1d35e2fSjeremylt       KSP ksp_coarse;
531d1d35e2fSjeremylt       PC  pc_coarse;
532ccaff030SJeremy L Thompson 
533ccaff030SJeremy L Thompson       // -------- Coarse KSP
534*2b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse));
535*2b730f8bSJeremy L Thompson       PetscCall(KSPSetType(ksp_coarse, KSPPREONLY));
536*2b730f8bSJeremy L Thompson       PetscCall(KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse));
537*2b730f8bSJeremy L Thompson       PetscCall(KSPSetOptionsPrefix(ksp_coarse, "coarse_"));
538ccaff030SJeremy L Thompson 
539ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
540*2b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(ksp_coarse, &pc_coarse));
541*2b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc_coarse, PCGAMG));
542*2b730f8bSJeremy L Thompson       PetscCall(PCSetOptionsPrefix(pc_coarse, "coarse_"));
543ccaff030SJeremy L Thompson 
544*2b730f8bSJeremy L Thompson       PetscCall(KSPSetFromOptions(ksp_coarse));
545*2b730f8bSJeremy L Thompson       PetscCall(PCSetFromOptions(pc_coarse));
546ccaff030SJeremy L Thompson 
547ccaff030SJeremy L Thompson       // ------ PCMG options
548*2b730f8bSJeremy L Thompson       PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE));
549*2b730f8bSJeremy L Thompson       PetscCall(PCMGSetNumberSmooth(pc, 3));
550*2b730f8bSJeremy L Thompson       PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type));
551ccaff030SJeremy L Thompson     }
552*2b730f8bSJeremy L Thompson     PetscCall(KSPSetFromOptions(ksp));
553*2b730f8bSJeremy L Thompson     PetscCall(PCSetFromOptions(pc));
554ccaff030SJeremy L Thompson   }
555038d0cd7Sjeremylt   {
556038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
557d1d35e2fSjeremylt     SNESLineSearch line_search;
558481a4840SJed Brown 
559*2b730f8bSJeremy L Thompson     PetscCall(SNESGetLineSearch(snes, &line_search));
560*2b730f8bSJeremy L Thompson     PetscCall(SNESLineSearchSetType(line_search, SNESLINESEARCHCP));
561481a4840SJed Brown   }
562481a4840SJed Brown 
563*2b730f8bSJeremy L Thompson   PetscCall(SNESSetFromOptions(snes));
564ccaff030SJeremy L Thompson 
565ccaff030SJeremy L Thompson   // Performance logging
566*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
567ccaff030SJeremy L Thompson 
568ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
569ccaff030SJeremy L Thompson   // Set initial guess
570ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
571*2b730f8bSJeremy L Thompson   PetscCall(PetscObjectSetName((PetscObject)U, ""));
572*2b730f8bSJeremy L Thompson   PetscCall(VecSet(U, 0.0));
573ccaff030SJeremy L Thompson 
574ccaff030SJeremy L Thompson   // View solution
575d1d35e2fSjeremylt   if (app_ctx->view_soln) {
576*2b730f8bSJeremy L Thompson     PetscCall(ViewSolution(comm, app_ctx, U, 0, 0.0));
577ccaff030SJeremy L Thompson   }
578ccaff030SJeremy L Thompson 
579ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
580ccaff030SJeremy L Thompson   // Solve SNES
581ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
582d1d35e2fSjeremylt   PetscBool snes_monitor = PETSC_FALSE;
583*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor));
5845f24f028Sjeremylt 
585ccaff030SJeremy L Thompson   // Performance logging
586*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve));
587*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(stage_snes_solve));
588ccaff030SJeremy L Thompson 
589ccaff030SJeremy L Thompson   // Timing
590*2b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)snes));
591d1d35e2fSjeremylt   start_time = MPI_Wtime();
592ccaff030SJeremy L Thompson 
593ccaff030SJeremy L Thompson   // Solve for each load increment
5945f24f028Sjeremylt   PetscInt increment;
595d1d35e2fSjeremylt   for (increment = 1; increment <= app_ctx->num_increments; increment++) {
5965f24f028Sjeremylt     // -- Log increment count
597d1d35e2fSjeremylt     if (snes_monitor) {
598*2b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "%" PetscInt_FMT " Load Increment\n", increment - 1));
5995f24f028Sjeremylt     }
6005f24f028Sjeremylt 
601ccaff030SJeremy L Thompson     // -- Scale the problem
602d1d35e2fSjeremylt     PetscScalar load_increment = 1.0 * increment / app_ctx->num_increments,
603*2b730f8bSJeremy L Thompson                 scalingFactor  = load_increment / (increment == 1 ? 1 : res_ctx->load_increment);
604d1d35e2fSjeremylt     res_ctx->load_increment    = load_increment;
605d1d35e2fSjeremylt     if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) {
606*2b730f8bSJeremy L Thompson       PetscCall(VecScale(F, scalingFactor));
607ccaff030SJeremy L Thompson     }
608ccaff030SJeremy L Thompson 
609ccaff030SJeremy L Thompson     // -- Solve
610*2b730f8bSJeremy L Thompson     PetscCall(SNESSolve(snes, F, U));
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson     // -- View solution
613d1d35e2fSjeremylt     if (app_ctx->view_soln) {
614*2b730f8bSJeremy L Thompson       PetscCall(ViewSolution(comm, app_ctx, U, increment, load_increment));
615ccaff030SJeremy L Thompson     }
616ccaff030SJeremy L Thompson 
617ccaff030SJeremy L Thompson     // -- Update SNES iteration count
618ccaff030SJeremy L Thompson     PetscInt its;
619*2b730f8bSJeremy L Thompson     PetscCall(SNESGetIterationNumber(snes, &its));
620bf0f51feSjeremylt     snes_its += its;
621*2b730f8bSJeremy L Thompson     PetscCall(SNESGetLinearSolveIterations(snes, &its));
622bf0f51feSjeremylt     ksp_its += its;
6233951731eSjeremylt 
6243951731eSjeremylt     // -- Check for divergence
6253951731eSjeremylt     SNESConvergedReason reason;
626*2b730f8bSJeremy L Thompson     PetscCall(SNESGetConvergedReason(snes, &reason));
627*2b730f8bSJeremy L Thompson     if (reason < 0) break;
6288355605fSKaren (Ren) Stengel     if (app_ctx->energy_viewer) {
6298355605fSKaren (Ren) Stengel       // -- Log strain energy for current load increment
6308355605fSKaren (Ren) Stengel       CeedScalar energy;
631*2b730f8bSJeremy L Thompson       PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy));
6328355605fSKaren (Ren) Stengel 
6338355605fSKaren (Ren) Stengel       if (!app_ctx->test_mode) {
6348355605fSKaren (Ren) Stengel         // -- Output
635*2b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm, "    Strain Energy                      : %.12e\n", energy));
6368355605fSKaren (Ren) Stengel       }
637*2b730f8bSJeremy L Thompson       PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment, energy));
6388355605fSKaren (Ren) Stengel     }
639ccaff030SJeremy L Thompson   }
640ccaff030SJeremy L Thompson 
641ccaff030SJeremy L Thompson   // Timing
642d1d35e2fSjeremylt   elapsed_time = MPI_Wtime() - start_time;
643ccaff030SJeremy L Thompson 
644ccaff030SJeremy L Thompson   // Performance logging
645*2b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
646ccaff030SJeremy L Thompson 
647ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
648ccaff030SJeremy L Thompson   // Output summary
649ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
650d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
651ccaff030SJeremy L Thompson     // -- SNES
652d1d35e2fSjeremylt     SNESType            snes_type;
653ccaff030SJeremy L Thompson     SNESConvergedReason reason;
654ccaff030SJeremy L Thompson     PetscReal           rnorm;
655*2b730f8bSJeremy L Thompson     PetscCall(SNESGetType(snes, &snes_type));
656*2b730f8bSJeremy L Thompson     PetscCall(SNESGetConvergedReason(snes, &reason));
657*2b730f8bSJeremy L Thompson     PetscCall(SNESGetFunctionNorm(snes, &rnorm));
658*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
659ccaff030SJeremy L Thompson                           "  SNES:\n"
660ccaff030SJeremy L Thompson                           "    SNES Type                          : %s\n"
661ccaff030SJeremy L Thompson                           "    SNES Convergence                   : %s\n"
662990fdeb6SJeremy L Thompson                           "    Number of Load Increments          : %" PetscInt_FMT "\n"
663990fdeb6SJeremy L Thompson                           "    Completed Load Increments          : %" PetscInt_FMT "\n"
66408140895SJed Brown                           "    Total SNES Iterations              : %" PetscInt_FMT "\n"
665ccaff030SJeremy L Thompson                           "    Final rnorm                        : %e\n",
666*2b730f8bSJeremy L Thompson                           snes_type, SNESConvergedReasons[reason], app_ctx->num_increments, increment - 1, snes_its, (double)rnorm));
667ccaff030SJeremy L Thompson 
668ccaff030SJeremy L Thompson     // -- KSP
669ccaff030SJeremy L Thompson     KSP     ksp;
670d1d35e2fSjeremylt     KSPType ksp_type;
671*2b730f8bSJeremy L Thompson     PetscCall(SNESGetKSP(snes, &ksp));
672*2b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
673*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
674ccaff030SJeremy L Thompson                           "  Linear Solver:\n"
6757418ca88SJeremy L Thompson                           "    KSP Type                           : %s\n"
67608140895SJed Brown                           "    Total KSP Iterations               : %" PetscInt_FMT "\n",
677*2b730f8bSJeremy L Thompson                           ksp_type, ksp_its));
678ccaff030SJeremy L Thompson 
679ccaff030SJeremy L Thompson     // -- PC
680ccaff030SJeremy L Thompson     PC     pc;
681d1d35e2fSjeremylt     PCType pc_type;
682*2b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
683*2b730f8bSJeremy L Thompson     PetscCall(PCGetType(pc, &pc_type));
684*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "    PC Type                            : %s\n", pc_type));
685e3e3df41Sjeremylt 
686d1d35e2fSjeremylt     if (!strcmp(pc_type, PCMG)) {
687d1d35e2fSjeremylt       PCMGType pcmg_type;
688*2b730f8bSJeremy L Thompson       PetscCall(PCMGGetType(pc, &pcmg_type));
689*2b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
690ccaff030SJeremy L Thompson                             "  P-Multigrid:\n"
691ccaff030SJeremy L Thompson                             "    PCMG Type                          : %s\n"
692ccaff030SJeremy L Thompson                             "    PCMG Cycle Type                    : %s\n",
693*2b730f8bSJeremy L Thompson                             PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type]));
694ccaff030SJeremy L Thompson 
695ccaff030SJeremy L Thompson       // -- Coarse Solve
696d1d35e2fSjeremylt       KSP    ksp_coarse;
697d1d35e2fSjeremylt       PC     pc_coarse;
698d1d35e2fSjeremylt       PCType pc_type;
699ccaff030SJeremy L Thompson 
700*2b730f8bSJeremy L Thompson       PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse));
701*2b730f8bSJeremy L Thompson       PetscCall(KSPGetType(ksp_coarse, &ksp_type));
702*2b730f8bSJeremy L Thompson       PetscCall(KSPGetPC(ksp_coarse, &pc_coarse));
703*2b730f8bSJeremy L Thompson       PetscCall(PCGetType(pc_coarse, &pc_type));
704*2b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
705ccaff030SJeremy L Thompson                             "    Coarse Solve:\n"
706ccaff030SJeremy L Thompson                             "      KSP Type                         : %s\n"
707ccaff030SJeremy L Thompson                             "      PC Type                          : %s\n",
708*2b730f8bSJeremy L Thompson                             ksp_type, pc_type));
709ccaff030SJeremy L Thompson     }
710ccaff030SJeremy L Thompson   }
711ccaff030SJeremy L Thompson 
712ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
713ccaff030SJeremy L Thompson   // Compute solve time
714ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
715d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
716*2b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm));
717*2b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm));
718*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
719ccaff030SJeremy L Thompson                           "  Performance:\n"
7207418ca88SJeremy L Thompson                           "    SNES Solve Time                    : %g (%g) sec\n"
7217418ca88SJeremy L Thompson                           "    DoFs/Sec in SNES                   : %g (%g) million\n",
722*2b730f8bSJeremy 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));
723ccaff030SJeremy L Thompson   }
724ccaff030SJeremy L Thompson 
725ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
726ccaff030SJeremy L Thompson   // Compute error
727ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
728d1d35e2fSjeremylt   if (app_ctx->forcing_choice == FORCE_MMS) {
729d1d35e2fSjeremylt     CeedScalar        l2_error = 1., l2_U_norm = 1.;
730d1d35e2fSjeremylt     const CeedScalar *true_array;
731d1d35e2fSjeremylt     Vec               error_vec, true_vec;
732ccaff030SJeremy L Thompson 
733ccaff030SJeremy L Thompson     // -- Work vectors
734*2b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &error_vec));
735*2b730f8bSJeremy L Thompson     PetscCall(VecSet(error_vec, 0.0));
736*2b730f8bSJeremy L Thompson     PetscCall(VecDuplicate(U, &true_vec));
737*2b730f8bSJeremy L Thompson     PetscCall(VecSet(true_vec, 0.0));
738ccaff030SJeremy L Thompson 
739ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
740*2b730f8bSJeremy L Thompson     CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln, CEED_MEM_HOST, &true_array);
741*2b730f8bSJeremy L Thompson     PetscCall(VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array));
742*2b730f8bSJeremy L Thompson     PetscCall(DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec));
743*2b730f8bSJeremy L Thompson     PetscCall(VecResetArray(res_ctx->Y_loc));
744d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array);
745ccaff030SJeremy L Thompson 
746ccaff030SJeremy L Thompson     // -- Compute L2 error
747*2b730f8bSJeremy L Thompson     PetscCall(VecWAXPY(error_vec, -1.0, U, true_vec));
748*2b730f8bSJeremy L Thompson     PetscCall(VecNorm(error_vec, NORM_2, &l2_error));
749*2b730f8bSJeremy L Thompson     PetscCall(VecNorm(U, NORM_2, &l2_U_norm));
750d1d35e2fSjeremylt     l2_error /= l2_U_norm;
751ccaff030SJeremy L Thompson 
752ccaff030SJeremy L Thompson     // -- Output
753d1d35e2fSjeremylt     if (!app_ctx->test_mode || l2_error > 0.05) {
754*2b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "    L2 Error                           : %e\n", l2_error));
755ccaff030SJeremy L Thompson     }
756ccaff030SJeremy L Thompson 
757ccaff030SJeremy L Thompson     // -- Cleanup
758*2b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&error_vec));
759*2b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&true_vec));
7602d93065eSjeremylt   }
7612d93065eSjeremylt 
7622d93065eSjeremylt   // ---------------------------------------------------------------------------
7632d93065eSjeremylt   // Compute energy
7642d93065eSjeremylt   // ---------------------------------------------------------------------------
7658355605fSKaren (Ren) Stengel   PetscReal energy;
766*2b730f8bSJeremy L Thompson   PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy));
7678355605fSKaren (Ren) Stengel   if (!app_ctx->test_mode) {
7682d93065eSjeremylt     // -- Output
769*2b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "    Strain Energy                      : %.12e\n", energy));
770ccaff030SJeremy L Thompson   }
771*2b730f8bSJeremy L Thompson   PetscCall(RegressionTests_solids(app_ctx, energy));
772ccaff030SJeremy L Thompson 
773ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
7745c25879aSJeremy L Thompson   // Output diagnostic quantities
7755c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
776d1d35e2fSjeremylt   if (app_ctx->view_soln || app_ctx->view_final_soln) {
7775c25879aSJeremy L Thompson     // -- Setup context
778d1d35e2fSjeremylt     UserMult diagnostic_ctx;
779*2b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(1, &diagnostic_ctx));
780*2b730f8bSJeremy L Thompson     PetscCall(PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)));
781d1d35e2fSjeremylt     diagnostic_ctx->dm = dm_diagnostic;
782d1d35e2fSjeremylt     diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic;
7835c25879aSJeremy L Thompson 
7845c25879aSJeremy L Thompson     // -- Compute and output
785*2b730f8bSJeremy L Thompson     PetscCall(ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx, app_ctx, U, ceed_data[fine_level]->elem_restr_diagnostic));
7865c25879aSJeremy L Thompson 
7875c25879aSJeremy L Thompson     // -- Cleanup
788*2b730f8bSJeremy L Thompson     PetscCall(PetscFree(diagnostic_ctx));
7895c25879aSJeremy L Thompson   }
7905c25879aSJeremy L Thompson 
7915c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
792ccaff030SJeremy L Thompson   // Free objects
793ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
794ccaff030SJeremy L Thompson   // Data in arrays per level
795d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
796ccaff030SJeremy L Thompson     // Vectors
797*2b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&U_g[level]));
798*2b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&U_loc[level]));
799ccaff030SJeremy L Thompson 
800ccaff030SJeremy L Thompson     // Jacobian matrix and data
801*2b730f8bSJeremy L Thompson     PetscCall(VecDestroy(&jacob_ctx[level]->Y_loc));
802*2b730f8bSJeremy L Thompson     PetscCall(MatDestroy(&jacob_mat[level]));
803*2b730f8bSJeremy L Thompson     PetscCall(PetscFree(jacob_ctx[level]));
804ccaff030SJeremy L Thompson 
805ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
806ccaff030SJeremy L Thompson     if (level > 0) {
807*2b730f8bSJeremy L Thompson       PetscCall(PetscFree(prolong_restr_ctx[level]));
808*2b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&prolong_restr_mat[level]));
809ccaff030SJeremy L Thompson     }
810ccaff030SJeremy L Thompson 
811ccaff030SJeremy L Thompson     // DM
812*2b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&level_dms[level]));
813ccaff030SJeremy L Thompson 
814ccaff030SJeremy L Thompson     // libCEED objects
815*2b730f8bSJeremy L Thompson     PetscCall(CeedDataDestroy(level, ceed_data[level]));
816ccaff030SJeremy L Thompson   }
817ccaff030SJeremy L Thompson 
818*2b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&app_ctx->energy_viewer));
8198355605fSKaren (Ren) Stengel 
820ccaff030SJeremy L Thompson   // Arrays
821*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_g));
822*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_loc));
823*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_g_size));
824*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_l_size));
825*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(U_loc_size));
826*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_ctx));
827*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_mat));
828*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(prolong_restr_ctx));
829*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(prolong_restr_mat));
830*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(app_ctx->level_degrees));
831*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(ceed_data));
832ccaff030SJeremy L Thompson 
833ccaff030SJeremy L Thompson   // libCEED objects
83417f0843fSJeremy L Thompson   CeedVectorDestroy(&form_jacob_ctx->coo_values);
835d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys);
836d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys_smoother);
837ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
838ccaff030SJeremy L Thompson 
839ccaff030SJeremy L Thompson   // PETSc objects
840*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&U));
841*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&R));
842*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&R_loc));
843*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&F));
844*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&F_loc));
845*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&neumann_bcs));
846*2b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&bcs_loc));
847*2b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&jacob_mat_coarse));
848*2b730f8bSJeremy L Thompson   PetscCall(SNESDestroy(&snes));
849*2b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_orig));
850*2b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_energy));
851*2b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm_diagnostic));
852*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(level_dms));
853ccaff030SJeremy L Thompson 
8545754ecacSJeremy L Thompson   // -- Function list
855*2b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupPhysics));
856*2b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupSmootherPhysics));
857*2b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel));
858*2b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedLevel));
8595754ecacSJeremy L Thompson 
860ccaff030SJeremy L Thompson   // Structs
861*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(res_ctx));
862*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(form_jacob_ctx));
863*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(jacob_coarse_ctx));
864*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(app_ctx));
865*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(problem_functions));
866*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(units));
867ccaff030SJeremy L Thompson 
868ccaff030SJeremy L Thompson   return PetscFinalize();
869ccaff030SJeremy L Thompson }
870