xref: /libCEED/examples/solids/elasticity.c (revision caae6498fc183b53f040f3f51aad74fe34626b10)
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 //
368355605fSKaren (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
378355605fSKaren (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
388355605fSKaren (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
525754ecacSJeremy L Thompson   ProblemFunctions problem_functions;      // Setup functions for each problem
53ccaff030SJeremy L Thompson   Units          units;                    // Contains units scaling
54ccaff030SJeremy L Thompson   // PETSc objects
55d1d35e2fSjeremylt   PetscLogStage  stage_dm_setup, stage_libceed_setup,
56d1d35e2fSjeremylt                  stage_snes_setup, stage_snes_solve;
57d1d35e2fSjeremylt   DM             dm_orig;                  // Distributed DM to clone
58d1d35e2fSjeremylt   DM             dm_energy, dm_diagnostic; // DMs for postprocessing
59d1d35e2fSjeremylt   DM             *level_dms;
60d1d35e2fSjeremylt   Vec            U, *U_g, *U_loc;          // U: solution, R: residual, F: forcing
61d1d35e2fSjeremylt   Vec            R, R_loc, F, F_loc;       // g: global, loc: local
62d1d35e2fSjeremylt   Vec            neumann_bcs = NULL, bcs_loc = NULL;
6317f0843fSJeremy L Thompson   SNES           snes;
64d1d35e2fSjeremylt   Mat            *jacob_mat, jacob_mat_coarse, *prolong_restr_mat;
65ccaff030SJeremy L Thompson   // PETSc data
66d1d35e2fSjeremylt   UserMult       res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx;
67d1d35e2fSjeremylt   FormJacobCtx   form_jacob_ctx;
68d1d35e2fSjeremylt   UserMultProlongRestr *prolong_restr_ctx;
69d1d35e2fSjeremylt   PCMGCycleType  pcmg_cycle_type = PC_MG_CYCLE_V;
70ccaff030SJeremy L Thompson   // libCEED objects
71d99fa3c5SJeremy L Thompson   Ceed           ceed;
72d1d35e2fSjeremylt   CeedData       *ceed_data;
73d1d35e2fSjeremylt   CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL;
74ccaff030SJeremy L Thompson   // Parameters
75d1d35e2fSjeremylt   PetscInt       num_comp_u = 3;                 // 3 DoFs in 3D
76d1d35e2fSjeremylt   PetscInt       num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic
77d1d35e2fSjeremylt   PetscInt       num_levels = 1, fine_level = 0;
78bf0f51feSjeremylt   PetscInt       *U_g_size, *U_l_size, *U_loc_size;
79bf0f51feSjeremylt   PetscInt       snes_its = 0, ksp_its = 0;
80d1d35e2fSjeremylt   double         start_time, elapsed_time, min_time, max_time;
81ccaff030SJeremy L Thompson 
82ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
83bf0f51feSjeremylt   if (ierr) return ierr;
84ccaff030SJeremy L Thompson 
85ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
86ccaff030SJeremy L Thompson   // Process command line options
87ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
88ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
89ccaff030SJeremy L Thompson 
90ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
91d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr);
92d1d35e2fSjeremylt   ierr = ProcessCommandLineOptions(comm, app_ctx); CHKERRQ(ierr);
935754ecacSJeremy L Thompson   ierr = PetscCalloc1(1, &problem_functions); CHKERRQ(ierr);
945754ecacSJeremy L Thompson   ierr = RegisterProblems(problem_functions); CHKERRQ(ierr);
95d1d35e2fSjeremylt   num_levels = app_ctx->num_levels;
96d1d35e2fSjeremylt   fine_level = num_levels - 1;
97ccaff030SJeremy L Thompson 
98ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
99d1d35e2fSjeremylt   // Initialize libCEED
10062e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
101d1d35e2fSjeremylt   // Initialize backend
102d1d35e2fSjeremylt   CeedInit(app_ctx->ceed_resource, &ceed);
10362e9c006SJeremy L Thompson 
10462e9c006SJeremy L Thompson   // Check preferred MemType
105d1d35e2fSjeremylt   CeedMemType mem_type_backend;
106d1d35e2fSjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1075754ecacSJeremy L Thompson   // Setup physics context and wrap in libCEED object
1085754ecacSJeremy L Thompson   {
1095754ecacSJeremy L Thompson     PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *);
1105754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name,
1115754ecacSJeremy L Thompson                                  &SetupPhysics); CHKERRQ(ierr);
1125754ecacSJeremy L Thompson     if (!SetupPhysics)
1135754ecacSJeremy L Thompson       SETERRQ1(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found",
1145754ecacSJeremy L Thompson                app_ctx->name);
1155754ecacSJeremy L Thompson     ierr = (*SetupPhysics)(comm, ceed, &units, &ctx_phys); CHKERRQ(ierr);
1165754ecacSJeremy L Thompson     PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext,
1175754ecacSJeremy L Thompson                                            CeedQFunctionContext *);
1185754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupSmootherPhysics,
1195754ecacSJeremy L Thompson                                  app_ctx->name, &SetupSmootherPhysics);
1205754ecacSJeremy L Thompson     CHKERRQ(ierr);
1215754ecacSJeremy L Thompson     if (!SetupSmootherPhysics)
1225754ecacSJeremy L Thompson       SETERRQ1(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found",
1235754ecacSJeremy L Thompson                app_ctx->name);
1245754ecacSJeremy L Thompson     ierr = (*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother);
1255754ecacSJeremy L Thompson     CHKERRQ(ierr);
126777ff853SJeremy L Thompson   }
127777ff853SJeremy L Thompson 
12862e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
129ccaff030SJeremy L Thompson   // Setup DM
130ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
131ccaff030SJeremy L Thompson   // Performance logging
132d1d35e2fSjeremylt   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup);
133ccaff030SJeremy L Thompson   CHKERRQ(ierr);
134d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_dm_setup); CHKERRQ(ierr);
135ccaff030SJeremy L Thompson 
136ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
137d1d35e2fSjeremylt   ierr = CreateDistributedDM(comm, app_ctx, &dm_orig); CHKERRQ(ierr);
138b68a8d79SJed Brown   VecType vectype;
139d1d35e2fSjeremylt   switch (mem_type_backend) {
140b68a8d79SJed Brown   case CEED_MEM_HOST: vectype = VECSTANDARD; break;
141b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
142b68a8d79SJed Brown     const char *resolved;
143b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
144b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
145b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
146b68a8d79SJed Brown     else vectype = VECSTANDARD;
147b68a8d79SJed Brown   }
148b68a8d79SJed Brown   }
149d1d35e2fSjeremylt   ierr = DMSetVecType(dm_orig, vectype); CHKERRQ(ierr);
150d1d35e2fSjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
151ccaff030SJeremy L Thompson 
152ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
153d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &level_dms); CHKERRQ(ierr);
154d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
155d1d35e2fSjeremylt     ierr = DMClone(dm_orig, &level_dms[level]); CHKERRQ(ierr);
156d1d35e2fSjeremylt     ierr = DMGetVecType(dm_orig, &vectype); CHKERRQ(ierr);
157d1d35e2fSjeremylt     ierr = DMSetVecType(level_dms[level], vectype); CHKERRQ(ierr);
158d1d35e2fSjeremylt     ierr = SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level],
159d1d35e2fSjeremylt                            PETSC_TRUE, num_comp_u); CHKERRQ(ierr);
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;
163d1d35e2fSjeremylt     ierr = DMGetLocalSection(level_dms[level], &section); CHKERRQ(ierr);
164a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr);
165a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
166a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
167a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
168a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
169a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
170a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
171a3c02c40SJeremy L Thompson   }
172a3c02c40SJeremy L Thompson 
173a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
174d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_energy); CHKERRQ(ierr);
175d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level],
176d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_e); CHKERRQ(ierr);
177d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_diagnostic); CHKERRQ(ierr);
178d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_diagnostic, app_ctx,
179d1d35e2fSjeremylt                          app_ctx->level_degrees[fine_level],
180d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_u + num_comp_d); CHKERRQ(ierr);
181d1d35e2fSjeremylt   ierr = DMSetVecType(dm_energy, vectype); CHKERRQ(ierr);
182d1d35e2fSjeremylt   ierr = DMSetVecType(dm_diagnostic, vectype); CHKERRQ(ierr);
183a3c02c40SJeremy L Thompson   {
184a3c02c40SJeremy L Thompson     // -- Label field components for viewing
185a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
186a3c02c40SJeremy L Thompson     PetscSection section;
187d1d35e2fSjeremylt     ierr = DMGetLocalSection(dm_diagnostic, &section); CHKERRQ(ierr);
188a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr);
1898ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
1908ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1918ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
1928ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1938ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
1948ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
195379387d4SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure");
1968ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1977ab18febSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain");
19813fdad4bSJeremy L Thompson     CHKERRQ(ierr);
19913fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2");
20013fdad4bSJeremy L Thompson     CHKERRQ(ierr);
20113fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 6, "detJ");
20213fdad4bSJeremy L Thompson     CHKERRQ(ierr);
20313fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity");
204a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
205ccaff030SJeremy L Thompson   }
206ccaff030SJeremy L Thompson 
207ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
208ccaff030SJeremy L Thompson   // Setup solution and work vectors
209ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
210ccaff030SJeremy L Thompson   // Allocate arrays
211d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g); CHKERRQ(ierr);
212d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc); CHKERRQ(ierr);
213d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g_size); CHKERRQ(ierr);
214d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_l_size); CHKERRQ(ierr);
215d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc_size); CHKERRQ(ierr);
216ccaff030SJeremy L Thompson 
217ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
218d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
219ccaff030SJeremy L Thompson     // -- Create global unknown vector U
220d1d35e2fSjeremylt     ierr = DMCreateGlobalVector(level_dms[level], &U_g[level]); CHKERRQ(ierr);
221d1d35e2fSjeremylt     ierr = VecGetSize(U_g[level], &U_g_size[level]); CHKERRQ(ierr);
222ccaff030SJeremy L Thompson     // Note: Local size for matShell
223d1d35e2fSjeremylt     ierr = VecGetLocalSize(U_g[level], &U_l_size[level]); CHKERRQ(ierr);
224ccaff030SJeremy L Thompson 
225d1d35e2fSjeremylt     // -- Create local unknown vector U_loc
226d1d35e2fSjeremylt     ierr = DMCreateLocalVector(level_dms[level], &U_loc[level]); CHKERRQ(ierr);
227ccaff030SJeremy L Thompson     // Note: local size for libCEED
228d1d35e2fSjeremylt     ierr = VecGetSize(U_loc[level], &U_loc_size[level]); CHKERRQ(ierr);
229ccaff030SJeremy L Thompson   }
230ccaff030SJeremy L Thompson 
231ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
232d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &U); CHKERRQ(ierr);
233d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &R); CHKERRQ(ierr);
234d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &F); CHKERRQ(ierr);
235d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &R_loc); CHKERRQ(ierr);
236d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &F_loc); CHKERRQ(ierr);
237ccaff030SJeremy L Thompson 
238ccaff030SJeremy L Thompson   // Performance logging
239ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
240ccaff030SJeremy L Thompson 
241ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
242ccaff030SJeremy L Thompson   // Set up libCEED
243ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
244ccaff030SJeremy L Thompson   // Performance logging
245d1d35e2fSjeremylt   ierr = PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup);
246ccaff030SJeremy L Thompson   CHKERRQ(ierr);
247d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_libceed_setup); CHKERRQ(ierr);
248ccaff030SJeremy L Thompson 
249ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
250d1d35e2fSjeremylt   CeedVector force_ceed;
251ccaff030SJeremy L Thompson   CeedScalar *f;
252d1d35e2fSjeremylt   PetscMemType force_mem_type;
253d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
254d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(F_loc, &f, &force_mem_type); CHKERRQ(ierr);
255d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed);
256d1d35e2fSjeremylt     CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f);
257ccaff030SJeremy L Thompson   }
258ccaff030SJeremy L Thompson 
259fe394131Sjeremylt   // -- Create libCEED local Neumann BCs vector
260d1d35e2fSjeremylt   CeedVector neumann_ceed;
261fe394131Sjeremylt   CeedScalar *n;
262d1d35e2fSjeremylt   PetscMemType nummann_mem_type;
263d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
264d1d35e2fSjeremylt     ierr = VecDuplicate(U, &neumann_bcs); CHKERRQ(ierr);
265d1d35e2fSjeremylt     ierr = VecDuplicate(U_loc[fine_level], &bcs_loc); CHKERRQ(ierr);
266d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type); CHKERRQ(ierr);
267d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed);
268d1d35e2fSjeremylt     CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type),
269fe394131Sjeremylt                        CEED_USE_POINTER, n);
270fe394131Sjeremylt   }
271fe394131Sjeremylt 
272ccaff030SJeremy L Thompson   // -- Setup libCEED objects
273d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
274d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
275d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &ceed_data[fine_level]); CHKERRQ(ierr);
276ccaff030SJeremy L Thompson   {
2775754ecacSJeremy L Thompson     PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx,
2785754ecacSJeremy L Thompson                                             CeedQFunctionContext, PetscInt,
2795754ecacSJeremy L Thompson                                             PetscInt, PetscInt, PetscInt,
2805754ecacSJeremy L Thompson                                             CeedVector, CeedVector, CeedData *);
2815754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupLibceedFineLevel,
2825754ecacSJeremy L Thompson                                  app_ctx->name, &SetupLibceedFineLevel);
2835754ecacSJeremy L Thompson     CHKERRQ(ierr);
2845754ecacSJeremy L Thompson     if (!SetupLibceedFineLevel)
2855754ecacSJeremy L Thompson       SETERRQ1(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found",
2865754ecacSJeremy L Thompson                app_ctx->name);
2875754ecacSJeremy L Thompson     ierr = (*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic,
2885754ecacSJeremy L Thompson                                     ceed, app_ctx, ctx_phys, fine_level,
2895754ecacSJeremy L Thompson                                     num_comp_u, U_g_size[fine_level],
2905754ecacSJeremy L Thompson                                     U_loc_size[fine_level],
2915754ecacSJeremy L Thompson                                     force_ceed, neumann_ceed, ceed_data);
292ccaff030SJeremy L Thompson     CHKERRQ(ierr);
293ccaff030SJeremy L Thompson   }
294d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
295d1d35e2fSjeremylt   for (PetscInt level = num_levels - 2; level >= 0; level--) {
296d1d35e2fSjeremylt     ierr = PetscCalloc1(1, &ceed_data[level]); CHKERRQ(ierr);
297d99fa3c5SJeremy L Thompson 
298d99fa3c5SJeremy L Thompson     // Get global communication restriction
299d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
300d1d35e2fSjeremylt     ierr = VecSet(U_loc[level+1], 1.0); CHKERRQ(ierr);
301d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[level+1], U_loc[level+1], ADD_VALUES,
302d1d35e2fSjeremylt                            U_g[level+1]); CHKERRQ(ierr);
303d1d35e2fSjeremylt     ierr = DMGlobalToLocal(level_dms[level+1], U_g[level+1], INSERT_VALUES,
304d1d35e2fSjeremylt                            U_loc[level+1]); CHKERRQ(ierr);
305d99fa3c5SJeremy L Thompson 
306d99fa3c5SJeremy L Thompson     // Place in libCEED array
307d99fa3c5SJeremy L Thompson     const PetscScalar *m;
308d1d35e2fSjeremylt     PetscMemType m_mem_type;
309d1d35e2fSjeremylt     ierr = VecGetArrayReadAndMemType(U_loc[level+1], &m, &m_mem_type);
310d1d35e2fSjeremylt     CHKERRQ(ierr);
311d1d35e2fSjeremylt     CeedVectorSetArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
312d99fa3c5SJeremy L Thompson                        CEED_USE_POINTER, (CeedScalar *)m);
313ccaff030SJeremy L Thompson 
3141c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
3155754ecacSJeremy L Thompson     PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt,
3165754ecacSJeremy L Thompson                                         PetscInt, PetscInt, PetscInt, CeedVector, CeedData *);
3175754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupLibceedLevel,
3185754ecacSJeremy L Thompson                                  app_ctx->name, &SetupLibceedLevel);
3195754ecacSJeremy L Thompson     CHKERRQ(ierr);
320*caae6498SJed Brown     if (!SetupLibceedLevel)
321*caae6498SJed Brown       SETERRQ1(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found",
322*caae6498SJed Brown                app_ctx->name);
3235754ecacSJeremy L Thompson     ierr = (*SetupLibceedLevel)(level_dms[level], ceed, app_ctx,
3245754ecacSJeremy L Thompson                                 level, num_comp_u, U_g_size[level],
3255754ecacSJeremy L Thompson                                 U_loc_size[level], ceed_data[level+1]->x_ceed,
3265754ecacSJeremy L Thompson                                 ceed_data);
327777ff853SJeremy L Thompson     CHKERRQ(ierr);
328d99fa3c5SJeremy L Thompson 
329d99fa3c5SJeremy L Thompson     // Restore PETSc vector
330d1d35e2fSjeremylt     CeedVectorTakeArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
331d99fa3c5SJeremy L Thompson                         (CeedScalar **)&m);
332d1d35e2fSjeremylt     ierr = VecRestoreArrayReadAndMemType(U_loc[level+1], &m); CHKERRQ(ierr);
333d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
334d1d35e2fSjeremylt     ierr = VecZeroEntries(U_loc[level+1]); CHKERRQ(ierr);
335ccaff030SJeremy L Thompson   }
336ccaff030SJeremy L Thompson 
337ccaff030SJeremy L Thompson   // Performance logging
338ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
339ccaff030SJeremy L Thompson 
340ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
341fe394131Sjeremylt   // Setup global forcing and Neumann BC vectors
342ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
343ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
344ccaff030SJeremy L Thompson 
345d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
346d1d35e2fSjeremylt     CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL);
347d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(F_loc, &f); CHKERRQ(ierr);
348d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F);
349ccaff030SJeremy L Thompson     CHKERRQ(ierr);
350d1d35e2fSjeremylt     CeedVectorDestroy(&force_ceed);
351ccaff030SJeremy L Thompson   }
352ccaff030SJeremy L Thompson 
353d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
354d1d35e2fSjeremylt     ierr = VecZeroEntries(neumann_bcs); CHKERRQ(ierr);
355d1d35e2fSjeremylt     CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL);
356d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(bcs_loc, &n); CHKERRQ(ierr);
357d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs);
358fe394131Sjeremylt     CHKERRQ(ierr);
359d1d35e2fSjeremylt     CeedVectorDestroy(&neumann_ceed);
360fe394131Sjeremylt   }
361fe394131Sjeremylt 
362ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
363ccaff030SJeremy L Thompson   // Print problem summary
364ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
365d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
366ccaff030SJeremy L Thompson     const char *usedresource;
367ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
36877841947SLeila Ghaffari     char hostname[PETSC_MAX_PATH_LEN];
36977841947SLeila Ghaffari     ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
37077841947SLeila Ghaffari     PetscInt comm_size;
37177841947SLeila Ghaffari     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson 
373ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
37417fd649aSJeremy L Thompson                        "\n-- Elasticity Example - libCEED + PETSc --\n"
37577841947SLeila Ghaffari                        "  MPI:\n"
37677841947SLeila Ghaffari                        "    Hostname                           : %s\n"
37777841947SLeila Ghaffari                        "    Total ranks                        : %d\n"
378ccaff030SJeremy L Thompson                        "  libCEED:\n"
37962e9c006SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
380b68a8d79SJed Brown                        "    libCEED Backend MemType            : %s\n",
38177841947SLeila Ghaffari                        hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend]);
38262e9c006SJeremy L Thompson     CHKERRQ(ierr);
383ccaff030SJeremy L Thompson 
38462e9c006SJeremy L Thompson     VecType vecType;
38562e9c006SJeremy L Thompson     ierr = VecGetType(U, &vecType); CHKERRQ(ierr);
38662e9c006SJeremy L Thompson     ierr = PetscPrintf(comm,
38762e9c006SJeremy L Thompson                        "  PETSc:\n"
38862e9c006SJeremy L Thompson                        "    PETSc Vec Type                     : %s\n",
38962e9c006SJeremy L Thompson                        vecType); CHKERRQ(ierr);
39062e9c006SJeremy L Thompson 
391ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
392ccaff030SJeremy L Thompson                        "  Problem:\n"
393ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
394ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
395ccaff030SJeremy L Thompson                        "  Mesh:\n"
396ccaff030SJeremy L Thompson                        "    File                               : %s\n"
397ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
398ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
399ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
400ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
401ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
402ccaff030SJeremy L Thompson                        "  Multigrid:\n"
403ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
404ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
4055754ecacSJeremy L Thompson                        app_ctx->name_for_disp,
406d1d35e2fSjeremylt                        forcing_types_for_disp[app_ctx->forcing_choice],
407d1d35e2fSjeremylt                        app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh",
408d1d35e2fSjeremylt                        app_ctx->degree + 1, app_ctx->degree + 1,
4095754ecacSJeremy L Thompson                        U_g_size[fine_level]/num_comp_u, U_l_size[fine_level]/num_comp_u,
4105754ecacSJeremy L Thompson                        num_comp_u,
411d1d35e2fSjeremylt                        (app_ctx->degree == 1 &&
412d1d35e2fSjeremylt                         app_ctx->multigrid_choice != MULTIGRID_NONE) ?
413f892e6d1Sjeremylt                        "Algebraic multigrid" :
414d1d35e2fSjeremylt                        multigrid_types_for_disp[app_ctx->multigrid_choice],
415d1d35e2fSjeremylt                        (app_ctx->degree == 1 ||
416d1d35e2fSjeremylt                         app_ctx->multigrid_choice == MULTIGRID_NONE) ?
417d1d35e2fSjeremylt                        0 : num_levels); CHKERRQ(ierr);
418ccaff030SJeremy L Thompson 
419d1d35e2fSjeremylt     if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
420e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
421d1d35e2fSjeremylt         CeedInt level = i ? fine_level : 0;
422ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
423ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
424ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
425ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
426ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
427ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
428d1d35e2fSjeremylt                            app_ctx->level_degrees[level] + 1,
429d1d35e2fSjeremylt                            U_g_size[level]/num_comp_u, U_l_size[level]/num_comp_u);
430ccaff030SJeremy L Thompson         CHKERRQ(ierr);
431ccaff030SJeremy L Thompson       }
432ccaff030SJeremy L Thompson     }
433ccaff030SJeremy L Thompson   }
434ccaff030SJeremy L Thompson 
435ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
436ccaff030SJeremy L Thompson   // Setup SNES
437ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
438ccaff030SJeremy L Thompson   // Performance logging
439d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup);
440ccaff030SJeremy L Thompson   CHKERRQ(ierr);
441d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_setup); CHKERRQ(ierr);
442ccaff030SJeremy L Thompson 
443ccaff030SJeremy L Thompson   // Create SNES
444ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
445d1d35e2fSjeremylt   ierr = SNESSetDM(snes, level_dms[fine_level]); CHKERRQ(ierr);
446ccaff030SJeremy L Thompson 
447ccaff030SJeremy L Thompson   // -- Jacobian evaluators
448d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_ctx); CHKERRQ(ierr);
449d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_mat); CHKERRQ(ierr);
450d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
451ccaff030SJeremy L Thompson     // -- Jacobian context for level
452d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &jacob_ctx[level]); CHKERRQ(ierr);
453d1d35e2fSjeremylt     ierr = SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level],
454d1d35e2fSjeremylt                             U_loc[level], ceed_data[level], ceed, ctx_phys,
455d1d35e2fSjeremylt                             ctx_phys_smoother, jacob_ctx[level]); CHKERRQ(ierr);
456ccaff030SJeremy L Thompson 
457ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
458d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level],
459d1d35e2fSjeremylt                           U_g_size[level], jacob_ctx[level], &jacob_mat[level]);
460ccaff030SJeremy L Thompson     CHKERRQ(ierr);
461d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_MULT,
462ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
463ccaff030SJeremy L Thompson     CHKERRQ(ierr);
464d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL,
465ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
466d1d35e2fSjeremylt     ierr = MatShellSetVecType(jacob_mat[level], vectype); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   }
468e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
469ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
470d1d35e2fSjeremylt   ierr = PetscMalloc1(1, &form_jacob_ctx); CHKERRQ(ierr);
471d1d35e2fSjeremylt   form_jacob_ctx->jacob_ctx = jacob_ctx;
472d1d35e2fSjeremylt   form_jacob_ctx->num_levels = num_levels;
473d1d35e2fSjeremylt   form_jacob_ctx->jacob_mat = jacob_mat;
474ccaff030SJeremy L Thompson 
475ccaff030SJeremy L Thompson   // -- Residual evaluation function
476d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &res_ctx); CHKERRQ(ierr);
477d1d35e2fSjeremylt   ierr = PetscMemcpy(res_ctx, jacob_ctx[fine_level],
478d1d35e2fSjeremylt                      sizeof(*jacob_ctx[fine_level])); CHKERRQ(ierr);
4795754ecacSJeremy L Thompson   res_ctx->op = ceed_data[fine_level]->op_residual;
4805754ecacSJeremy L Thompson   res_ctx->qf = ceed_data[fine_level]->qf_residual;
481d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0)
482d1d35e2fSjeremylt     res_ctx->neumann_bcs = neumann_bcs;
483fe394131Sjeremylt   else
484d1d35e2fSjeremylt     res_ctx->neumann_bcs = NULL;
485d1d35e2fSjeremylt   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx); CHKERRQ(ierr);
486ccaff030SJeremy L Thompson 
487ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
488d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_ctx); CHKERRQ(ierr);
489d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_mat); CHKERRQ(ierr);
490d1d35e2fSjeremylt   for (PetscInt level = 1; level < num_levels; level++) {
491ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
492d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &prolong_restr_ctx[level]); CHKERRQ(ierr);
493d1d35e2fSjeremylt     ierr = SetupProlongRestrictCtx(comm, app_ctx, level_dms[level-1],
494d1d35e2fSjeremylt                                    level_dms[level], U_g[level], U_loc[level-1],
495d1d35e2fSjeremylt                                    U_loc[level], ceed_data[level-1],
496d1d35e2fSjeremylt                                    ceed_data[level], ceed,
497d1d35e2fSjeremylt                                    prolong_restr_ctx[level]); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson 
499ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
500d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level-1], U_g_size[level],
501d1d35e2fSjeremylt                           U_g_size[level-1], prolong_restr_ctx[level],
502d1d35e2fSjeremylt                           &prolong_restr_mat[level]); CHKERRQ(ierr);
503ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
504d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT,
505ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
506d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE,
507ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
508ccaff030SJeremy L Thompson     CHKERRQ(ierr);
509d1d35e2fSjeremylt     ierr = MatShellSetVecType(prolong_restr_mat[level], vectype); CHKERRQ(ierr);
510ccaff030SJeremy L Thompson   }
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
51317f0843fSJeremy L Thompson   // Setup for AMG coarse solve
514ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
515d1d35e2fSjeremylt   if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
516e3e3df41Sjeremylt     // -- Jacobian Matrix
517d1d35e2fSjeremylt     ierr = DMCreateMatrix(level_dms[0], &jacob_mat_coarse); CHKERRQ(ierr);
518e3e3df41Sjeremylt 
519d1d35e2fSjeremylt     if (app_ctx->degree > 1) {
52017f0843fSJeremy L Thompson       // -- Assemble sparsity pattern
52117f0843fSJeremy L Thompson       CeedInt num_entries, *rows, *cols;
52217f0843fSJeremy L Thompson       CeedVector coo_values;
52317f0843fSJeremy L Thompson       CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries,
52417f0843fSJeremy L Thompson                                          &rows, &cols);
52517f0843fSJeremy L Thompson       ISLocalToGlobalMapping ltog_row, ltog_col;
52617f0843fSJeremy L Thompson       ierr = MatGetLocalToGlobalMapping(jacob_mat_coarse, &ltog_row, &ltog_col);
527ccaff030SJeremy L Thompson       CHKERRQ(ierr);
52817f0843fSJeremy L Thompson       ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows);
52917f0843fSJeremy L Thompson       CHKERRQ(ierr);
53017f0843fSJeremy L Thompson       ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols);
53117f0843fSJeremy L Thompson       CHKERRQ(ierr);
53217f0843fSJeremy L Thompson       ierr = MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows, cols);
53317f0843fSJeremy L Thompson       CHKERRQ(ierr);
53417f0843fSJeremy L Thompson       free(rows);
53517f0843fSJeremy L Thompson       free(cols);
53617f0843fSJeremy L Thompson       CeedVectorCreate(ceed, num_entries, &coo_values);
537ccaff030SJeremy L Thompson 
538d1d35e2fSjeremylt       // -- Update form_jacob_ctx
53917f0843fSJeremy L Thompson       form_jacob_ctx->coo_values = coo_values;
54017f0843fSJeremy L Thompson       form_jacob_ctx->op_coarse = ceed_data[0]->op_jacobian;
541d1d35e2fSjeremylt       form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse;
542e3e3df41Sjeremylt     }
543e3e3df41Sjeremylt   }
544e3e3df41Sjeremylt 
545e3e3df41Sjeremylt   // Set Jacobian function
546d1d35e2fSjeremylt   if (app_ctx->degree > 1) {
547d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level],
548d1d35e2fSjeremylt                            FormJacobian, form_jacob_ctx); CHKERRQ(ierr);
549e3e3df41Sjeremylt   } else {
550d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse,
551e3e3df41Sjeremylt                            SNESComputeJacobianDefaultColor, NULL);
552e3e3df41Sjeremylt     CHKERRQ(ierr);
553e3e3df41Sjeremylt   }
554ccaff030SJeremy L Thompson 
555ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
556ccaff030SJeremy L Thompson   // Setup KSP
557ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
558ccaff030SJeremy L Thompson   {
559ccaff030SJeremy L Thompson     PC pc;
560ccaff030SJeremy L Thompson     KSP ksp;
561ccaff030SJeremy L Thompson 
562ccaff030SJeremy L Thompson     // -- KSP
563ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
565ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
567ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
568ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
569ccaff030SJeremy L Thompson 
570ccaff030SJeremy L Thompson     // -- Preconditioning
571ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
572d1d35e2fSjeremylt     ierr = PCSetDM(pc, level_dms[fine_level]); CHKERRQ(ierr);
573ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
574ccaff030SJeremy L Thompson 
575d1d35e2fSjeremylt     if (app_ctx->multigrid_choice == MULTIGRID_NONE) {
576ccaff030SJeremy L Thompson       // ---- No Multigrid
577ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
578ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
579d1d35e2fSjeremylt     } else if (app_ctx->degree == 1) {
580f892e6d1Sjeremylt       // ---- AMG for degree 1
581f892e6d1Sjeremylt       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
582ccaff030SJeremy L Thompson     } else {
583ccaff030SJeremy L Thompson       // ---- PCMG
584ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
585ccaff030SJeremy L Thompson 
586ccaff030SJeremy L Thompson       // ------ PCMG levels
587d1d35e2fSjeremylt       ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
588d1d35e2fSjeremylt       for (PetscInt level = 0; level < num_levels; level++) {
589ccaff030SJeremy L Thompson         // -------- Smoother
590d1d35e2fSjeremylt         KSP ksp_smoother, ksp_est;
591d1d35e2fSjeremylt         PC pc_smoother;
592ccaff030SJeremy L Thompson 
593ccaff030SJeremy L Thompson         // ---------- Smoother KSP
594d1d35e2fSjeremylt         ierr = PCMGGetSmoother(pc, level, &ksp_smoother); CHKERRQ(ierr);
595d1d35e2fSjeremylt         ierr = KSPSetDM(ksp_smoother, level_dms[level]); CHKERRQ(ierr);
596d1d35e2fSjeremylt         ierr = KSPSetDMActive(ksp_smoother, PETSC_FALSE); CHKERRQ(ierr);
597ccaff030SJeremy L Thompson 
598ccaff030SJeremy L Thompson         // ---------- Chebyshev options
599d1d35e2fSjeremylt         ierr = KSPSetType(ksp_smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
600d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1);
601ccaff030SJeremy L Thompson         CHKERRQ(ierr);
602d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est); CHKERRQ(ierr);
603d1d35e2fSjeremylt         ierr = KSPSetType(ksp_est, KSPCG); CHKERRQ(ierr);
604d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE);
605ccaff030SJeremy L Thompson         CHKERRQ(ierr);
606d1d35e2fSjeremylt         ierr = KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]);
607ccaff030SJeremy L Thompson         CHKERRQ(ierr);
608ccaff030SJeremy L Thompson 
609ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
610d1d35e2fSjeremylt         ierr = KSPGetPC(ksp_smoother, &pc_smoother); CHKERRQ(ierr);
611d1d35e2fSjeremylt         ierr = PCSetType(pc_smoother, PCJACOBI); CHKERRQ(ierr);
612d1d35e2fSjeremylt         ierr = PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
613ccaff030SJeremy L Thompson 
614ccaff030SJeremy L Thompson         // -------- Work vector
615d1d35e2fSjeremylt         if (level != fine_level) {
616d1d35e2fSjeremylt           ierr = PCMGSetX(pc, level, U_g[level]); CHKERRQ(ierr);
617ccaff030SJeremy L Thompson         }
618ccaff030SJeremy L Thompson 
619ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
620ccaff030SJeremy L Thompson         if (level > 0) {
621d1d35e2fSjeremylt           ierr = PCMGSetInterpolation(pc, level, prolong_restr_mat[level]);
622ccaff030SJeremy L Thompson           CHKERRQ(ierr);
623d1d35e2fSjeremylt           ierr = PCMGSetRestriction(pc, level, prolong_restr_mat[level]);
624ccaff030SJeremy L Thompson           CHKERRQ(ierr);
625ccaff030SJeremy L Thompson         }
626ccaff030SJeremy L Thompson       }
627ccaff030SJeremy L Thompson 
628ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
629d1d35e2fSjeremylt       KSP ksp_coarse;
630d1d35e2fSjeremylt       PC pc_coarse;
631ccaff030SJeremy L Thompson 
632ccaff030SJeremy L Thompson       // -------- Coarse KSP
633d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
634d1d35e2fSjeremylt       ierr = KSPSetType(ksp_coarse, KSPPREONLY); CHKERRQ(ierr);
635d1d35e2fSjeremylt       ierr = KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse);
636ccaff030SJeremy L Thompson       CHKERRQ(ierr);
637d1d35e2fSjeremylt       ierr = KSPSetOptionsPrefix(ksp_coarse, "coarse_"); CHKERRQ(ierr);
638ccaff030SJeremy L Thompson 
639ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
640d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
641d1d35e2fSjeremylt       ierr = PCSetType(pc_coarse, PCGAMG); CHKERRQ(ierr);
642d1d35e2fSjeremylt       ierr = PCSetOptionsPrefix(pc_coarse, "coarse_"); CHKERRQ(ierr);
643ccaff030SJeremy L Thompson 
644d1d35e2fSjeremylt       ierr = KSPSetFromOptions(ksp_coarse); CHKERRQ(ierr);
645d1d35e2fSjeremylt       ierr = PCSetFromOptions(pc_coarse); CHKERRQ(ierr);
646ccaff030SJeremy L Thompson 
647ccaff030SJeremy L Thompson       // ------ PCMG options
648ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
649ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
650d1d35e2fSjeremylt       ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
651ccaff030SJeremy L Thompson     }
652ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
653ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
654ccaff030SJeremy L Thompson   }
655038d0cd7Sjeremylt   {
656038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
657d1d35e2fSjeremylt     SNESLineSearch line_search;
658481a4840SJed Brown 
659d1d35e2fSjeremylt     ierr = SNESGetLineSearch(snes, &line_search); CHKERRQ(ierr);
660d1d35e2fSjeremylt     ierr = SNESLineSearchSetType(line_search, SNESLINESEARCHCP); CHKERRQ(ierr);
661481a4840SJed Brown   }
662481a4840SJed Brown 
663ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
664ccaff030SJeremy L Thompson 
665ccaff030SJeremy L Thompson   // Performance logging
666ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
667ccaff030SJeremy L Thompson 
668ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
669ccaff030SJeremy L Thompson   // Set initial guess
670ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
671f81c27eaSJed Brown   ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr);
672ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
673ccaff030SJeremy L Thompson 
674ccaff030SJeremy L Thompson   // View solution
675d1d35e2fSjeremylt   if (app_ctx->view_soln) {
676d1d35e2fSjeremylt     ierr = ViewSolution(comm, app_ctx, U, 0, 0.0); CHKERRQ(ierr);
677ccaff030SJeremy L Thompson   }
678ccaff030SJeremy L Thompson 
679ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
680ccaff030SJeremy L Thompson   // Solve SNES
681ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
682d1d35e2fSjeremylt   PetscBool snes_monitor = PETSC_FALSE;
683d1d35e2fSjeremylt   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor);
6845f24f028Sjeremylt   CHKERRQ(ierr);
6855f24f028Sjeremylt 
686ccaff030SJeremy L Thompson   // Performance logging
687d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve);
688ccaff030SJeremy L Thompson   CHKERRQ(ierr);
689d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_solve); CHKERRQ(ierr);
690ccaff030SJeremy L Thompson 
691ccaff030SJeremy L Thompson   // Timing
692ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
693d1d35e2fSjeremylt   start_time = MPI_Wtime();
694ccaff030SJeremy L Thompson 
695ccaff030SJeremy L Thompson   // Solve for each load increment
6965f24f028Sjeremylt   PetscInt increment;
697d1d35e2fSjeremylt   for (increment = 1; increment <= app_ctx->num_increments; increment++) {
6985f24f028Sjeremylt     // -- Log increment count
699d1d35e2fSjeremylt     if (snes_monitor) {
7005f24f028Sjeremylt       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
7015f24f028Sjeremylt       CHKERRQ(ierr);
7025f24f028Sjeremylt     }
7035f24f028Sjeremylt 
704ccaff030SJeremy L Thompson     // -- Scale the problem
705d1d35e2fSjeremylt     PetscScalar load_increment = 1.0*increment / app_ctx->num_increments,
706d1d35e2fSjeremylt                 scalingFactor = load_increment /
707d1d35e2fSjeremylt                                 (increment == 1 ? 1 : res_ctx->load_increment);
708d1d35e2fSjeremylt     res_ctx->load_increment = load_increment;
709d1d35e2fSjeremylt     if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) {
710ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
711ccaff030SJeremy L Thompson     }
712ccaff030SJeremy L Thompson 
713ccaff030SJeremy L Thompson     // -- Solve
714ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
715ccaff030SJeremy L Thompson 
716ccaff030SJeremy L Thompson     // -- View solution
717d1d35e2fSjeremylt     if (app_ctx->view_soln) {
718d1d35e2fSjeremylt       ierr = ViewSolution(comm, app_ctx, U, increment, load_increment); CHKERRQ(ierr);
719ccaff030SJeremy L Thompson     }
720ccaff030SJeremy L Thompson 
721ccaff030SJeremy L Thompson     // -- Update SNES iteration count
722ccaff030SJeremy L Thompson     PetscInt its;
723ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
724bf0f51feSjeremylt     snes_its += its;
7257418ca88SJeremy L Thompson     ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr);
726bf0f51feSjeremylt     ksp_its += its;
7273951731eSjeremylt 
7283951731eSjeremylt     // -- Check for divergence
7293951731eSjeremylt     SNESConvergedReason reason;
7303951731eSjeremylt     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
7313951731eSjeremylt     if (reason < 0)
7323951731eSjeremylt       break;
7338355605fSKaren (Ren) Stengel     if (app_ctx->energy_viewer) {
7348355605fSKaren (Ren) Stengel       // -- Log strain energy for current load increment
7358355605fSKaren (Ren) Stengel       CeedScalar energy;
7368355605fSKaren (Ren) Stengel       ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
7378355605fSKaren (Ren) Stengel                                  U, &energy); CHKERRQ(ierr);
7388355605fSKaren (Ren) Stengel 
7398355605fSKaren (Ren) Stengel       if (!app_ctx->test_mode) {
7408355605fSKaren (Ren) Stengel         // -- Output
7418355605fSKaren (Ren) Stengel         ierr = PetscPrintf(comm,
7428355605fSKaren (Ren) Stengel                            "    Strain Energy                      : %.12e\n",
7438355605fSKaren (Ren) Stengel                            energy); CHKERRQ(ierr);
7448355605fSKaren (Ren) Stengel       }
7458355605fSKaren (Ren) Stengel       ierr = PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment,
7468355605fSKaren (Ren) Stengel                                     energy); CHKERRQ(ierr);
7478355605fSKaren (Ren) Stengel     }
748ccaff030SJeremy L Thompson   }
749ccaff030SJeremy L Thompson 
750ccaff030SJeremy L Thompson   // Timing
751d1d35e2fSjeremylt   elapsed_time = MPI_Wtime() - start_time;
752ccaff030SJeremy L Thompson 
753ccaff030SJeremy L Thompson   // Performance logging
754ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
755ccaff030SJeremy L Thompson 
756ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
757ccaff030SJeremy L Thompson   // Output summary
758ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
759d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
760ccaff030SJeremy L Thompson     // -- SNES
761d1d35e2fSjeremylt     SNESType snes_type;
762ccaff030SJeremy L Thompson     SNESConvergedReason reason;
763ccaff030SJeremy L Thompson     PetscReal rnorm;
764d1d35e2fSjeremylt     ierr = SNESGetType(snes, &snes_type); CHKERRQ(ierr);
765ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
767ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
768ccaff030SJeremy L Thompson                        "  SNES:\n"
769ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
770ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
771ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
7725f24f028Sjeremylt                        "    Completed Load Increments          : %d\n"
773ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
774ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
775d1d35e2fSjeremylt                        snes_type, SNESConvergedReasons[reason],
776d1d35e2fSjeremylt                        app_ctx->num_increments, increment - 1,
777bf0f51feSjeremylt                        snes_its, (double)rnorm); CHKERRQ(ierr);
778ccaff030SJeremy L Thompson 
779ccaff030SJeremy L Thompson     // -- KSP
780ccaff030SJeremy L Thompson     KSP ksp;
781d1d35e2fSjeremylt     KSPType ksp_type;
782ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
783d1d35e2fSjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
784ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
785ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
7867418ca88SJeremy L Thompson                        "    KSP Type                           : %s\n"
7877418ca88SJeremy L Thompson                        "    Total KSP Iterations               : %D\n",
788bf0f51feSjeremylt                        ksp_type, ksp_its); CHKERRQ(ierr);
789ccaff030SJeremy L Thompson 
790ccaff030SJeremy L Thompson     // -- PC
791ccaff030SJeremy L Thompson     PC pc;
792d1d35e2fSjeremylt     PCType pc_type;
793ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
794d1d35e2fSjeremylt     ierr = PCGetType(pc, &pc_type); CHKERRQ(ierr);
795e3e3df41Sjeremylt     ierr = PetscPrintf(comm,
796e3e3df41Sjeremylt                        "    PC Type                            : %s\n",
797d1d35e2fSjeremylt                        pc_type); CHKERRQ(ierr);
798e3e3df41Sjeremylt 
799d1d35e2fSjeremylt     if (!strcmp(pc_type, PCMG)) {
800d1d35e2fSjeremylt       PCMGType pcmg_type;
801d1d35e2fSjeremylt       ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
802ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
803ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
804ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
805ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
806d1d35e2fSjeremylt                          PCMGTypes[pcmg_type],
807d1d35e2fSjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
808ccaff030SJeremy L Thompson 
809ccaff030SJeremy L Thompson       // -- Coarse Solve
810d1d35e2fSjeremylt       KSP ksp_coarse;
811d1d35e2fSjeremylt       PC pc_coarse;
812d1d35e2fSjeremylt       PCType pc_type;
813ccaff030SJeremy L Thompson 
814d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
815d1d35e2fSjeremylt       ierr = KSPGetType(ksp_coarse, &ksp_type); CHKERRQ(ierr);
816d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
817d1d35e2fSjeremylt       ierr = PCGetType(pc_coarse, &pc_type); CHKERRQ(ierr);
818ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
819ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
820ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
821ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
822d1d35e2fSjeremylt                          ksp_type, pc_type); CHKERRQ(ierr);
823ccaff030SJeremy L Thompson     }
824ccaff030SJeremy L Thompson   }
825ccaff030SJeremy L Thompson 
826ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
827ccaff030SJeremy L Thompson   // Compute solve time
828ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
829d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
830d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm);
831ccaff030SJeremy L Thompson     CHKERRQ(ierr);
832d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm);
833ccaff030SJeremy L Thompson     CHKERRQ(ierr);
834ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
835ccaff030SJeremy L Thompson                        "  Performance:\n"
8367418ca88SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n"
8377418ca88SJeremy L Thompson                        "    DoFs/Sec in SNES                   : %g (%g) million\n",
838bf0f51feSjeremylt                        max_time, min_time, 1e-6*U_g_size[fine_level]*ksp_its/max_time,
839bf0f51feSjeremylt                        1e-6*U_g_size[fine_level]*ksp_its/min_time); CHKERRQ(ierr);
840ccaff030SJeremy L Thompson   }
841ccaff030SJeremy L Thompson 
842ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
843ccaff030SJeremy L Thompson   // Compute error
844ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
845d1d35e2fSjeremylt   if (app_ctx->forcing_choice == FORCE_MMS) {
846d1d35e2fSjeremylt     CeedScalar l2_error = 1., l2_U_norm = 1.;
847d1d35e2fSjeremylt     const CeedScalar *true_array;
848d1d35e2fSjeremylt     Vec error_vec, true_vec;
849ccaff030SJeremy L Thompson 
850ccaff030SJeremy L Thompson     // -- Work vectors
851d1d35e2fSjeremylt     ierr = VecDuplicate(U, &error_vec); CHKERRQ(ierr);
852d1d35e2fSjeremylt     ierr = VecSet(error_vec, 0.0); CHKERRQ(ierr);
853d1d35e2fSjeremylt     ierr = VecDuplicate(U, &true_vec); CHKERRQ(ierr);
854d1d35e2fSjeremylt     ierr = VecSet(true_vec, 0.0); CHKERRQ(ierr);
855ccaff030SJeremy L Thompson 
856ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
857d1d35e2fSjeremylt     CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln,
858d1d35e2fSjeremylt                            CEED_MEM_HOST, &true_array);
859d1d35e2fSjeremylt     ierr = VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array);
86062e9c006SJeremy L Thompson     CHKERRQ(ierr);
861d1d35e2fSjeremylt     ierr = DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec);
862ccaff030SJeremy L Thompson     CHKERRQ(ierr);
863d1d35e2fSjeremylt     ierr = VecResetArray(res_ctx->Y_loc); CHKERRQ(ierr);
864d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array);
865ccaff030SJeremy L Thompson 
866ccaff030SJeremy L Thompson     // -- Compute L2 error
867d1d35e2fSjeremylt     ierr = VecWAXPY(error_vec, -1.0, U, true_vec); CHKERRQ(ierr);
868d1d35e2fSjeremylt     ierr = VecNorm(error_vec, NORM_2, &l2_error); CHKERRQ(ierr);
869d1d35e2fSjeremylt     ierr = VecNorm(U, NORM_2, &l2_U_norm); CHKERRQ(ierr);
870d1d35e2fSjeremylt     l2_error /= l2_U_norm;
871ccaff030SJeremy L Thompson 
872ccaff030SJeremy L Thompson     // -- Output
873d1d35e2fSjeremylt     if (!app_ctx->test_mode || l2_error > 0.05) {
874ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
875ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
876d1d35e2fSjeremylt                          l2_error); CHKERRQ(ierr);
877ccaff030SJeremy L Thompson     }
878ccaff030SJeremy L Thompson 
879ccaff030SJeremy L Thompson     // -- Cleanup
880d1d35e2fSjeremylt     ierr = VecDestroy(&error_vec); CHKERRQ(ierr);
881d1d35e2fSjeremylt     ierr = VecDestroy(&true_vec); CHKERRQ(ierr);
8822d93065eSjeremylt   }
8832d93065eSjeremylt 
8842d93065eSjeremylt   // ---------------------------------------------------------------------------
8852d93065eSjeremylt   // Compute energy
8862d93065eSjeremylt   // ---------------------------------------------------------------------------
8878355605fSKaren (Ren) Stengel   PetscReal energy;
888d1d35e2fSjeremylt   ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
889a3c02c40SJeremy L Thompson                              U, &energy); CHKERRQ(ierr);
8908355605fSKaren (Ren) Stengel   if (!app_ctx->test_mode) {
8912d93065eSjeremylt     // -- Output
8922d93065eSjeremylt     ierr = PetscPrintf(comm,
89372d03b64SArash Mehraban                        "    Strain Energy                      : %.12e\n",
8942d93065eSjeremylt                        energy); CHKERRQ(ierr);
895ccaff030SJeremy L Thompson   }
8968355605fSKaren (Ren) Stengel   ierr = RegressionTests_solids(app_ctx, energy); CHKERRQ(ierr);
897ccaff030SJeremy L Thompson 
898ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
8995c25879aSJeremy L Thompson   // Output diagnostic quantities
9005c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
901d1d35e2fSjeremylt   if (app_ctx->view_soln || app_ctx->view_final_soln) {
9025c25879aSJeremy L Thompson     // -- Setup context
903d1d35e2fSjeremylt     UserMult diagnostic_ctx;
904d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &diagnostic_ctx); CHKERRQ(ierr);
905d1d35e2fSjeremylt     ierr = PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)); CHKERRQ(ierr);
906d1d35e2fSjeremylt     diagnostic_ctx->dm = dm_diagnostic;
907d1d35e2fSjeremylt     diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic;
9085c25879aSJeremy L Thompson 
9095c25879aSJeremy L Thompson     // -- Compute and output
910d1d35e2fSjeremylt     ierr = ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx,
911d1d35e2fSjeremylt                                     app_ctx, U,
912d1d35e2fSjeremylt                                     ceed_data[fine_level]->elem_restr_diagnostic);
9135c25879aSJeremy L Thompson     CHKERRQ(ierr);
9145c25879aSJeremy L Thompson 
9155c25879aSJeremy L Thompson     // -- Cleanup
916d1d35e2fSjeremylt     ierr = PetscFree(diagnostic_ctx); CHKERRQ(ierr);
9175c25879aSJeremy L Thompson   }
9185c25879aSJeremy L Thompson 
9195c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
920ccaff030SJeremy L Thompson   // Free objects
921ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
922ccaff030SJeremy L Thompson   // Data in arrays per level
923d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
924ccaff030SJeremy L Thompson     // Vectors
925d1d35e2fSjeremylt     ierr = VecDestroy(&U_g[level]); CHKERRQ(ierr);
926d1d35e2fSjeremylt     ierr = VecDestroy(&U_loc[level]); CHKERRQ(ierr);
927ccaff030SJeremy L Thompson 
928ccaff030SJeremy L Thompson     // Jacobian matrix and data
929d1d35e2fSjeremylt     ierr = VecDestroy(&jacob_ctx[level]->Y_loc); CHKERRQ(ierr);
930d1d35e2fSjeremylt     ierr = MatDestroy(&jacob_mat[level]); CHKERRQ(ierr);
931d1d35e2fSjeremylt     ierr = PetscFree(jacob_ctx[level]); CHKERRQ(ierr);
932ccaff030SJeremy L Thompson 
933ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
934ccaff030SJeremy L Thompson     if (level > 0) {
935d1d35e2fSjeremylt       ierr = PetscFree(prolong_restr_ctx[level]); CHKERRQ(ierr);
936d1d35e2fSjeremylt       ierr = MatDestroy(&prolong_restr_mat[level]); CHKERRQ(ierr);
937ccaff030SJeremy L Thompson     }
938ccaff030SJeremy L Thompson 
939ccaff030SJeremy L Thompson     // DM
940d1d35e2fSjeremylt     ierr = DMDestroy(&level_dms[level]); CHKERRQ(ierr);
941ccaff030SJeremy L Thompson 
942ccaff030SJeremy L Thompson     // libCEED objects
943d1d35e2fSjeremylt     ierr = CeedDataDestroy(level, ceed_data[level]); CHKERRQ(ierr);
944ccaff030SJeremy L Thompson   }
945ccaff030SJeremy L Thompson 
9468355605fSKaren (Ren) Stengel   ierr = PetscViewerDestroy(&app_ctx->energy_viewer); CHKERRQ(ierr);
9478355605fSKaren (Ren) Stengel 
948ccaff030SJeremy L Thompson   // Arrays
949d1d35e2fSjeremylt   ierr = PetscFree(U_g); CHKERRQ(ierr);
950d1d35e2fSjeremylt   ierr = PetscFree(U_loc); CHKERRQ(ierr);
951d1d35e2fSjeremylt   ierr = PetscFree(U_g_size); CHKERRQ(ierr);
952d1d35e2fSjeremylt   ierr = PetscFree(U_l_size); CHKERRQ(ierr);
953d1d35e2fSjeremylt   ierr = PetscFree(U_loc_size); CHKERRQ(ierr);
954d1d35e2fSjeremylt   ierr = PetscFree(jacob_ctx); CHKERRQ(ierr);
955d1d35e2fSjeremylt   ierr = PetscFree(jacob_mat); CHKERRQ(ierr);
956d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_ctx); CHKERRQ(ierr);
957d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_mat); CHKERRQ(ierr);
958d1d35e2fSjeremylt   ierr = PetscFree(app_ctx->level_degrees); CHKERRQ(ierr);
959d1d35e2fSjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
960ccaff030SJeremy L Thompson 
961ccaff030SJeremy L Thompson   // libCEED objects
96217f0843fSJeremy L Thompson   CeedVectorDestroy(&form_jacob_ctx->coo_values);
963d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys);
964d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys_smoother);
965ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
966ccaff030SJeremy L Thompson 
967ccaff030SJeremy L Thompson   // PETSc objects
968ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
969ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
970d1d35e2fSjeremylt   ierr = VecDestroy(&R_loc); CHKERRQ(ierr);
971ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
972d1d35e2fSjeremylt   ierr = VecDestroy(&F_loc); CHKERRQ(ierr);
973d1d35e2fSjeremylt   ierr = VecDestroy(&neumann_bcs); CHKERRQ(ierr);
974d1d35e2fSjeremylt   ierr = VecDestroy(&bcs_loc); CHKERRQ(ierr);
975d1d35e2fSjeremylt   ierr = MatDestroy(&jacob_mat_coarse); CHKERRQ(ierr);
976ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
977d1d35e2fSjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
978d1d35e2fSjeremylt   ierr = DMDestroy(&dm_energy); CHKERRQ(ierr);
979d1d35e2fSjeremylt   ierr = DMDestroy(&dm_diagnostic); CHKERRQ(ierr);
980d1d35e2fSjeremylt   ierr = PetscFree(level_dms); CHKERRQ(ierr);
981ccaff030SJeremy L Thompson 
9825754ecacSJeremy L Thompson   // -- Function list
9835754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupPhysics);
9845754ecacSJeremy L Thompson   CHKERRQ(ierr);
9855754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel);
9865754ecacSJeremy L Thompson   CHKERRQ(ierr);
9875754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedLevel);
9885754ecacSJeremy L Thompson   CHKERRQ(ierr);
9895754ecacSJeremy L Thompson 
990ccaff030SJeremy L Thompson   // Structs
991d1d35e2fSjeremylt   ierr = PetscFree(res_ctx); CHKERRQ(ierr);
992d1d35e2fSjeremylt   ierr = PetscFree(form_jacob_ctx); CHKERRQ(ierr);
993d1d35e2fSjeremylt   ierr = PetscFree(jacob_coarse_ctx); CHKERRQ(ierr);
994d1d35e2fSjeremylt   ierr = PetscFree(app_ctx); CHKERRQ(ierr);
9955754ecacSJeremy L Thompson   ierr = PetscFree(problem_functions); CHKERRQ(ierr);
996ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
997ccaff030SJeremy L Thompson 
998ccaff030SJeremy L Thompson   return PetscFinalize();
999ccaff030SJeremy L Thompson }
1000