xref: /libCEED/examples/solids/elasticity.c (revision 8355605fb12655f53893674f518ac7968aa830ae)
1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4ccaff030SJeremy L Thompson //
5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Elasticity
18ccaff030SJeremy L Thompson //
19ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve
20ccaff030SJeremy L Thompson //   elasticity problems.
21ccaff030SJeremy L Thompson //
22ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex.
23ccaff030SJeremy L Thompson //
24ccaff030SJeremy L Thompson // Build with:
25ccaff030SJeremy L Thompson //
26ccaff030SJeremy L Thompson //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27ccaff030SJeremy L Thompson //
28ccaff030SJeremy L Thompson // Sample runs:
29ccaff030SJeremy L Thompson //
30eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms
31eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_clamp 998,999 -bc_clamp_998_translate 0.1,0.2,0.3 -problem SS-NH -forcing none -ceed /cpu/self
32eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_clamp 998,999 -bc_clamp_998_rotate 1,0,0,0.2 -problem FSInitial-NH1 -forcing none -ceed /gpu/cuda
33ccaff030SJeremy L Thompson //
34ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
35ccaff030SJeremy L Thompson //
36*8355605fSKaren (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
37*8355605fSKaren (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
38*8355605fSKaren (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
39ccaff030SJeremy L Thompson 
40ccaff030SJeremy L Thompson /// @file
41ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
44ccaff030SJeremy L Thompson 
45ccaff030SJeremy L Thompson #include "elasticity.h"
46ccaff030SJeremy L Thompson 
47ccaff030SJeremy L Thompson int main(int argc, char **argv) {
48ccaff030SJeremy L Thompson   PetscInt       ierr;
49ccaff030SJeremy L Thompson   MPI_Comm       comm;
50ccaff030SJeremy L Thompson   // Context structs
51d1d35e2fSjeremylt   AppCtx         app_ctx;                  // Contains problem options
52*8355605fSKaren (Ren) Stengel   Physics        phys = NULL;              // Physical constants for Neo-Hookean
53*8355605fSKaren (Ren) Stengel   Physics_MR     phys_MR = NULL;           // Physical constants for Mooney-Rivlin
54d1d35e2fSjeremylt   Physics        phys_smoother = NULL;     // Separate context if nu_smoother set
55ccaff030SJeremy L Thompson   Units          units;                    // Contains units scaling
56ccaff030SJeremy L Thompson   // PETSc objects
57d1d35e2fSjeremylt   PetscLogStage  stage_dm_setup, stage_libceed_setup,
58d1d35e2fSjeremylt                  stage_snes_setup, stage_snes_solve;
59d1d35e2fSjeremylt   DM             dm_orig;                  // Distributed DM to clone
60d1d35e2fSjeremylt   DM             dm_energy, dm_diagnostic; // DMs for postprocessing
61d1d35e2fSjeremylt   DM             *level_dms;
62d1d35e2fSjeremylt   Vec            U, *U_g, *U_loc;          // U: solution, R: residual, F: forcing
63d1d35e2fSjeremylt   Vec            R, R_loc, F, F_loc;       // g: global, loc: local
64d1d35e2fSjeremylt   Vec            neumann_bcs = NULL, bcs_loc = NULL;
65d1d35e2fSjeremylt   SNES           snes, snes_coarse = NULL;
66d1d35e2fSjeremylt   Mat            *jacob_mat, jacob_mat_coarse, *prolong_restr_mat;
67ccaff030SJeremy L Thompson   // PETSc data
68d1d35e2fSjeremylt   UserMult       res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx;
69d1d35e2fSjeremylt   FormJacobCtx   form_jacob_ctx;
70d1d35e2fSjeremylt   UserMultProlongRestr *prolong_restr_ctx;
71d1d35e2fSjeremylt   PCMGCycleType  pcmg_cycle_type = PC_MG_CYCLE_V;
72ccaff030SJeremy L Thompson   // libCEED objects
73d99fa3c5SJeremy L Thompson   Ceed           ceed;
74d1d35e2fSjeremylt   CeedData       *ceed_data;
75d1d35e2fSjeremylt   CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL;
76ccaff030SJeremy L Thompson   // Parameters
77d1d35e2fSjeremylt   PetscInt       num_comp_u = 3;                 // 3 DoFs in 3D
78d1d35e2fSjeremylt   PetscInt       num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic
79d1d35e2fSjeremylt   PetscInt       num_levels = 1, fine_level = 0;
80bf0f51feSjeremylt   PetscInt       *U_g_size, *U_l_size, *U_loc_size;
81bf0f51feSjeremylt   PetscInt       snes_its = 0, ksp_its = 0;
82d1d35e2fSjeremylt   double         start_time, elapsed_time, min_time, max_time;
83ccaff030SJeremy L Thompson 
84ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
85bf0f51feSjeremylt   if (ierr) return ierr;
86ccaff030SJeremy L Thompson 
87ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
88ccaff030SJeremy L Thompson   // Process command line options
89ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
90ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
91ccaff030SJeremy L Thompson 
92ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
93d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr);
94d1d35e2fSjeremylt   ierr = ProcessCommandLineOptions(comm, app_ctx); CHKERRQ(ierr);
95d1d35e2fSjeremylt   num_levels = app_ctx->num_levels;
96d1d35e2fSjeremylt   fine_level = num_levels - 1;
97ccaff030SJeremy L Thompson 
98*8355605fSKaren (Ren) Stengel   if (app_ctx->problem_choice != ELAS_FSInitial_MR1) {
99ccaff030SJeremy L Thompson     // -- Set Poison's ratio, Young's Modulus
100ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
101*8355605fSKaren (Ren) Stengel     ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
102*8355605fSKaren (Ren) Stengel     ierr = ProcessPhysics_General(comm, app_ctx, phys, phys_MR, units);
103*8355605fSKaren (Ren) Stengel     CHKERRQ(ierr);
104d1d35e2fSjeremylt     if (fabs(app_ctx->nu_smoother) > 1E-14) {
105d1d35e2fSjeremylt       ierr = PetscMalloc1(1, &phys_smoother); CHKERRQ(ierr);
106d1d35e2fSjeremylt       ierr = PetscMemcpy(phys_smoother, phys, sizeof(*phys)); CHKERRQ(ierr);
107d1d35e2fSjeremylt       phys_smoother->nu = app_ctx->nu_smoother;
108f7b4142eSJeremy L Thompson     }
109*8355605fSKaren (Ren) Stengel   } else {
110*8355605fSKaren (Ren) Stengel     // -- Set Mooney-Rivlin parameters
111*8355605fSKaren (Ren) Stengel     ierr = PetscMalloc1(1, &phys_MR); CHKERRQ(ierr);
112*8355605fSKaren (Ren) Stengel     ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
113*8355605fSKaren (Ren) Stengel     ierr = ProcessPhysics_General(comm, app_ctx, phys, phys_MR, units);
114*8355605fSKaren (Ren) Stengel     CHKERRQ(ierr);
115*8355605fSKaren (Ren) Stengel   }
116ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
117d1d35e2fSjeremylt   // Initialize libCEED
11862e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
119d1d35e2fSjeremylt   // Initialize backend
120d1d35e2fSjeremylt   CeedInit(app_ctx->ceed_resource, &ceed);
12162e9c006SJeremy L Thompson 
12262e9c006SJeremy L Thompson   // Check preferred MemType
123d1d35e2fSjeremylt   CeedMemType mem_type_backend;
124d1d35e2fSjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
12562e9c006SJeremy L Thompson 
126777ff853SJeremy L Thompson   // Wrap context in libCEED objects
127d1d35e2fSjeremylt   CeedQFunctionContextCreate(ceed, &ctx_phys);
128*8355605fSKaren (Ren) Stengel   switch (app_ctx->problem_choice) {
129*8355605fSKaren (Ren) Stengel   case ELAS_LINEAR:
130d1d35e2fSjeremylt     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
131777ff853SJeremy L Thompson                                 sizeof(*phys), phys);
132*8355605fSKaren (Ren) Stengel     break;
133*8355605fSKaren (Ren) Stengel   case ELAS_SS_NH:
134*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
135*8355605fSKaren (Ren) Stengel                                 sizeof(*phys), phys);
136*8355605fSKaren (Ren) Stengel     break;
137*8355605fSKaren (Ren) Stengel   case ELAS_FSInitial_NH1:
138*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
139*8355605fSKaren (Ren) Stengel                                 sizeof(*phys), phys);
140*8355605fSKaren (Ren) Stengel     break;
141*8355605fSKaren (Ren) Stengel   case ELAS_FSInitial_NH2:
142*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
143*8355605fSKaren (Ren) Stengel                                 sizeof(*phys), phys);
144*8355605fSKaren (Ren) Stengel     break;
145*8355605fSKaren (Ren) Stengel   case ELAS_FSCurrent_NH1:
146*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
147*8355605fSKaren (Ren) Stengel                                 sizeof(*phys), phys);
148*8355605fSKaren (Ren) Stengel     break;
149*8355605fSKaren (Ren) Stengel   case ELAS_FSCurrent_NH2:
150*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
151*8355605fSKaren (Ren) Stengel                                 sizeof(*phys), phys);
152*8355605fSKaren (Ren) Stengel     break;
153*8355605fSKaren (Ren) Stengel   case ELAS_FSInitial_MR1:
154*8355605fSKaren (Ren) Stengel     CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER,
155*8355605fSKaren (Ren) Stengel                                 sizeof(*phys_MR), phys_MR);
156*8355605fSKaren (Ren) Stengel     break;
157*8355605fSKaren (Ren) Stengel   }
158*8355605fSKaren (Ren) Stengel 
159d1d35e2fSjeremylt   if (phys_smoother) {
160d1d35e2fSjeremylt     CeedQFunctionContextCreate(ceed, &ctx_phys_smoother);
161d1d35e2fSjeremylt     CeedQFunctionContextSetData(ctx_phys_smoother, CEED_MEM_HOST, CEED_USE_POINTER,
162d1d35e2fSjeremylt                                 sizeof(*phys_smoother), phys_smoother);
163777ff853SJeremy L Thompson   }
164777ff853SJeremy L Thompson 
16562e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
166ccaff030SJeremy L Thompson   // Setup DM
167ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
168ccaff030SJeremy L Thompson   // Performance logging
169d1d35e2fSjeremylt   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup);
170ccaff030SJeremy L Thompson   CHKERRQ(ierr);
171d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_dm_setup); CHKERRQ(ierr);
172ccaff030SJeremy L Thompson 
173ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
174d1d35e2fSjeremylt   ierr = CreateDistributedDM(comm, app_ctx, &dm_orig); CHKERRQ(ierr);
175b68a8d79SJed Brown   VecType vectype;
176d1d35e2fSjeremylt   switch (mem_type_backend) {
177b68a8d79SJed Brown   case CEED_MEM_HOST: vectype = VECSTANDARD; break;
178b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
179b68a8d79SJed Brown     const char *resolved;
180b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
181b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
182b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
183b68a8d79SJed Brown     else vectype = VECSTANDARD;
184b68a8d79SJed Brown   }
185b68a8d79SJed Brown   }
186d1d35e2fSjeremylt   ierr = DMSetVecType(dm_orig, vectype); CHKERRQ(ierr);
187d1d35e2fSjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
188ccaff030SJeremy L Thompson 
189ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
190d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &level_dms); CHKERRQ(ierr);
191d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
192d1d35e2fSjeremylt     ierr = DMClone(dm_orig, &level_dms[level]); CHKERRQ(ierr);
193d1d35e2fSjeremylt     ierr = DMGetVecType(dm_orig, &vectype); CHKERRQ(ierr);
194d1d35e2fSjeremylt     ierr = DMSetVecType(level_dms[level], vectype); CHKERRQ(ierr);
195d1d35e2fSjeremylt     ierr = SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level],
196d1d35e2fSjeremylt                            PETSC_TRUE, num_comp_u); CHKERRQ(ierr);
197a3c02c40SJeremy L Thompson     // -- Label field components for viewing
198a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
199a3c02c40SJeremy L Thompson     PetscSection section;
200d1d35e2fSjeremylt     ierr = DMGetLocalSection(level_dms[level], &section); CHKERRQ(ierr);
201a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr);
202a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
203a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
204a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
205a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
206a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
207a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
208a3c02c40SJeremy L Thompson   }
209a3c02c40SJeremy L Thompson 
210a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
211d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_energy); CHKERRQ(ierr);
212d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level],
213d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_e); CHKERRQ(ierr);
214d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_diagnostic); CHKERRQ(ierr);
215d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_diagnostic, app_ctx,
216d1d35e2fSjeremylt                          app_ctx->level_degrees[fine_level],
217d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_u + num_comp_d); CHKERRQ(ierr);
218d1d35e2fSjeremylt   ierr = DMSetVecType(dm_energy, vectype); CHKERRQ(ierr);
219d1d35e2fSjeremylt   ierr = DMSetVecType(dm_diagnostic, vectype); CHKERRQ(ierr);
220a3c02c40SJeremy L Thompson   {
221a3c02c40SJeremy L Thompson     // -- Label field components for viewing
222a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
223a3c02c40SJeremy L Thompson     PetscSection section;
224d1d35e2fSjeremylt     ierr = DMGetLocalSection(dm_diagnostic, &section); CHKERRQ(ierr);
225a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr);
2268ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
2278ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
2288ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
2298ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
2308ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
2318ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
232379387d4SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure");
2338ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
2347ab18febSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain");
23513fdad4bSJeremy L Thompson     CHKERRQ(ierr);
23613fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2");
23713fdad4bSJeremy L Thompson     CHKERRQ(ierr);
23813fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 6, "detJ");
23913fdad4bSJeremy L Thompson     CHKERRQ(ierr);
24013fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity");
241a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
242ccaff030SJeremy L Thompson   }
243ccaff030SJeremy L Thompson 
244ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
245ccaff030SJeremy L Thompson   // Setup solution and work vectors
246ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
247ccaff030SJeremy L Thompson   // Allocate arrays
248d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g); CHKERRQ(ierr);
249d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc); CHKERRQ(ierr);
250d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g_size); CHKERRQ(ierr);
251d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_l_size); CHKERRQ(ierr);
252d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc_size); CHKERRQ(ierr);
253ccaff030SJeremy L Thompson 
254ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
255d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
256ccaff030SJeremy L Thompson     // -- Create global unknown vector U
257d1d35e2fSjeremylt     ierr = DMCreateGlobalVector(level_dms[level], &U_g[level]); CHKERRQ(ierr);
258d1d35e2fSjeremylt     ierr = VecGetSize(U_g[level], &U_g_size[level]); CHKERRQ(ierr);
259ccaff030SJeremy L Thompson     // Note: Local size for matShell
260d1d35e2fSjeremylt     ierr = VecGetLocalSize(U_g[level], &U_l_size[level]); CHKERRQ(ierr);
261ccaff030SJeremy L Thompson 
262d1d35e2fSjeremylt     // -- Create local unknown vector U_loc
263d1d35e2fSjeremylt     ierr = DMCreateLocalVector(level_dms[level], &U_loc[level]); CHKERRQ(ierr);
264ccaff030SJeremy L Thompson     // Note: local size for libCEED
265d1d35e2fSjeremylt     ierr = VecGetSize(U_loc[level], &U_loc_size[level]); CHKERRQ(ierr);
266ccaff030SJeremy L Thompson   }
267ccaff030SJeremy L Thompson 
268ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
269d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &U); CHKERRQ(ierr);
270d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &R); CHKERRQ(ierr);
271d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &F); CHKERRQ(ierr);
272d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &R_loc); CHKERRQ(ierr);
273d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &F_loc); CHKERRQ(ierr);
274ccaff030SJeremy L Thompson 
275ccaff030SJeremy L Thompson   // Performance logging
276ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
277ccaff030SJeremy L Thompson 
278ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
279ccaff030SJeremy L Thompson   // Set up libCEED
280ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
281ccaff030SJeremy L Thompson   // Performance logging
282d1d35e2fSjeremylt   ierr = PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup);
283ccaff030SJeremy L Thompson   CHKERRQ(ierr);
284d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_libceed_setup); CHKERRQ(ierr);
285ccaff030SJeremy L Thompson 
286ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
287d1d35e2fSjeremylt   CeedVector force_ceed;
288ccaff030SJeremy L Thompson   CeedScalar *f;
289d1d35e2fSjeremylt   PetscMemType force_mem_type;
290d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
291d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(F_loc, &f, &force_mem_type); CHKERRQ(ierr);
292d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed);
293d1d35e2fSjeremylt     CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f);
294ccaff030SJeremy L Thompson   }
295ccaff030SJeremy L Thompson 
296fe394131Sjeremylt   // -- Create libCEED local Neumann BCs vector
297d1d35e2fSjeremylt   CeedVector neumann_ceed;
298fe394131Sjeremylt   CeedScalar *n;
299d1d35e2fSjeremylt   PetscMemType nummann_mem_type;
300d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
301d1d35e2fSjeremylt     ierr = VecDuplicate(U, &neumann_bcs); CHKERRQ(ierr);
302d1d35e2fSjeremylt     ierr = VecDuplicate(U_loc[fine_level], &bcs_loc); CHKERRQ(ierr);
303d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type); CHKERRQ(ierr);
304d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed);
305d1d35e2fSjeremylt     CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type),
306fe394131Sjeremylt                        CEED_USE_POINTER, n);
307fe394131Sjeremylt   }
308fe394131Sjeremylt 
309ccaff030SJeremy L Thompson   // -- Setup libCEED objects
310d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
311d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
312d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &ceed_data[fine_level]); CHKERRQ(ierr);
313ccaff030SJeremy L Thompson   {
314d1d35e2fSjeremylt     ierr = SetupLibceedFineLevel(level_dms[fine_level], dm_energy, dm_diagnostic,
315d1d35e2fSjeremylt                                  ceed, app_ctx, ctx_phys, ceed_data, fine_level,
316d1d35e2fSjeremylt                                  num_comp_u, U_g_size[fine_level], U_loc_size[fine_level],
317d1d35e2fSjeremylt                                  force_ceed, neumann_ceed);
318ccaff030SJeremy L Thompson     CHKERRQ(ierr);
319ccaff030SJeremy L Thompson   }
320d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
321d1d35e2fSjeremylt   for (PetscInt level = num_levels - 2; level >= 0; level--) {
322d1d35e2fSjeremylt     ierr = PetscCalloc1(1, &ceed_data[level]); CHKERRQ(ierr);
323d99fa3c5SJeremy L Thompson 
324d99fa3c5SJeremy L Thompson     // Get global communication restriction
325d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
326d1d35e2fSjeremylt     ierr = VecSet(U_loc[level+1], 1.0); CHKERRQ(ierr);
327d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[level+1], U_loc[level+1], ADD_VALUES,
328d1d35e2fSjeremylt                            U_g[level+1]); CHKERRQ(ierr);
329d1d35e2fSjeremylt     ierr = DMGlobalToLocal(level_dms[level+1], U_g[level+1], INSERT_VALUES,
330d1d35e2fSjeremylt                            U_loc[level+1]); CHKERRQ(ierr);
331d99fa3c5SJeremy L Thompson 
332d99fa3c5SJeremy L Thompson     // Place in libCEED array
333d99fa3c5SJeremy L Thompson     const PetscScalar *m;
334d1d35e2fSjeremylt     PetscMemType m_mem_type;
335d1d35e2fSjeremylt     ierr = VecGetArrayReadAndMemType(U_loc[level+1], &m, &m_mem_type);
336d1d35e2fSjeremylt     CHKERRQ(ierr);
337d1d35e2fSjeremylt     CeedVectorSetArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
338d99fa3c5SJeremy L Thompson                        CEED_USE_POINTER, (CeedScalar *)m);
339ccaff030SJeremy L Thompson 
3401c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
341d1d35e2fSjeremylt     ierr = SetupLibceedLevel(level_dms[level], ceed, app_ctx,
342d1d35e2fSjeremylt                              ceed_data, level, num_comp_u, U_g_size[level],
343d1d35e2fSjeremylt                              U_loc_size[level], ceed_data[level+1]->x_ceed);
344777ff853SJeremy L Thompson     CHKERRQ(ierr);
345d99fa3c5SJeremy L Thompson 
346d99fa3c5SJeremy L Thompson     // Restore PETSc vector
347d1d35e2fSjeremylt     CeedVectorTakeArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
348d99fa3c5SJeremy L Thompson                         (CeedScalar **)&m);
349d1d35e2fSjeremylt     ierr = VecRestoreArrayReadAndMemType(U_loc[level+1], &m); CHKERRQ(ierr);
350d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
351d1d35e2fSjeremylt     ierr = VecZeroEntries(U_loc[level+1]); CHKERRQ(ierr);
352ccaff030SJeremy L Thompson   }
353ccaff030SJeremy L Thompson 
354ccaff030SJeremy L Thompson   // Performance logging
355ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
356ccaff030SJeremy L Thompson 
357ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
358fe394131Sjeremylt   // Setup global forcing and Neumann BC vectors
359ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
360ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
361ccaff030SJeremy L Thompson 
362d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
363d1d35e2fSjeremylt     CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL);
364d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(F_loc, &f); CHKERRQ(ierr);
365d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F);
366ccaff030SJeremy L Thompson     CHKERRQ(ierr);
367d1d35e2fSjeremylt     CeedVectorDestroy(&force_ceed);
368ccaff030SJeremy L Thompson   }
369ccaff030SJeremy L Thompson 
370d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
371d1d35e2fSjeremylt     ierr = VecZeroEntries(neumann_bcs); CHKERRQ(ierr);
372d1d35e2fSjeremylt     CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL);
373d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(bcs_loc, &n); CHKERRQ(ierr);
374d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs);
375fe394131Sjeremylt     CHKERRQ(ierr);
376d1d35e2fSjeremylt     CeedVectorDestroy(&neumann_ceed);
377fe394131Sjeremylt   }
378fe394131Sjeremylt 
379ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
380ccaff030SJeremy L Thompson   // Print problem summary
381ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
382d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
383ccaff030SJeremy L Thompson     const char *usedresource;
384ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
38577841947SLeila Ghaffari     char hostname[PETSC_MAX_PATH_LEN];
38677841947SLeila Ghaffari     ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
38777841947SLeila Ghaffari     PetscInt comm_size;
38877841947SLeila Ghaffari     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
389ccaff030SJeremy L Thompson 
390ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
39117fd649aSJeremy L Thompson                        "\n-- Elasticity Example - libCEED + PETSc --\n"
39277841947SLeila Ghaffari                        "  MPI:\n"
39377841947SLeila Ghaffari                        "    Hostname                           : %s\n"
39477841947SLeila Ghaffari                        "    Total ranks                        : %d\n"
395ccaff030SJeremy L Thompson                        "  libCEED:\n"
39662e9c006SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
397b68a8d79SJed Brown                        "    libCEED Backend MemType            : %s\n",
39877841947SLeila Ghaffari                        hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend]);
39962e9c006SJeremy L Thompson     CHKERRQ(ierr);
400ccaff030SJeremy L Thompson 
40162e9c006SJeremy L Thompson     VecType vecType;
40262e9c006SJeremy L Thompson     ierr = VecGetType(U, &vecType); CHKERRQ(ierr);
40362e9c006SJeremy L Thompson     ierr = PetscPrintf(comm,
40462e9c006SJeremy L Thompson                        "  PETSc:\n"
40562e9c006SJeremy L Thompson                        "    PETSc Vec Type                     : %s\n",
40662e9c006SJeremy L Thompson                        vecType); CHKERRQ(ierr);
40762e9c006SJeremy L Thompson 
408ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
409ccaff030SJeremy L Thompson                        "  Problem:\n"
410ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
411ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
412ccaff030SJeremy L Thompson                        "  Mesh:\n"
413ccaff030SJeremy L Thompson                        "    File                               : %s\n"
414ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
415ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
416ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
417ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
418ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
419ccaff030SJeremy L Thompson                        "  Multigrid:\n"
420ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
421ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
422d1d35e2fSjeremylt                        problemTypesForDisp[app_ctx->problem_choice],
423d1d35e2fSjeremylt                        forcing_types_for_disp[app_ctx->forcing_choice],
424d1d35e2fSjeremylt                        app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh",
425d1d35e2fSjeremylt                        app_ctx->degree + 1, app_ctx->degree + 1,
426d1d35e2fSjeremylt                        U_g_size[fine_level]/num_comp_u, U_l_size[fine_level]/num_comp_u, num_comp_u,
427d1d35e2fSjeremylt                        (app_ctx->degree == 1 &&
428d1d35e2fSjeremylt                         app_ctx->multigrid_choice != MULTIGRID_NONE) ?
429f892e6d1Sjeremylt                        "Algebraic multigrid" :
430d1d35e2fSjeremylt                        multigrid_types_for_disp[app_ctx->multigrid_choice],
431d1d35e2fSjeremylt                        (app_ctx->degree == 1 ||
432d1d35e2fSjeremylt                         app_ctx->multigrid_choice == MULTIGRID_NONE) ?
433d1d35e2fSjeremylt                        0 : num_levels); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson 
435d1d35e2fSjeremylt     if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
436e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
437d1d35e2fSjeremylt         CeedInt level = i ? fine_level : 0;
438ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
439ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
440ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
441ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
442ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
443ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
444d1d35e2fSjeremylt                            app_ctx->level_degrees[level] + 1,
445d1d35e2fSjeremylt                            U_g_size[level]/num_comp_u, U_l_size[level]/num_comp_u);
446ccaff030SJeremy L Thompson         CHKERRQ(ierr);
447ccaff030SJeremy L Thompson       }
448ccaff030SJeremy L Thompson     }
449ccaff030SJeremy L Thompson   }
450ccaff030SJeremy L Thompson 
451ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
452ccaff030SJeremy L Thompson   // Setup SNES
453ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
454ccaff030SJeremy L Thompson   // Performance logging
455d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup);
456ccaff030SJeremy L Thompson   CHKERRQ(ierr);
457d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_setup); CHKERRQ(ierr);
458ccaff030SJeremy L Thompson 
459ccaff030SJeremy L Thompson   // Create SNES
460ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
461d1d35e2fSjeremylt   ierr = SNESSetDM(snes, level_dms[fine_level]); CHKERRQ(ierr);
462ccaff030SJeremy L Thompson 
463ccaff030SJeremy L Thompson   // -- Jacobian evaluators
464d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_ctx); CHKERRQ(ierr);
465d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_mat); CHKERRQ(ierr);
466d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
467ccaff030SJeremy L Thompson     // -- Jacobian context for level
468d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &jacob_ctx[level]); CHKERRQ(ierr);
469d1d35e2fSjeremylt     ierr = SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level],
470d1d35e2fSjeremylt                             U_loc[level], ceed_data[level], ceed, ctx_phys,
471d1d35e2fSjeremylt                             ctx_phys_smoother, jacob_ctx[level]); CHKERRQ(ierr);
472ccaff030SJeremy L Thompson 
473ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
474d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level],
475d1d35e2fSjeremylt                           U_g_size[level], jacob_ctx[level], &jacob_mat[level]);
476ccaff030SJeremy L Thompson     CHKERRQ(ierr);
477d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_MULT,
478ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
479ccaff030SJeremy L Thompson     CHKERRQ(ierr);
480d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL,
481ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
482d1d35e2fSjeremylt     ierr = MatShellSetVecType(jacob_mat[level], vectype); CHKERRQ(ierr);
483ccaff030SJeremy L Thompson   }
484e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
485ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
486d1d35e2fSjeremylt   ierr = PetscMalloc1(1, &form_jacob_ctx); CHKERRQ(ierr);
487d1d35e2fSjeremylt   form_jacob_ctx->jacob_ctx = jacob_ctx;
488d1d35e2fSjeremylt   form_jacob_ctx->num_levels = num_levels;
489d1d35e2fSjeremylt   form_jacob_ctx->jacob_mat = jacob_mat;
490ccaff030SJeremy L Thompson 
491ccaff030SJeremy L Thompson   // -- Residual evaluation function
492d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &res_ctx); CHKERRQ(ierr);
493d1d35e2fSjeremylt   ierr = PetscMemcpy(res_ctx, jacob_ctx[fine_level],
494d1d35e2fSjeremylt                      sizeof(*jacob_ctx[fine_level])); CHKERRQ(ierr);
495d1d35e2fSjeremylt   res_ctx->op = ceed_data[fine_level]->op_apply;
496d1d35e2fSjeremylt   res_ctx->qf = ceed_data[fine_level]->qf_apply;
497d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0)
498d1d35e2fSjeremylt     res_ctx->neumann_bcs = neumann_bcs;
499fe394131Sjeremylt   else
500d1d35e2fSjeremylt     res_ctx->neumann_bcs = NULL;
501d1d35e2fSjeremylt   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx); CHKERRQ(ierr);
502ccaff030SJeremy L Thompson 
503ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
504d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_ctx); CHKERRQ(ierr);
505d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_mat); CHKERRQ(ierr);
506d1d35e2fSjeremylt   for (PetscInt level = 1; level < num_levels; level++) {
507ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
508d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &prolong_restr_ctx[level]); CHKERRQ(ierr);
509d1d35e2fSjeremylt     ierr = SetupProlongRestrictCtx(comm, app_ctx, level_dms[level-1],
510d1d35e2fSjeremylt                                    level_dms[level], U_g[level], U_loc[level-1],
511d1d35e2fSjeremylt                                    U_loc[level], ceed_data[level-1],
512d1d35e2fSjeremylt                                    ceed_data[level], ceed,
513d1d35e2fSjeremylt                                    prolong_restr_ctx[level]); CHKERRQ(ierr);
514ccaff030SJeremy L Thompson 
515ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
516d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level-1], U_g_size[level],
517d1d35e2fSjeremylt                           U_g_size[level-1], prolong_restr_ctx[level],
518d1d35e2fSjeremylt                           &prolong_restr_mat[level]); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
520d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT,
521ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
522d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE,
523ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
524ccaff030SJeremy L Thompson     CHKERRQ(ierr);
525d1d35e2fSjeremylt     ierr = MatShellSetVecType(prolong_restr_mat[level], vectype); CHKERRQ(ierr);
526ccaff030SJeremy L Thompson   }
527ccaff030SJeremy L Thompson 
528ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
529ccaff030SJeremy L Thompson   // Setup dummy SNES for AMG coarse solve
530ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
531d1d35e2fSjeremylt   if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
532e3e3df41Sjeremylt     // -- Jacobian Matrix
533d1d35e2fSjeremylt     ierr = DMSetMatType(level_dms[0], MATAIJ); CHKERRQ(ierr);
534d1d35e2fSjeremylt     ierr = DMCreateMatrix(level_dms[0], &jacob_mat_coarse); CHKERRQ(ierr);
535e3e3df41Sjeremylt 
536d1d35e2fSjeremylt     if (app_ctx->degree > 1) {
537d1d35e2fSjeremylt       ierr = SNESCreate(comm, &snes_coarse); CHKERRQ(ierr);
538d1d35e2fSjeremylt       ierr = SNESSetDM(snes_coarse, level_dms[0]); CHKERRQ(ierr);
539d1d35e2fSjeremylt       ierr = SNESSetSolution(snes_coarse, U_g[0]); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson 
541e3e3df41Sjeremylt       // -- Jacobian function
542d1d35e2fSjeremylt       ierr = SNESSetJacobian(snes_coarse, jacob_mat_coarse, jacob_mat_coarse, NULL,
543ccaff030SJeremy L Thompson                              NULL); CHKERRQ(ierr);
544ccaff030SJeremy L Thompson 
545ccaff030SJeremy L Thompson       // -- Residual evaluation function
546d1d35e2fSjeremylt       ierr = PetscMalloc1(1, &jacob_coarse_ctx); CHKERRQ(ierr);
547d1d35e2fSjeremylt       ierr = PetscMemcpy(jacob_coarse_ctx, jacob_ctx[0], sizeof(*jacob_ctx[0]));
548ccaff030SJeremy L Thompson       CHKERRQ(ierr);
549d1d35e2fSjeremylt       ierr = SNESSetFunction(snes_coarse, U_g[0], ApplyJacobianCoarse_Ceed,
550d1d35e2fSjeremylt                              jacob_coarse_ctx); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson 
552d1d35e2fSjeremylt       // -- Update form_jacob_ctx
553d1d35e2fSjeremylt       form_jacob_ctx->u_coarse = U_g[0];
554d1d35e2fSjeremylt       form_jacob_ctx->snes_coarse = snes_coarse;
555d1d35e2fSjeremylt       form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse;
556e3e3df41Sjeremylt     }
557e3e3df41Sjeremylt   }
558e3e3df41Sjeremylt 
559e3e3df41Sjeremylt   // Set Jacobian function
560d1d35e2fSjeremylt   if (app_ctx->degree > 1) {
561d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level],
562d1d35e2fSjeremylt                            FormJacobian, form_jacob_ctx); CHKERRQ(ierr);
563e3e3df41Sjeremylt   } else {
564d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse,
565e3e3df41Sjeremylt                            SNESComputeJacobianDefaultColor, NULL);
566e3e3df41Sjeremylt     CHKERRQ(ierr);
567e3e3df41Sjeremylt   }
568ccaff030SJeremy L Thompson 
569ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
570ccaff030SJeremy L Thompson   // Setup KSP
571ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
572ccaff030SJeremy L Thompson   {
573ccaff030SJeremy L Thompson     PC pc;
574ccaff030SJeremy L Thompson     KSP ksp;
575ccaff030SJeremy L Thompson 
576ccaff030SJeremy L Thompson     // -- KSP
577ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
578ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
579ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
580ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
581ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
582ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson 
584ccaff030SJeremy L Thompson     // -- Preconditioning
585ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
586d1d35e2fSjeremylt     ierr = PCSetDM(pc, level_dms[fine_level]); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson 
589d1d35e2fSjeremylt     if (app_ctx->multigrid_choice == MULTIGRID_NONE) {
590ccaff030SJeremy L Thompson       // ---- No Multigrid
591ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
593d1d35e2fSjeremylt     } else if (app_ctx->degree == 1) {
594f892e6d1Sjeremylt       // ---- AMG for degree 1
595f892e6d1Sjeremylt       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
596ccaff030SJeremy L Thompson     } else {
597ccaff030SJeremy L Thompson       // ---- PCMG
598ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
599ccaff030SJeremy L Thompson 
600ccaff030SJeremy L Thompson       // ------ PCMG levels
601d1d35e2fSjeremylt       ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
602d1d35e2fSjeremylt       for (PetscInt level = 0; level < num_levels; level++) {
603ccaff030SJeremy L Thompson         // -------- Smoother
604d1d35e2fSjeremylt         KSP ksp_smoother, ksp_est;
605d1d35e2fSjeremylt         PC pc_smoother;
606ccaff030SJeremy L Thompson 
607ccaff030SJeremy L Thompson         // ---------- Smoother KSP
608d1d35e2fSjeremylt         ierr = PCMGGetSmoother(pc, level, &ksp_smoother); CHKERRQ(ierr);
609d1d35e2fSjeremylt         ierr = KSPSetDM(ksp_smoother, level_dms[level]); CHKERRQ(ierr);
610d1d35e2fSjeremylt         ierr = KSPSetDMActive(ksp_smoother, PETSC_FALSE); CHKERRQ(ierr);
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson         // ---------- Chebyshev options
613d1d35e2fSjeremylt         ierr = KSPSetType(ksp_smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
614d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1);
615ccaff030SJeremy L Thompson         CHKERRQ(ierr);
616d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est); CHKERRQ(ierr);
617d1d35e2fSjeremylt         ierr = KSPSetType(ksp_est, KSPCG); CHKERRQ(ierr);
618d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE);
619ccaff030SJeremy L Thompson         CHKERRQ(ierr);
620d1d35e2fSjeremylt         ierr = KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]);
621ccaff030SJeremy L Thompson         CHKERRQ(ierr);
622ccaff030SJeremy L Thompson 
623ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
624d1d35e2fSjeremylt         ierr = KSPGetPC(ksp_smoother, &pc_smoother); CHKERRQ(ierr);
625d1d35e2fSjeremylt         ierr = PCSetType(pc_smoother, PCJACOBI); CHKERRQ(ierr);
626d1d35e2fSjeremylt         ierr = PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
627ccaff030SJeremy L Thompson 
628ccaff030SJeremy L Thompson         // -------- Work vector
629d1d35e2fSjeremylt         if (level != fine_level) {
630d1d35e2fSjeremylt           ierr = PCMGSetX(pc, level, U_g[level]); CHKERRQ(ierr);
631ccaff030SJeremy L Thompson         }
632ccaff030SJeremy L Thompson 
633ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
634ccaff030SJeremy L Thompson         if (level > 0) {
635d1d35e2fSjeremylt           ierr = PCMGSetInterpolation(pc, level, prolong_restr_mat[level]);
636ccaff030SJeremy L Thompson           CHKERRQ(ierr);
637d1d35e2fSjeremylt           ierr = PCMGSetRestriction(pc, level, prolong_restr_mat[level]);
638ccaff030SJeremy L Thompson           CHKERRQ(ierr);
639ccaff030SJeremy L Thompson         }
640ccaff030SJeremy L Thompson       }
641ccaff030SJeremy L Thompson 
642ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
643d1d35e2fSjeremylt       KSP ksp_coarse;
644d1d35e2fSjeremylt       PC pc_coarse;
645ccaff030SJeremy L Thompson 
646ccaff030SJeremy L Thompson       // -------- Coarse KSP
647d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
648d1d35e2fSjeremylt       ierr = KSPSetType(ksp_coarse, KSPPREONLY); CHKERRQ(ierr);
649d1d35e2fSjeremylt       ierr = KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse);
650ccaff030SJeremy L Thompson       CHKERRQ(ierr);
651d1d35e2fSjeremylt       ierr = KSPSetOptionsPrefix(ksp_coarse, "coarse_"); CHKERRQ(ierr);
652ccaff030SJeremy L Thompson 
653ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
654d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
655d1d35e2fSjeremylt       ierr = PCSetType(pc_coarse, PCGAMG); CHKERRQ(ierr);
656d1d35e2fSjeremylt       ierr = PCSetOptionsPrefix(pc_coarse, "coarse_"); CHKERRQ(ierr);
657ccaff030SJeremy L Thompson 
658d1d35e2fSjeremylt       ierr = KSPSetFromOptions(ksp_coarse); CHKERRQ(ierr);
659d1d35e2fSjeremylt       ierr = PCSetFromOptions(pc_coarse); CHKERRQ(ierr);
660ccaff030SJeremy L Thompson 
661ccaff030SJeremy L Thompson       // ------ PCMG options
662ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
663ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
664d1d35e2fSjeremylt       ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
665ccaff030SJeremy L Thompson     }
666ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
667ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
668ccaff030SJeremy L Thompson   }
669038d0cd7Sjeremylt   {
670038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
671d1d35e2fSjeremylt     SNESLineSearch line_search;
672481a4840SJed Brown 
673d1d35e2fSjeremylt     ierr = SNESGetLineSearch(snes, &line_search); CHKERRQ(ierr);
674d1d35e2fSjeremylt     ierr = SNESLineSearchSetType(line_search, SNESLINESEARCHCP); CHKERRQ(ierr);
675481a4840SJed Brown   }
676481a4840SJed Brown 
677ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
678ccaff030SJeremy L Thompson 
679ccaff030SJeremy L Thompson   // Performance logging
680ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
681ccaff030SJeremy L Thompson 
682ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
683ccaff030SJeremy L Thompson   // Set initial guess
684ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
685f81c27eaSJed Brown   ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr);
686ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
687ccaff030SJeremy L Thompson 
688ccaff030SJeremy L Thompson   // View solution
689d1d35e2fSjeremylt   if (app_ctx->view_soln) {
690d1d35e2fSjeremylt     ierr = ViewSolution(comm, app_ctx, U, 0, 0.0); CHKERRQ(ierr);
691ccaff030SJeremy L Thompson   }
692ccaff030SJeremy L Thompson 
693ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
694ccaff030SJeremy L Thompson   // Solve SNES
695ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
696d1d35e2fSjeremylt   PetscBool snes_monitor = PETSC_FALSE;
697d1d35e2fSjeremylt   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor);
6985f24f028Sjeremylt   CHKERRQ(ierr);
6995f24f028Sjeremylt 
700ccaff030SJeremy L Thompson   // Performance logging
701d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve);
702ccaff030SJeremy L Thompson   CHKERRQ(ierr);
703d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_solve); CHKERRQ(ierr);
704ccaff030SJeremy L Thompson 
705ccaff030SJeremy L Thompson   // Timing
706ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
707d1d35e2fSjeremylt   start_time = MPI_Wtime();
708ccaff030SJeremy L Thompson 
709ccaff030SJeremy L Thompson   // Solve for each load increment
7105f24f028Sjeremylt   PetscInt increment;
711d1d35e2fSjeremylt   for (increment = 1; increment <= app_ctx->num_increments; increment++) {
7125f24f028Sjeremylt     // -- Log increment count
713d1d35e2fSjeremylt     if (snes_monitor) {
7145f24f028Sjeremylt       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
7155f24f028Sjeremylt       CHKERRQ(ierr);
7165f24f028Sjeremylt     }
7175f24f028Sjeremylt 
718ccaff030SJeremy L Thompson     // -- Scale the problem
719d1d35e2fSjeremylt     PetscScalar load_increment = 1.0*increment / app_ctx->num_increments,
720d1d35e2fSjeremylt                 scalingFactor = load_increment /
721d1d35e2fSjeremylt                                 (increment == 1 ? 1 : res_ctx->load_increment);
722d1d35e2fSjeremylt     res_ctx->load_increment = load_increment;
723d1d35e2fSjeremylt     if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) {
724ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
725ccaff030SJeremy L Thompson     }
726ccaff030SJeremy L Thompson 
727ccaff030SJeremy L Thompson     // -- Solve
728ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
729ccaff030SJeremy L Thompson 
730ccaff030SJeremy L Thompson     // -- View solution
731d1d35e2fSjeremylt     if (app_ctx->view_soln) {
732d1d35e2fSjeremylt       ierr = ViewSolution(comm, app_ctx, U, increment, load_increment); CHKERRQ(ierr);
733ccaff030SJeremy L Thompson     }
734ccaff030SJeremy L Thompson 
735ccaff030SJeremy L Thompson     // -- Update SNES iteration count
736ccaff030SJeremy L Thompson     PetscInt its;
737ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
738bf0f51feSjeremylt     snes_its += its;
7397418ca88SJeremy L Thompson     ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr);
740bf0f51feSjeremylt     ksp_its += its;
7413951731eSjeremylt 
7423951731eSjeremylt     // -- Check for divergence
7433951731eSjeremylt     SNESConvergedReason reason;
7443951731eSjeremylt     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
7453951731eSjeremylt     if (reason < 0)
7463951731eSjeremylt       break;
747*8355605fSKaren (Ren) Stengel     if (app_ctx->energy_viewer) {
748*8355605fSKaren (Ren) Stengel       // -- Log strain energy for current load increment
749*8355605fSKaren (Ren) Stengel       CeedScalar energy;
750*8355605fSKaren (Ren) Stengel       ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
751*8355605fSKaren (Ren) Stengel                                  U, &energy); CHKERRQ(ierr);
752*8355605fSKaren (Ren) Stengel 
753*8355605fSKaren (Ren) Stengel       if (!app_ctx->test_mode) {
754*8355605fSKaren (Ren) Stengel         // -- Output
755*8355605fSKaren (Ren) Stengel         ierr = PetscPrintf(comm,
756*8355605fSKaren (Ren) Stengel                            "    Strain Energy                      : %.12e\n",
757*8355605fSKaren (Ren) Stengel                            energy); CHKERRQ(ierr);
758*8355605fSKaren (Ren) Stengel       }
759*8355605fSKaren (Ren) Stengel       ierr = PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment,
760*8355605fSKaren (Ren) Stengel                                     energy); CHKERRQ(ierr);
761*8355605fSKaren (Ren) Stengel     }
762ccaff030SJeremy L Thompson   }
763ccaff030SJeremy L Thompson 
764ccaff030SJeremy L Thompson   // Timing
765d1d35e2fSjeremylt   elapsed_time = MPI_Wtime() - start_time;
766ccaff030SJeremy L Thompson 
767ccaff030SJeremy L Thompson   // Performance logging
768ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
769ccaff030SJeremy L Thompson 
770ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
771ccaff030SJeremy L Thompson   // Output summary
772ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
773d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
774ccaff030SJeremy L Thompson     // -- SNES
775d1d35e2fSjeremylt     SNESType snes_type;
776ccaff030SJeremy L Thompson     SNESConvergedReason reason;
777ccaff030SJeremy L Thompson     PetscReal rnorm;
778d1d35e2fSjeremylt     ierr = SNESGetType(snes, &snes_type); CHKERRQ(ierr);
779ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
780ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
781ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
782ccaff030SJeremy L Thompson                        "  SNES:\n"
783ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
784ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
785ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
7865f24f028Sjeremylt                        "    Completed Load Increments          : %d\n"
787ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
788ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
789d1d35e2fSjeremylt                        snes_type, SNESConvergedReasons[reason],
790d1d35e2fSjeremylt                        app_ctx->num_increments, increment - 1,
791bf0f51feSjeremylt                        snes_its, (double)rnorm); CHKERRQ(ierr);
792ccaff030SJeremy L Thompson 
793ccaff030SJeremy L Thompson     // -- KSP
794ccaff030SJeremy L Thompson     KSP ksp;
795d1d35e2fSjeremylt     KSPType ksp_type;
796ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
797d1d35e2fSjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
798ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
799ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
8007418ca88SJeremy L Thompson                        "    KSP Type                           : %s\n"
8017418ca88SJeremy L Thompson                        "    Total KSP Iterations               : %D\n",
802bf0f51feSjeremylt                        ksp_type, ksp_its); CHKERRQ(ierr);
803ccaff030SJeremy L Thompson 
804ccaff030SJeremy L Thompson     // -- PC
805ccaff030SJeremy L Thompson     PC pc;
806d1d35e2fSjeremylt     PCType pc_type;
807ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
808d1d35e2fSjeremylt     ierr = PCGetType(pc, &pc_type); CHKERRQ(ierr);
809e3e3df41Sjeremylt     ierr = PetscPrintf(comm,
810e3e3df41Sjeremylt                        "    PC Type                            : %s\n",
811d1d35e2fSjeremylt                        pc_type); CHKERRQ(ierr);
812e3e3df41Sjeremylt 
813d1d35e2fSjeremylt     if (!strcmp(pc_type, PCMG)) {
814d1d35e2fSjeremylt       PCMGType pcmg_type;
815d1d35e2fSjeremylt       ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
816ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
817ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
818ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
819ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
820d1d35e2fSjeremylt                          PCMGTypes[pcmg_type],
821d1d35e2fSjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
822ccaff030SJeremy L Thompson 
823ccaff030SJeremy L Thompson       // -- Coarse Solve
824d1d35e2fSjeremylt       KSP ksp_coarse;
825d1d35e2fSjeremylt       PC pc_coarse;
826d1d35e2fSjeremylt       PCType pc_type;
827ccaff030SJeremy L Thompson 
828d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
829d1d35e2fSjeremylt       ierr = KSPGetType(ksp_coarse, &ksp_type); CHKERRQ(ierr);
830d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
831d1d35e2fSjeremylt       ierr = PCGetType(pc_coarse, &pc_type); CHKERRQ(ierr);
832ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
833ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
834ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
835ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
836d1d35e2fSjeremylt                          ksp_type, pc_type); CHKERRQ(ierr);
837ccaff030SJeremy L Thompson     }
838ccaff030SJeremy L Thompson   }
839ccaff030SJeremy L Thompson 
840ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
841ccaff030SJeremy L Thompson   // Compute solve time
842ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
843d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
844d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm);
845ccaff030SJeremy L Thompson     CHKERRQ(ierr);
846d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm);
847ccaff030SJeremy L Thompson     CHKERRQ(ierr);
848ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
849ccaff030SJeremy L Thompson                        "  Performance:\n"
8507418ca88SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n"
8517418ca88SJeremy L Thompson                        "    DoFs/Sec in SNES                   : %g (%g) million\n",
852bf0f51feSjeremylt                        max_time, min_time, 1e-6*U_g_size[fine_level]*ksp_its/max_time,
853bf0f51feSjeremylt                        1e-6*U_g_size[fine_level]*ksp_its/min_time); CHKERRQ(ierr);
854ccaff030SJeremy L Thompson   }
855ccaff030SJeremy L Thompson 
856ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
857ccaff030SJeremy L Thompson   // Compute error
858ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
859d1d35e2fSjeremylt   if (app_ctx->forcing_choice == FORCE_MMS) {
860d1d35e2fSjeremylt     CeedScalar l2_error = 1., l2_U_norm = 1.;
861d1d35e2fSjeremylt     const CeedScalar *true_array;
862d1d35e2fSjeremylt     Vec error_vec, true_vec;
863ccaff030SJeremy L Thompson 
864ccaff030SJeremy L Thompson     // -- Work vectors
865d1d35e2fSjeremylt     ierr = VecDuplicate(U, &error_vec); CHKERRQ(ierr);
866d1d35e2fSjeremylt     ierr = VecSet(error_vec, 0.0); CHKERRQ(ierr);
867d1d35e2fSjeremylt     ierr = VecDuplicate(U, &true_vec); CHKERRQ(ierr);
868d1d35e2fSjeremylt     ierr = VecSet(true_vec, 0.0); CHKERRQ(ierr);
869ccaff030SJeremy L Thompson 
870ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
871d1d35e2fSjeremylt     CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln,
872d1d35e2fSjeremylt                            CEED_MEM_HOST, &true_array);
873d1d35e2fSjeremylt     ierr = VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array);
87462e9c006SJeremy L Thompson     CHKERRQ(ierr);
875d1d35e2fSjeremylt     ierr = DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec);
876ccaff030SJeremy L Thompson     CHKERRQ(ierr);
877d1d35e2fSjeremylt     ierr = VecResetArray(res_ctx->Y_loc); CHKERRQ(ierr);
878d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array);
879ccaff030SJeremy L Thompson 
880ccaff030SJeremy L Thompson     // -- Compute L2 error
881d1d35e2fSjeremylt     ierr = VecWAXPY(error_vec, -1.0, U, true_vec); CHKERRQ(ierr);
882d1d35e2fSjeremylt     ierr = VecNorm(error_vec, NORM_2, &l2_error); CHKERRQ(ierr);
883d1d35e2fSjeremylt     ierr = VecNorm(U, NORM_2, &l2_U_norm); CHKERRQ(ierr);
884d1d35e2fSjeremylt     l2_error /= l2_U_norm;
885ccaff030SJeremy L Thompson 
886ccaff030SJeremy L Thompson     // -- Output
887d1d35e2fSjeremylt     if (!app_ctx->test_mode || l2_error > 0.05) {
888ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
889ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
890d1d35e2fSjeremylt                          l2_error); CHKERRQ(ierr);
891ccaff030SJeremy L Thompson     }
892ccaff030SJeremy L Thompson 
893ccaff030SJeremy L Thompson     // -- Cleanup
894d1d35e2fSjeremylt     ierr = VecDestroy(&error_vec); CHKERRQ(ierr);
895d1d35e2fSjeremylt     ierr = VecDestroy(&true_vec); CHKERRQ(ierr);
8962d93065eSjeremylt   }
8972d93065eSjeremylt 
8982d93065eSjeremylt   // ---------------------------------------------------------------------------
8992d93065eSjeremylt   // Compute energy
9002d93065eSjeremylt   // ---------------------------------------------------------------------------
901*8355605fSKaren (Ren) Stengel   PetscReal energy;
902d1d35e2fSjeremylt   ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
903a3c02c40SJeremy L Thompson                              U, &energy); CHKERRQ(ierr);
904*8355605fSKaren (Ren) Stengel   if (!app_ctx->test_mode) {
9052d93065eSjeremylt     // -- Output
9062d93065eSjeremylt     ierr = PetscPrintf(comm,
90772d03b64SArash Mehraban                        "    Strain Energy                      : %.12e\n",
9082d93065eSjeremylt                        energy); CHKERRQ(ierr);
909ccaff030SJeremy L Thompson   }
910*8355605fSKaren (Ren) Stengel   ierr = RegressionTests_solids(app_ctx, energy); CHKERRQ(ierr);
911ccaff030SJeremy L Thompson 
912ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
9135c25879aSJeremy L Thompson   // Output diagnostic quantities
9145c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
915d1d35e2fSjeremylt   if (app_ctx->view_soln || app_ctx->view_final_soln) {
9165c25879aSJeremy L Thompson     // -- Setup context
917d1d35e2fSjeremylt     UserMult diagnostic_ctx;
918d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &diagnostic_ctx); CHKERRQ(ierr);
919d1d35e2fSjeremylt     ierr = PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)); CHKERRQ(ierr);
920d1d35e2fSjeremylt     diagnostic_ctx->dm = dm_diagnostic;
921d1d35e2fSjeremylt     diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic;
9225c25879aSJeremy L Thompson 
9235c25879aSJeremy L Thompson     // -- Compute and output
924d1d35e2fSjeremylt     ierr = ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx,
925d1d35e2fSjeremylt                                     app_ctx, U,
926d1d35e2fSjeremylt                                     ceed_data[fine_level]->elem_restr_diagnostic);
9275c25879aSJeremy L Thompson     CHKERRQ(ierr);
9285c25879aSJeremy L Thompson 
9295c25879aSJeremy L Thompson     // -- Cleanup
930d1d35e2fSjeremylt     ierr = PetscFree(diagnostic_ctx); CHKERRQ(ierr);
9315c25879aSJeremy L Thompson   }
9325c25879aSJeremy L Thompson 
9335c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
934ccaff030SJeremy L Thompson   // Free objects
935ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
936ccaff030SJeremy L Thompson   // Data in arrays per level
937d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
938ccaff030SJeremy L Thompson     // Vectors
939d1d35e2fSjeremylt     ierr = VecDestroy(&U_g[level]); CHKERRQ(ierr);
940d1d35e2fSjeremylt     ierr = VecDestroy(&U_loc[level]); CHKERRQ(ierr);
941ccaff030SJeremy L Thompson 
942ccaff030SJeremy L Thompson     // Jacobian matrix and data
943d1d35e2fSjeremylt     ierr = VecDestroy(&jacob_ctx[level]->Y_loc); CHKERRQ(ierr);
944d1d35e2fSjeremylt     ierr = MatDestroy(&jacob_mat[level]); CHKERRQ(ierr);
945d1d35e2fSjeremylt     ierr = PetscFree(jacob_ctx[level]); CHKERRQ(ierr);
946ccaff030SJeremy L Thompson 
947ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
948ccaff030SJeremy L Thompson     if (level > 0) {
949d1d35e2fSjeremylt       ierr = PetscFree(prolong_restr_ctx[level]); CHKERRQ(ierr);
950d1d35e2fSjeremylt       ierr = MatDestroy(&prolong_restr_mat[level]); CHKERRQ(ierr);
951ccaff030SJeremy L Thompson     }
952ccaff030SJeremy L Thompson 
953ccaff030SJeremy L Thompson     // DM
954d1d35e2fSjeremylt     ierr = DMDestroy(&level_dms[level]); CHKERRQ(ierr);
955ccaff030SJeremy L Thompson 
956ccaff030SJeremy L Thompson     // libCEED objects
957d1d35e2fSjeremylt     ierr = CeedDataDestroy(level, ceed_data[level]); CHKERRQ(ierr);
958ccaff030SJeremy L Thompson   }
959ccaff030SJeremy L Thompson 
960*8355605fSKaren (Ren) Stengel   ierr = PetscViewerDestroy(&app_ctx->energy_viewer); CHKERRQ(ierr);
961*8355605fSKaren (Ren) Stengel 
962ccaff030SJeremy L Thompson   // Arrays
963d1d35e2fSjeremylt   ierr = PetscFree(U_g); CHKERRQ(ierr);
964d1d35e2fSjeremylt   ierr = PetscFree(U_loc); CHKERRQ(ierr);
965d1d35e2fSjeremylt   ierr = PetscFree(U_g_size); CHKERRQ(ierr);
966d1d35e2fSjeremylt   ierr = PetscFree(U_l_size); CHKERRQ(ierr);
967d1d35e2fSjeremylt   ierr = PetscFree(U_loc_size); CHKERRQ(ierr);
968d1d35e2fSjeremylt   ierr = PetscFree(jacob_ctx); CHKERRQ(ierr);
969d1d35e2fSjeremylt   ierr = PetscFree(jacob_mat); CHKERRQ(ierr);
970d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_ctx); CHKERRQ(ierr);
971d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_mat); CHKERRQ(ierr);
972d1d35e2fSjeremylt   ierr = PetscFree(app_ctx->level_degrees); CHKERRQ(ierr);
973d1d35e2fSjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
974ccaff030SJeremy L Thompson 
975ccaff030SJeremy L Thompson   // libCEED objects
976d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys);
977d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys_smoother);
978ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
979ccaff030SJeremy L Thompson 
980ccaff030SJeremy L Thompson   // PETSc objects
981ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
982ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
983d1d35e2fSjeremylt   ierr = VecDestroy(&R_loc); CHKERRQ(ierr);
984ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
985d1d35e2fSjeremylt   ierr = VecDestroy(&F_loc); CHKERRQ(ierr);
986d1d35e2fSjeremylt   ierr = VecDestroy(&neumann_bcs); CHKERRQ(ierr);
987d1d35e2fSjeremylt   ierr = VecDestroy(&bcs_loc); CHKERRQ(ierr);
988d1d35e2fSjeremylt   ierr = MatDestroy(&jacob_mat_coarse); CHKERRQ(ierr);
989ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
990d1d35e2fSjeremylt   ierr = SNESDestroy(&snes_coarse); CHKERRQ(ierr);
991d1d35e2fSjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
992d1d35e2fSjeremylt   ierr = DMDestroy(&dm_energy); CHKERRQ(ierr);
993d1d35e2fSjeremylt   ierr = DMDestroy(&dm_diagnostic); CHKERRQ(ierr);
994d1d35e2fSjeremylt   ierr = PetscFree(level_dms); CHKERRQ(ierr);
995ccaff030SJeremy L Thompson 
996ccaff030SJeremy L Thompson   // Structs
997d1d35e2fSjeremylt   ierr = PetscFree(res_ctx); CHKERRQ(ierr);
998d1d35e2fSjeremylt   ierr = PetscFree(form_jacob_ctx); CHKERRQ(ierr);
999d1d35e2fSjeremylt   ierr = PetscFree(jacob_coarse_ctx); CHKERRQ(ierr);
1000d1d35e2fSjeremylt   ierr = PetscFree(app_ctx); CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   ierr = PetscFree(phys); CHKERRQ(ierr);
1002*8355605fSKaren (Ren) Stengel   ierr = PetscFree(phys_MR); CHKERRQ(ierr);
1003d1d35e2fSjeremylt   ierr = PetscFree(phys_smoother); CHKERRQ(ierr);
1004ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1005ccaff030SJeremy L Thompson 
1006ccaff030SJeremy L Thompson   return PetscFinalize();
1007ccaff030SJeremy L Thompson }
1008