xref: /libCEED/examples/solids/elasticity.c (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ccaff030SJeremy L Thompson //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ccaff030SJeremy L Thompson //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7ccaff030SJeremy L Thompson 
8ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Elasticity
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve
11ccaff030SJeremy L Thompson //   elasticity problems.
12ccaff030SJeremy L Thompson //
13ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex.
14ccaff030SJeremy L Thompson //
15ccaff030SJeremy L Thompson // Build with:
16ccaff030SJeremy L Thompson //
17ccaff030SJeremy L Thompson //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18ccaff030SJeremy L Thompson //
19ccaff030SJeremy L Thompson // Sample runs:
20ccaff030SJeremy L Thompson //
21eccc2849SRezgar Shakeri //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms
22eccc2849SRezgar 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
23eccc2849SRezgar 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
24ccaff030SJeremy L Thompson //
25ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
26ccaff030SJeremy L Thompson //
278355605fSKaren (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
288355605fSKaren (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
298355605fSKaren (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
30ccaff030SJeremy L Thompson 
31ccaff030SJeremy L Thompson /// @file
32ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
33ccaff030SJeremy L Thompson 
34ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson #include "elasticity.h"
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson int main(int argc, char **argv) {
39ccaff030SJeremy L Thompson   PetscInt       ierr;
40ccaff030SJeremy L Thompson   MPI_Comm       comm;
41ccaff030SJeremy L Thompson   // Context structs
42d1d35e2fSjeremylt   AppCtx         app_ctx;                  // Contains problem options
435754ecacSJeremy L Thompson   ProblemFunctions problem_functions;      // Setup functions for each problem
44ccaff030SJeremy L Thompson   Units          units;                    // Contains units scaling
45ccaff030SJeremy L Thompson   // PETSc objects
46d1d35e2fSjeremylt   PetscLogStage  stage_dm_setup, stage_libceed_setup,
47d1d35e2fSjeremylt                  stage_snes_setup, stage_snes_solve;
48d1d35e2fSjeremylt   DM             dm_orig;                  // Distributed DM to clone
49d1d35e2fSjeremylt   DM             dm_energy, dm_diagnostic; // DMs for postprocessing
50d1d35e2fSjeremylt   DM             *level_dms;
51d1d35e2fSjeremylt   Vec            U, *U_g, *U_loc;          // U: solution, R: residual, F: forcing
52d1d35e2fSjeremylt   Vec            R, R_loc, F, F_loc;       // g: global, loc: local
53d1d35e2fSjeremylt   Vec            neumann_bcs = NULL, bcs_loc = NULL;
5417f0843fSJeremy L Thompson   SNES           snes;
55d1d35e2fSjeremylt   Mat            *jacob_mat, jacob_mat_coarse, *prolong_restr_mat;
56ccaff030SJeremy L Thompson   // PETSc data
57d1d35e2fSjeremylt   UserMult       res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx;
58d1d35e2fSjeremylt   FormJacobCtx   form_jacob_ctx;
59d1d35e2fSjeremylt   UserMultProlongRestr *prolong_restr_ctx;
60d1d35e2fSjeremylt   PCMGCycleType  pcmg_cycle_type = PC_MG_CYCLE_V;
61ccaff030SJeremy L Thompson   // libCEED objects
62d99fa3c5SJeremy L Thompson   Ceed           ceed;
63d1d35e2fSjeremylt   CeedData       *ceed_data;
64d1d35e2fSjeremylt   CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL;
65ccaff030SJeremy L Thompson   // Parameters
66d1d35e2fSjeremylt   PetscInt       num_comp_u = 3;                 // 3 DoFs in 3D
67d1d35e2fSjeremylt   PetscInt       num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic
68d1d35e2fSjeremylt   PetscInt       num_levels = 1, fine_level = 0;
69bf0f51feSjeremylt   PetscInt       *U_g_size, *U_l_size, *U_loc_size;
70bf0f51feSjeremylt   PetscInt       snes_its = 0, ksp_its = 0;
71d1d35e2fSjeremylt   double         start_time, elapsed_time, min_time, max_time;
72ccaff030SJeremy L Thompson 
73ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
74bf0f51feSjeremylt   if (ierr) return ierr;
75ccaff030SJeremy L Thompson 
76ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
77ccaff030SJeremy L Thompson   // Process command line options
78ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
79ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
80ccaff030SJeremy L Thompson 
81ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
82d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr);
83d1d35e2fSjeremylt   ierr = ProcessCommandLineOptions(comm, app_ctx); CHKERRQ(ierr);
845754ecacSJeremy L Thompson   ierr = PetscCalloc1(1, &problem_functions); CHKERRQ(ierr);
855754ecacSJeremy L Thompson   ierr = RegisterProblems(problem_functions); CHKERRQ(ierr);
86d1d35e2fSjeremylt   num_levels = app_ctx->num_levels;
87d1d35e2fSjeremylt   fine_level = num_levels - 1;
88ccaff030SJeremy L Thompson 
89ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
90d1d35e2fSjeremylt   // Initialize libCEED
9162e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
92d1d35e2fSjeremylt   // Initialize backend
93d1d35e2fSjeremylt   CeedInit(app_ctx->ceed_resource, &ceed);
9462e9c006SJeremy L Thompson 
9562e9c006SJeremy L Thompson   // Check preferred MemType
96d1d35e2fSjeremylt   CeedMemType mem_type_backend;
97d1d35e2fSjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
985754ecacSJeremy L Thompson   // Setup physics context and wrap in libCEED object
995754ecacSJeremy L Thompson   {
1005754ecacSJeremy L Thompson     PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *);
1015754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name,
1025754ecacSJeremy L Thompson                                  &SetupPhysics); CHKERRQ(ierr);
1035754ecacSJeremy L Thompson     if (!SetupPhysics)
1047578c821SJeremy L Thompson       SETERRQ(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found",
1055754ecacSJeremy L Thompson               app_ctx->name);
1065754ecacSJeremy L Thompson     ierr = (*SetupPhysics)(comm, ceed, &units, &ctx_phys); CHKERRQ(ierr);
1075754ecacSJeremy L Thompson     PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext,
1085754ecacSJeremy L Thompson                                            CeedQFunctionContext *);
1095754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupSmootherPhysics,
1105754ecacSJeremy L Thompson                                  app_ctx->name, &SetupSmootherPhysics);
1115754ecacSJeremy L Thompson     CHKERRQ(ierr);
1125754ecacSJeremy L Thompson     if (!SetupSmootherPhysics)
1137578c821SJeremy L Thompson       SETERRQ(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found",
1145754ecacSJeremy L Thompson               app_ctx->name);
1155754ecacSJeremy L Thompson     ierr = (*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother);
1165754ecacSJeremy L Thompson     CHKERRQ(ierr);
117777ff853SJeremy L Thompson   }
118777ff853SJeremy L Thompson 
11962e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
120ccaff030SJeremy L Thompson   // Setup DM
121ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
122ccaff030SJeremy L Thompson   // Performance logging
123d1d35e2fSjeremylt   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup);
124ccaff030SJeremy L Thompson   CHKERRQ(ierr);
125d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_dm_setup); CHKERRQ(ierr);
126ccaff030SJeremy L Thompson 
127ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
128d1d35e2fSjeremylt   ierr = CreateDistributedDM(comm, app_ctx, &dm_orig); CHKERRQ(ierr);
129b68a8d79SJed Brown   VecType vectype;
130d1d35e2fSjeremylt   switch (mem_type_backend) {
131b68a8d79SJed Brown   case CEED_MEM_HOST: vectype = VECSTANDARD; break;
132b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
133b68a8d79SJed Brown     const char *resolved;
134b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
135b68a8d79SJed Brown     if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA;
136b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP;
137b68a8d79SJed Brown     else vectype = VECSTANDARD;
138b68a8d79SJed Brown   }
139b68a8d79SJed Brown   }
140d1d35e2fSjeremylt   ierr = DMSetVecType(dm_orig, vectype); CHKERRQ(ierr);
1413fc8a154SJed Brown   ierr = DMPlexDistributeSetDefault(dm_orig, PETSC_FALSE); CHKERRQ(ierr);
142d1d35e2fSjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
143ccaff030SJeremy L Thompson 
144ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
145d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &level_dms); CHKERRQ(ierr);
146d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
147d1d35e2fSjeremylt     ierr = DMClone(dm_orig, &level_dms[level]); CHKERRQ(ierr);
148d1d35e2fSjeremylt     ierr = DMGetVecType(dm_orig, &vectype); CHKERRQ(ierr);
149d1d35e2fSjeremylt     ierr = DMSetVecType(level_dms[level], vectype); CHKERRQ(ierr);
150d1d35e2fSjeremylt     ierr = SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level],
151d1d35e2fSjeremylt                            PETSC_TRUE, num_comp_u); CHKERRQ(ierr);
152a3c02c40SJeremy L Thompson     // -- Label field components for viewing
153a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
154a3c02c40SJeremy L Thompson     PetscSection section;
155d1d35e2fSjeremylt     ierr = DMGetLocalSection(level_dms[level], &section); CHKERRQ(ierr);
156a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr);
157a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
158a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
159a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
160a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
161a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
162a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
163a3c02c40SJeremy L Thompson   }
164a3c02c40SJeremy L Thompson 
165a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
166d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_energy); CHKERRQ(ierr);
167d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level],
168d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_e); CHKERRQ(ierr);
169d1d35e2fSjeremylt   ierr = DMClone(dm_orig, &dm_diagnostic); CHKERRQ(ierr);
170d1d35e2fSjeremylt   ierr = SetupDMByDegree(dm_diagnostic, app_ctx,
171d1d35e2fSjeremylt                          app_ctx->level_degrees[fine_level],
172d1d35e2fSjeremylt                          PETSC_FALSE, num_comp_u + num_comp_d); CHKERRQ(ierr);
173d1d35e2fSjeremylt   ierr = DMSetVecType(dm_energy, vectype); CHKERRQ(ierr);
174d1d35e2fSjeremylt   ierr = DMSetVecType(dm_diagnostic, vectype); CHKERRQ(ierr);
175a3c02c40SJeremy L Thompson   {
176a3c02c40SJeremy L Thompson     // -- Label field components for viewing
177a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
178a3c02c40SJeremy L Thompson     PetscSection section;
179d1d35e2fSjeremylt     ierr = DMGetLocalSection(dm_diagnostic, &section); CHKERRQ(ierr);
180a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr);
1818ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
1828ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1838ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
1848ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1858ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
1868ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
187379387d4SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure");
1888ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1897ab18febSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain");
19013fdad4bSJeremy L Thompson     CHKERRQ(ierr);
19113fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2");
19213fdad4bSJeremy L Thompson     CHKERRQ(ierr);
19313fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 6, "detJ");
19413fdad4bSJeremy L Thompson     CHKERRQ(ierr);
19513fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity");
196a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
197ccaff030SJeremy L Thompson   }
198ccaff030SJeremy L Thompson 
199ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
200ccaff030SJeremy L Thompson   // Setup solution and work vectors
201ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
202ccaff030SJeremy L Thompson   // Allocate arrays
203d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g); CHKERRQ(ierr);
204d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc); CHKERRQ(ierr);
205d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_g_size); CHKERRQ(ierr);
206d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_l_size); CHKERRQ(ierr);
207d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &U_loc_size); CHKERRQ(ierr);
208ccaff030SJeremy L Thompson 
209ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
210d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
211ccaff030SJeremy L Thompson     // -- Create global unknown vector U
212d1d35e2fSjeremylt     ierr = DMCreateGlobalVector(level_dms[level], &U_g[level]); CHKERRQ(ierr);
213d1d35e2fSjeremylt     ierr = VecGetSize(U_g[level], &U_g_size[level]); CHKERRQ(ierr);
214ccaff030SJeremy L Thompson     // Note: Local size for matShell
215d1d35e2fSjeremylt     ierr = VecGetLocalSize(U_g[level], &U_l_size[level]); CHKERRQ(ierr);
216ccaff030SJeremy L Thompson 
217d1d35e2fSjeremylt     // -- Create local unknown vector U_loc
218d1d35e2fSjeremylt     ierr = DMCreateLocalVector(level_dms[level], &U_loc[level]); CHKERRQ(ierr);
219ccaff030SJeremy L Thompson     // Note: local size for libCEED
220d1d35e2fSjeremylt     ierr = VecGetSize(U_loc[level], &U_loc_size[level]); CHKERRQ(ierr);
221ccaff030SJeremy L Thompson   }
222ccaff030SJeremy L Thompson 
223ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
224d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &U); CHKERRQ(ierr);
225d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &R); CHKERRQ(ierr);
226d1d35e2fSjeremylt   ierr = VecDuplicate(U_g[fine_level], &F); CHKERRQ(ierr);
227d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &R_loc); CHKERRQ(ierr);
228d1d35e2fSjeremylt   ierr = VecDuplicate(U_loc[fine_level], &F_loc); CHKERRQ(ierr);
229ccaff030SJeremy L Thompson 
230ccaff030SJeremy L Thompson   // Performance logging
231ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
232ccaff030SJeremy L Thompson 
233ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
234ccaff030SJeremy L Thompson   // Set up libCEED
235ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
236ccaff030SJeremy L Thompson   // Performance logging
237d1d35e2fSjeremylt   ierr = PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup);
238ccaff030SJeremy L Thompson   CHKERRQ(ierr);
239d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_libceed_setup); CHKERRQ(ierr);
240ccaff030SJeremy L Thompson 
241ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
242d1d35e2fSjeremylt   CeedVector force_ceed;
243ccaff030SJeremy L Thompson   CeedScalar *f;
244d1d35e2fSjeremylt   PetscMemType force_mem_type;
245d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
246d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(F_loc, &f, &force_mem_type); CHKERRQ(ierr);
247d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed);
248d1d35e2fSjeremylt     CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f);
249ccaff030SJeremy L Thompson   }
250ccaff030SJeremy L Thompson 
251fe394131Sjeremylt   // -- Create libCEED local Neumann BCs vector
252d1d35e2fSjeremylt   CeedVector neumann_ceed;
253fe394131Sjeremylt   CeedScalar *n;
254d1d35e2fSjeremylt   PetscMemType nummann_mem_type;
255d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
256d1d35e2fSjeremylt     ierr = VecDuplicate(U, &neumann_bcs); CHKERRQ(ierr);
257d1d35e2fSjeremylt     ierr = VecDuplicate(U_loc[fine_level], &bcs_loc); CHKERRQ(ierr);
258d1d35e2fSjeremylt     ierr = VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type); CHKERRQ(ierr);
259d1d35e2fSjeremylt     CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed);
260d1d35e2fSjeremylt     CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type),
261fe394131Sjeremylt                        CEED_USE_POINTER, n);
262fe394131Sjeremylt   }
263fe394131Sjeremylt 
264ccaff030SJeremy L Thompson   // -- Setup libCEED objects
265d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
266d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
267d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &ceed_data[fine_level]); CHKERRQ(ierr);
268ccaff030SJeremy L Thompson   {
2695754ecacSJeremy L Thompson     PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx,
2705754ecacSJeremy L Thompson                                             CeedQFunctionContext, PetscInt,
2715754ecacSJeremy L Thompson                                             PetscInt, PetscInt, PetscInt,
2725754ecacSJeremy L Thompson                                             CeedVector, CeedVector, CeedData *);
2735754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupLibceedFineLevel,
2745754ecacSJeremy L Thompson                                  app_ctx->name, &SetupLibceedFineLevel);
2755754ecacSJeremy L Thompson     CHKERRQ(ierr);
2765754ecacSJeremy L Thompson     if (!SetupLibceedFineLevel)
2777578c821SJeremy L Thompson       SETERRQ(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found",
2785754ecacSJeremy L Thompson               app_ctx->name);
2795754ecacSJeremy L Thompson     ierr = (*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic,
2805754ecacSJeremy L Thompson                                     ceed, app_ctx, ctx_phys, fine_level,
2815754ecacSJeremy L Thompson                                     num_comp_u, U_g_size[fine_level],
2825754ecacSJeremy L Thompson                                     U_loc_size[fine_level],
2835754ecacSJeremy L Thompson                                     force_ceed, neumann_ceed, ceed_data);
284ccaff030SJeremy L Thompson     CHKERRQ(ierr);
285ccaff030SJeremy L Thompson   }
286d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
287d1d35e2fSjeremylt   for (PetscInt level = num_levels - 2; level >= 0; level--) {
288d1d35e2fSjeremylt     ierr = PetscCalloc1(1, &ceed_data[level]); CHKERRQ(ierr);
289d99fa3c5SJeremy L Thompson 
290d99fa3c5SJeremy L Thompson     // Get global communication restriction
291d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
292d1d35e2fSjeremylt     ierr = VecSet(U_loc[level+1], 1.0); CHKERRQ(ierr);
293d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[level+1], U_loc[level+1], ADD_VALUES,
294d1d35e2fSjeremylt                            U_g[level+1]); CHKERRQ(ierr);
295d1d35e2fSjeremylt     ierr = DMGlobalToLocal(level_dms[level+1], U_g[level+1], INSERT_VALUES,
296d1d35e2fSjeremylt                            U_loc[level+1]); CHKERRQ(ierr);
297d99fa3c5SJeremy L Thompson 
298d99fa3c5SJeremy L Thompson     // Place in libCEED array
299d99fa3c5SJeremy L Thompson     const PetscScalar *m;
300d1d35e2fSjeremylt     PetscMemType m_mem_type;
301d1d35e2fSjeremylt     ierr = VecGetArrayReadAndMemType(U_loc[level+1], &m, &m_mem_type);
302d1d35e2fSjeremylt     CHKERRQ(ierr);
303d1d35e2fSjeremylt     CeedVectorSetArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
304d99fa3c5SJeremy L Thompson                        CEED_USE_POINTER, (CeedScalar *)m);
305ccaff030SJeremy L Thompson 
3061c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
3075754ecacSJeremy L Thompson     PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt,
3085754ecacSJeremy L Thompson                                         PetscInt, PetscInt, PetscInt, CeedVector, CeedData *);
3095754ecacSJeremy L Thompson     ierr = PetscFunctionListFind(problem_functions->setupLibceedLevel,
3105754ecacSJeremy L Thompson                                  app_ctx->name, &SetupLibceedLevel);
3115754ecacSJeremy L Thompson     CHKERRQ(ierr);
312caae6498SJed Brown     if (!SetupLibceedLevel)
3137578c821SJeremy L Thompson       SETERRQ(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found",
314caae6498SJed Brown               app_ctx->name);
3155754ecacSJeremy L Thompson     ierr = (*SetupLibceedLevel)(level_dms[level], ceed, app_ctx,
3165754ecacSJeremy L Thompson                                 level, num_comp_u, U_g_size[level],
3175754ecacSJeremy L Thompson                                 U_loc_size[level], ceed_data[level+1]->x_ceed,
3185754ecacSJeremy L Thompson                                 ceed_data);
319777ff853SJeremy L Thompson     CHKERRQ(ierr);
320d99fa3c5SJeremy L Thompson 
321d99fa3c5SJeremy L Thompson     // Restore PETSc vector
322d1d35e2fSjeremylt     CeedVectorTakeArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type),
323d99fa3c5SJeremy L Thompson                         (CeedScalar **)&m);
324d1d35e2fSjeremylt     ierr = VecRestoreArrayReadAndMemType(U_loc[level+1], &m); CHKERRQ(ierr);
325d1d35e2fSjeremylt     ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr);
326d1d35e2fSjeremylt     ierr = VecZeroEntries(U_loc[level+1]); CHKERRQ(ierr);
327ccaff030SJeremy L Thompson   }
328ccaff030SJeremy L Thompson 
329ccaff030SJeremy L Thompson   // Performance logging
330ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
331ccaff030SJeremy L Thompson 
332ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
333fe394131Sjeremylt   // Setup global forcing and Neumann BC vectors
334ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
335ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
336ccaff030SJeremy L Thompson 
337d1d35e2fSjeremylt   if (app_ctx->forcing_choice != FORCE_NONE) {
338d1d35e2fSjeremylt     CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL);
339d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(F_loc, &f); CHKERRQ(ierr);
340d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F);
341ccaff030SJeremy L Thompson     CHKERRQ(ierr);
342d1d35e2fSjeremylt     CeedVectorDestroy(&force_ceed);
343ccaff030SJeremy L Thompson   }
344ccaff030SJeremy L Thompson 
345d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0) {
346d1d35e2fSjeremylt     ierr = VecZeroEntries(neumann_bcs); CHKERRQ(ierr);
347d1d35e2fSjeremylt     CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL);
348d1d35e2fSjeremylt     ierr = VecRestoreArrayAndMemType(bcs_loc, &n); CHKERRQ(ierr);
349d1d35e2fSjeremylt     ierr = DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs);
350fe394131Sjeremylt     CHKERRQ(ierr);
351d1d35e2fSjeremylt     CeedVectorDestroy(&neumann_ceed);
352fe394131Sjeremylt   }
353fe394131Sjeremylt 
354ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
355ccaff030SJeremy L Thompson   // Print problem summary
356ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
357d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
358ccaff030SJeremy L Thompson     const char *usedresource;
359ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
36077841947SLeila Ghaffari     char hostname[PETSC_MAX_PATH_LEN];
36177841947SLeila Ghaffari     ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
36277841947SLeila Ghaffari     PetscInt comm_size;
36377841947SLeila Ghaffari     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson 
365ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
36617fd649aSJeremy L Thompson                        "\n-- Elasticity Example - libCEED + PETSc --\n"
36777841947SLeila Ghaffari                        "  MPI:\n"
36877841947SLeila Ghaffari                        "    Hostname                           : %s\n"
36977841947SLeila Ghaffari                        "    Total ranks                        : %d\n"
370ccaff030SJeremy L Thompson                        "  libCEED:\n"
37162e9c006SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
372b68a8d79SJed Brown                        "    libCEED Backend MemType            : %s\n",
37377841947SLeila Ghaffari                        hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend]);
37462e9c006SJeremy L Thompson     CHKERRQ(ierr);
375ccaff030SJeremy L Thompson 
37662e9c006SJeremy L Thompson     VecType vecType;
37762e9c006SJeremy L Thompson     ierr = VecGetType(U, &vecType); CHKERRQ(ierr);
37862e9c006SJeremy L Thompson     ierr = PetscPrintf(comm,
37962e9c006SJeremy L Thompson                        "  PETSc:\n"
38062e9c006SJeremy L Thompson                        "    PETSc Vec Type                     : %s\n",
38162e9c006SJeremy L Thompson                        vecType); CHKERRQ(ierr);
38262e9c006SJeremy L Thompson 
383ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
384ccaff030SJeremy L Thompson                        "  Problem:\n"
385ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
386ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
387ccaff030SJeremy L Thompson                        "  Mesh:\n"
388ccaff030SJeremy L Thompson                        "    File                               : %s\n"
389ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
390ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
391ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
392ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
393ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
394ccaff030SJeremy L Thompson                        "  Multigrid:\n"
395ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
396ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
3975754ecacSJeremy L Thompson                        app_ctx->name_for_disp,
398d1d35e2fSjeremylt                        forcing_types_for_disp[app_ctx->forcing_choice],
399d1d35e2fSjeremylt                        app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh",
400d1d35e2fSjeremylt                        app_ctx->degree + 1, app_ctx->degree + 1,
4015754ecacSJeremy L Thompson                        U_g_size[fine_level]/num_comp_u, U_l_size[fine_level]/num_comp_u,
4025754ecacSJeremy L Thompson                        num_comp_u,
403d1d35e2fSjeremylt                        (app_ctx->degree == 1 &&
404d1d35e2fSjeremylt                         app_ctx->multigrid_choice != MULTIGRID_NONE) ?
405f892e6d1Sjeremylt                        "Algebraic multigrid" :
406d1d35e2fSjeremylt                        multigrid_types_for_disp[app_ctx->multigrid_choice],
407d1d35e2fSjeremylt                        (app_ctx->degree == 1 ||
408d1d35e2fSjeremylt                         app_ctx->multigrid_choice == MULTIGRID_NONE) ?
409d1d35e2fSjeremylt                        0 : num_levels); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson 
411d1d35e2fSjeremylt     if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
412e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
413d1d35e2fSjeremylt         CeedInt level = i ? fine_level : 0;
414ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
415ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
416ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
417ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
418ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
419ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
420d1d35e2fSjeremylt                            app_ctx->level_degrees[level] + 1,
421d1d35e2fSjeremylt                            U_g_size[level]/num_comp_u, U_l_size[level]/num_comp_u);
422ccaff030SJeremy L Thompson         CHKERRQ(ierr);
423ccaff030SJeremy L Thompson       }
424ccaff030SJeremy L Thompson     }
425ccaff030SJeremy L Thompson   }
426ccaff030SJeremy L Thompson 
427ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
428ccaff030SJeremy L Thompson   // Setup SNES
429ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
430ccaff030SJeremy L Thompson   // Performance logging
431d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup);
432ccaff030SJeremy L Thompson   CHKERRQ(ierr);
433d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_setup); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson 
435ccaff030SJeremy L Thompson   // Create SNES
436ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
437d1d35e2fSjeremylt   ierr = SNESSetDM(snes, level_dms[fine_level]); CHKERRQ(ierr);
438ccaff030SJeremy L Thompson 
439ccaff030SJeremy L Thompson   // -- Jacobian evaluators
440d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_ctx); CHKERRQ(ierr);
441d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &jacob_mat); CHKERRQ(ierr);
442d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
443ccaff030SJeremy L Thompson     // -- Jacobian context for level
444d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &jacob_ctx[level]); CHKERRQ(ierr);
445d1d35e2fSjeremylt     ierr = SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level],
446d1d35e2fSjeremylt                             U_loc[level], ceed_data[level], ceed, ctx_phys,
447d1d35e2fSjeremylt                             ctx_phys_smoother, jacob_ctx[level]); CHKERRQ(ierr);
448ccaff030SJeremy L Thompson 
449ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
450d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level],
451d1d35e2fSjeremylt                           U_g_size[level], jacob_ctx[level], &jacob_mat[level]);
452ccaff030SJeremy L Thompson     CHKERRQ(ierr);
453d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_MULT,
454ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
455ccaff030SJeremy L Thompson     CHKERRQ(ierr);
456d1d35e2fSjeremylt     ierr = MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL,
457ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
458d1d35e2fSjeremylt     ierr = MatShellSetVecType(jacob_mat[level], vectype); CHKERRQ(ierr);
459ccaff030SJeremy L Thompson   }
460e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
461ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
462d1d35e2fSjeremylt   ierr = PetscMalloc1(1, &form_jacob_ctx); CHKERRQ(ierr);
463d1d35e2fSjeremylt   form_jacob_ctx->jacob_ctx = jacob_ctx;
464d1d35e2fSjeremylt   form_jacob_ctx->num_levels = num_levels;
465d1d35e2fSjeremylt   form_jacob_ctx->jacob_mat = jacob_mat;
466ccaff030SJeremy L Thompson 
467ccaff030SJeremy L Thompson   // -- Residual evaluation function
468d1d35e2fSjeremylt   ierr = PetscCalloc1(1, &res_ctx); CHKERRQ(ierr);
469d1d35e2fSjeremylt   ierr = PetscMemcpy(res_ctx, jacob_ctx[fine_level],
470d1d35e2fSjeremylt                      sizeof(*jacob_ctx[fine_level])); CHKERRQ(ierr);
4715754ecacSJeremy L Thompson   res_ctx->op = ceed_data[fine_level]->op_residual;
4725754ecacSJeremy L Thompson   res_ctx->qf = ceed_data[fine_level]->qf_residual;
473d1d35e2fSjeremylt   if (app_ctx->bc_traction_count > 0)
474d1d35e2fSjeremylt     res_ctx->neumann_bcs = neumann_bcs;
475fe394131Sjeremylt   else
476d1d35e2fSjeremylt     res_ctx->neumann_bcs = NULL;
477d1d35e2fSjeremylt   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson 
479ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
480d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_ctx); CHKERRQ(ierr);
481d1d35e2fSjeremylt   ierr = PetscMalloc1(num_levels, &prolong_restr_mat); CHKERRQ(ierr);
482d1d35e2fSjeremylt   for (PetscInt level = 1; level < num_levels; level++) {
483ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
484d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &prolong_restr_ctx[level]); CHKERRQ(ierr);
485d1d35e2fSjeremylt     ierr = SetupProlongRestrictCtx(comm, app_ctx, level_dms[level-1],
486d1d35e2fSjeremylt                                    level_dms[level], U_g[level], U_loc[level-1],
487d1d35e2fSjeremylt                                    U_loc[level], ceed_data[level-1],
488d1d35e2fSjeremylt                                    ceed_data[level], ceed,
489d1d35e2fSjeremylt                                    prolong_restr_ctx[level]); CHKERRQ(ierr);
490ccaff030SJeremy L Thompson 
491ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
492d1d35e2fSjeremylt     ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level-1], U_g_size[level],
493d1d35e2fSjeremylt                           U_g_size[level-1], prolong_restr_ctx[level],
494d1d35e2fSjeremylt                           &prolong_restr_mat[level]); CHKERRQ(ierr);
495ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
496d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT,
497ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
498d1d35e2fSjeremylt     ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE,
499ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
500ccaff030SJeremy L Thompson     CHKERRQ(ierr);
501d1d35e2fSjeremylt     ierr = MatShellSetVecType(prolong_restr_mat[level], vectype); CHKERRQ(ierr);
502ccaff030SJeremy L Thompson   }
503ccaff030SJeremy L Thompson 
504ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
50517f0843fSJeremy L Thompson   // Setup for AMG coarse solve
506ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
507d1d35e2fSjeremylt   if (app_ctx->multigrid_choice != MULTIGRID_NONE) {
508e3e3df41Sjeremylt     // -- Jacobian Matrix
509d1d35e2fSjeremylt     ierr = DMCreateMatrix(level_dms[0], &jacob_mat_coarse); CHKERRQ(ierr);
510e3e3df41Sjeremylt 
511d1d35e2fSjeremylt     if (app_ctx->degree > 1) {
51217f0843fSJeremy L Thompson       // -- Assemble sparsity pattern
5133047f789SJeremy L Thompson       PetscCount num_entries;
5143047f789SJeremy L Thompson       CeedInt *rows, *cols;
51517f0843fSJeremy L Thompson       CeedVector coo_values;
51617f0843fSJeremy L Thompson       CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries,
51717f0843fSJeremy L Thompson                                          &rows, &cols);
51817f0843fSJeremy L Thompson       ISLocalToGlobalMapping ltog_row, ltog_col;
51917f0843fSJeremy L Thompson       ierr = MatGetLocalToGlobalMapping(jacob_mat_coarse, &ltog_row, &ltog_col);
520ccaff030SJeremy L Thompson       CHKERRQ(ierr);
52117f0843fSJeremy L Thompson       ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows);
52217f0843fSJeremy L Thompson       CHKERRQ(ierr);
52317f0843fSJeremy L Thompson       ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols);
52417f0843fSJeremy L Thompson       CHKERRQ(ierr);
52517f0843fSJeremy L Thompson       ierr = MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows, cols);
52617f0843fSJeremy L Thompson       CHKERRQ(ierr);
52717f0843fSJeremy L Thompson       free(rows);
52817f0843fSJeremy L Thompson       free(cols);
52917f0843fSJeremy L Thompson       CeedVectorCreate(ceed, num_entries, &coo_values);
530ccaff030SJeremy L Thompson 
531d1d35e2fSjeremylt       // -- Update form_jacob_ctx
53217f0843fSJeremy L Thompson       form_jacob_ctx->coo_values = coo_values;
53317f0843fSJeremy L Thompson       form_jacob_ctx->op_coarse = ceed_data[0]->op_jacobian;
534d1d35e2fSjeremylt       form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse;
535e3e3df41Sjeremylt     }
536e3e3df41Sjeremylt   }
537e3e3df41Sjeremylt 
538e3e3df41Sjeremylt   // Set Jacobian function
539d1d35e2fSjeremylt   if (app_ctx->degree > 1) {
540d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level],
541d1d35e2fSjeremylt                            FormJacobian, form_jacob_ctx); CHKERRQ(ierr);
542e3e3df41Sjeremylt   } else {
543d1d35e2fSjeremylt     ierr = SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse,
544e3e3df41Sjeremylt                            SNESComputeJacobianDefaultColor, NULL);
545e3e3df41Sjeremylt     CHKERRQ(ierr);
546e3e3df41Sjeremylt   }
547ccaff030SJeremy L Thompson 
548ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
549ccaff030SJeremy L Thompson   // Setup KSP
550ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
551ccaff030SJeremy L Thompson   {
552ccaff030SJeremy L Thompson     PC pc;
553ccaff030SJeremy L Thompson     KSP ksp;
554ccaff030SJeremy L Thompson 
555ccaff030SJeremy L Thompson     // -- KSP
556ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
557ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
560ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
562ccaff030SJeremy L Thompson 
563ccaff030SJeremy L Thompson     // -- Preconditioning
564ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
565d1d35e2fSjeremylt     ierr = PCSetDM(pc, level_dms[fine_level]); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
567ccaff030SJeremy L Thompson 
568d1d35e2fSjeremylt     if (app_ctx->multigrid_choice == MULTIGRID_NONE) {
569ccaff030SJeremy L Thompson       // ---- No Multigrid
570ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
571ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
572d1d35e2fSjeremylt     } else if (app_ctx->degree == 1) {
573f892e6d1Sjeremylt       // ---- AMG for degree 1
574f892e6d1Sjeremylt       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
575ccaff030SJeremy L Thompson     } else {
576ccaff030SJeremy L Thompson       // ---- PCMG
577ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
578ccaff030SJeremy L Thompson 
579ccaff030SJeremy L Thompson       // ------ PCMG levels
580d1d35e2fSjeremylt       ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
581d1d35e2fSjeremylt       for (PetscInt level = 0; level < num_levels; level++) {
582ccaff030SJeremy L Thompson         // -------- Smoother
583d1d35e2fSjeremylt         KSP ksp_smoother, ksp_est;
584d1d35e2fSjeremylt         PC pc_smoother;
585ccaff030SJeremy L Thompson 
586ccaff030SJeremy L Thompson         // ---------- Smoother KSP
587d1d35e2fSjeremylt         ierr = PCMGGetSmoother(pc, level, &ksp_smoother); CHKERRQ(ierr);
588d1d35e2fSjeremylt         ierr = KSPSetDM(ksp_smoother, level_dms[level]); CHKERRQ(ierr);
589d1d35e2fSjeremylt         ierr = KSPSetDMActive(ksp_smoother, PETSC_FALSE); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson 
591ccaff030SJeremy L Thompson         // ---------- Chebyshev options
592d1d35e2fSjeremylt         ierr = KSPSetType(ksp_smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
593d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1);
594ccaff030SJeremy L Thompson         CHKERRQ(ierr);
595d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est); CHKERRQ(ierr);
596d1d35e2fSjeremylt         ierr = KSPSetType(ksp_est, KSPCG); CHKERRQ(ierr);
597d1d35e2fSjeremylt         ierr = KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE);
598ccaff030SJeremy L Thompson         CHKERRQ(ierr);
599d1d35e2fSjeremylt         ierr = KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]);
600ccaff030SJeremy L Thompson         CHKERRQ(ierr);
601ccaff030SJeremy L Thompson 
602ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
603d1d35e2fSjeremylt         ierr = KSPGetPC(ksp_smoother, &pc_smoother); CHKERRQ(ierr);
604d1d35e2fSjeremylt         ierr = PCSetType(pc_smoother, PCJACOBI); CHKERRQ(ierr);
605d1d35e2fSjeremylt         ierr = PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
606ccaff030SJeremy L Thompson 
607ccaff030SJeremy L Thompson         // -------- Work vector
608d1d35e2fSjeremylt         if (level != fine_level) {
609d1d35e2fSjeremylt           ierr = PCMGSetX(pc, level, U_g[level]); CHKERRQ(ierr);
610ccaff030SJeremy L Thompson         }
611ccaff030SJeremy L Thompson 
612ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
613ccaff030SJeremy L Thompson         if (level > 0) {
614d1d35e2fSjeremylt           ierr = PCMGSetInterpolation(pc, level, prolong_restr_mat[level]);
615ccaff030SJeremy L Thompson           CHKERRQ(ierr);
616d1d35e2fSjeremylt           ierr = PCMGSetRestriction(pc, level, prolong_restr_mat[level]);
617ccaff030SJeremy L Thompson           CHKERRQ(ierr);
618ccaff030SJeremy L Thompson         }
619ccaff030SJeremy L Thompson       }
620ccaff030SJeremy L Thompson 
621ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
622d1d35e2fSjeremylt       KSP ksp_coarse;
623d1d35e2fSjeremylt       PC pc_coarse;
624ccaff030SJeremy L Thompson 
625ccaff030SJeremy L Thompson       // -------- Coarse KSP
626d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
627d1d35e2fSjeremylt       ierr = KSPSetType(ksp_coarse, KSPPREONLY); CHKERRQ(ierr);
628d1d35e2fSjeremylt       ierr = KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse);
629ccaff030SJeremy L Thompson       CHKERRQ(ierr);
630d1d35e2fSjeremylt       ierr = KSPSetOptionsPrefix(ksp_coarse, "coarse_"); CHKERRQ(ierr);
631ccaff030SJeremy L Thompson 
632ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
633d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
634d1d35e2fSjeremylt       ierr = PCSetType(pc_coarse, PCGAMG); CHKERRQ(ierr);
635d1d35e2fSjeremylt       ierr = PCSetOptionsPrefix(pc_coarse, "coarse_"); CHKERRQ(ierr);
636ccaff030SJeremy L Thompson 
637d1d35e2fSjeremylt       ierr = KSPSetFromOptions(ksp_coarse); CHKERRQ(ierr);
638d1d35e2fSjeremylt       ierr = PCSetFromOptions(pc_coarse); CHKERRQ(ierr);
639ccaff030SJeremy L Thompson 
640ccaff030SJeremy L Thompson       // ------ PCMG options
641ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
643d1d35e2fSjeremylt       ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
644ccaff030SJeremy L Thompson     }
645ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
646ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
647ccaff030SJeremy L Thompson   }
648038d0cd7Sjeremylt   {
649038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
650d1d35e2fSjeremylt     SNESLineSearch line_search;
651481a4840SJed Brown 
652d1d35e2fSjeremylt     ierr = SNESGetLineSearch(snes, &line_search); CHKERRQ(ierr);
653d1d35e2fSjeremylt     ierr = SNESLineSearchSetType(line_search, SNESLINESEARCHCP); CHKERRQ(ierr);
654481a4840SJed Brown   }
655481a4840SJed Brown 
656ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
657ccaff030SJeremy L Thompson 
658ccaff030SJeremy L Thompson   // Performance logging
659ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
660ccaff030SJeremy L Thompson 
661ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
662ccaff030SJeremy L Thompson   // Set initial guess
663ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
664f81c27eaSJed Brown   ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr);
665ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
666ccaff030SJeremy L Thompson 
667ccaff030SJeremy L Thompson   // View solution
668d1d35e2fSjeremylt   if (app_ctx->view_soln) {
669d1d35e2fSjeremylt     ierr = ViewSolution(comm, app_ctx, U, 0, 0.0); CHKERRQ(ierr);
670ccaff030SJeremy L Thompson   }
671ccaff030SJeremy L Thompson 
672ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
673ccaff030SJeremy L Thompson   // Solve SNES
674ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
675d1d35e2fSjeremylt   PetscBool snes_monitor = PETSC_FALSE;
676d1d35e2fSjeremylt   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor);
6775f24f028Sjeremylt   CHKERRQ(ierr);
6785f24f028Sjeremylt 
679ccaff030SJeremy L Thompson   // Performance logging
680d1d35e2fSjeremylt   ierr = PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve);
681ccaff030SJeremy L Thompson   CHKERRQ(ierr);
682d1d35e2fSjeremylt   ierr = PetscLogStagePush(stage_snes_solve); CHKERRQ(ierr);
683ccaff030SJeremy L Thompson 
684ccaff030SJeremy L Thompson   // Timing
685ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
686d1d35e2fSjeremylt   start_time = MPI_Wtime();
687ccaff030SJeremy L Thompson 
688ccaff030SJeremy L Thompson   // Solve for each load increment
6895f24f028Sjeremylt   PetscInt increment;
690d1d35e2fSjeremylt   for (increment = 1; increment <= app_ctx->num_increments; increment++) {
6915f24f028Sjeremylt     // -- Log increment count
692d1d35e2fSjeremylt     if (snes_monitor) {
6935f24f028Sjeremylt       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
6945f24f028Sjeremylt       CHKERRQ(ierr);
6955f24f028Sjeremylt     }
6965f24f028Sjeremylt 
697ccaff030SJeremy L Thompson     // -- Scale the problem
698d1d35e2fSjeremylt     PetscScalar load_increment = 1.0*increment / app_ctx->num_increments,
699d1d35e2fSjeremylt                 scalingFactor = load_increment /
700d1d35e2fSjeremylt                                 (increment == 1 ? 1 : res_ctx->load_increment);
701d1d35e2fSjeremylt     res_ctx->load_increment = load_increment;
702d1d35e2fSjeremylt     if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) {
703ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
704ccaff030SJeremy L Thompson     }
705ccaff030SJeremy L Thompson 
706ccaff030SJeremy L Thompson     // -- Solve
707ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
708ccaff030SJeremy L Thompson 
709ccaff030SJeremy L Thompson     // -- View solution
710d1d35e2fSjeremylt     if (app_ctx->view_soln) {
711d1d35e2fSjeremylt       ierr = ViewSolution(comm, app_ctx, U, increment, load_increment); CHKERRQ(ierr);
712ccaff030SJeremy L Thompson     }
713ccaff030SJeremy L Thompson 
714ccaff030SJeremy L Thompson     // -- Update SNES iteration count
715ccaff030SJeremy L Thompson     PetscInt its;
716ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
717bf0f51feSjeremylt     snes_its += its;
7187418ca88SJeremy L Thompson     ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr);
719bf0f51feSjeremylt     ksp_its += its;
7203951731eSjeremylt 
7213951731eSjeremylt     // -- Check for divergence
7223951731eSjeremylt     SNESConvergedReason reason;
7233951731eSjeremylt     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
7243951731eSjeremylt     if (reason < 0)
7253951731eSjeremylt       break;
7268355605fSKaren (Ren) Stengel     if (app_ctx->energy_viewer) {
7278355605fSKaren (Ren) Stengel       // -- Log strain energy for current load increment
7288355605fSKaren (Ren) Stengel       CeedScalar energy;
7298355605fSKaren (Ren) Stengel       ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
7308355605fSKaren (Ren) Stengel                                  U, &energy); CHKERRQ(ierr);
7318355605fSKaren (Ren) Stengel 
7328355605fSKaren (Ren) Stengel       if (!app_ctx->test_mode) {
7338355605fSKaren (Ren) Stengel         // -- Output
7348355605fSKaren (Ren) Stengel         ierr = PetscPrintf(comm,
7358355605fSKaren (Ren) Stengel                            "    Strain Energy                      : %.12e\n",
7368355605fSKaren (Ren) Stengel                            energy); CHKERRQ(ierr);
7378355605fSKaren (Ren) Stengel       }
7388355605fSKaren (Ren) Stengel       ierr = PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment,
7398355605fSKaren (Ren) Stengel                                     energy); CHKERRQ(ierr);
7408355605fSKaren (Ren) Stengel     }
741ccaff030SJeremy L Thompson   }
742ccaff030SJeremy L Thompson 
743ccaff030SJeremy L Thompson   // Timing
744d1d35e2fSjeremylt   elapsed_time = MPI_Wtime() - start_time;
745ccaff030SJeremy L Thompson 
746ccaff030SJeremy L Thompson   // Performance logging
747ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
748ccaff030SJeremy L Thompson 
749ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
750ccaff030SJeremy L Thompson   // Output summary
751ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
752d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
753ccaff030SJeremy L Thompson     // -- SNES
754d1d35e2fSjeremylt     SNESType snes_type;
755ccaff030SJeremy L Thompson     SNESConvergedReason reason;
756ccaff030SJeremy L Thompson     PetscReal rnorm;
757d1d35e2fSjeremylt     ierr = SNESGetType(snes, &snes_type); CHKERRQ(ierr);
758ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
759ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
760ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
761ccaff030SJeremy L Thompson                        "  SNES:\n"
762ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
763ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
764ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
7655f24f028Sjeremylt                        "    Completed Load Increments          : %d\n"
766ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
767ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
768d1d35e2fSjeremylt                        snes_type, SNESConvergedReasons[reason],
769d1d35e2fSjeremylt                        app_ctx->num_increments, increment - 1,
770bf0f51feSjeremylt                        snes_its, (double)rnorm); CHKERRQ(ierr);
771ccaff030SJeremy L Thompson 
772ccaff030SJeremy L Thompson     // -- KSP
773ccaff030SJeremy L Thompson     KSP ksp;
774d1d35e2fSjeremylt     KSPType ksp_type;
775ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
776d1d35e2fSjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
777ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
778ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
7797418ca88SJeremy L Thompson                        "    KSP Type                           : %s\n"
7807418ca88SJeremy L Thompson                        "    Total KSP Iterations               : %D\n",
781bf0f51feSjeremylt                        ksp_type, ksp_its); CHKERRQ(ierr);
782ccaff030SJeremy L Thompson 
783ccaff030SJeremy L Thompson     // -- PC
784ccaff030SJeremy L Thompson     PC pc;
785d1d35e2fSjeremylt     PCType pc_type;
786ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
787d1d35e2fSjeremylt     ierr = PCGetType(pc, &pc_type); CHKERRQ(ierr);
788e3e3df41Sjeremylt     ierr = PetscPrintf(comm,
789e3e3df41Sjeremylt                        "    PC Type                            : %s\n",
790d1d35e2fSjeremylt                        pc_type); CHKERRQ(ierr);
791e3e3df41Sjeremylt 
792d1d35e2fSjeremylt     if (!strcmp(pc_type, PCMG)) {
793d1d35e2fSjeremylt       PCMGType pcmg_type;
794d1d35e2fSjeremylt       ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
795ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
796ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
797ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
798ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
799d1d35e2fSjeremylt                          PCMGTypes[pcmg_type],
800d1d35e2fSjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
801ccaff030SJeremy L Thompson 
802ccaff030SJeremy L Thompson       // -- Coarse Solve
803d1d35e2fSjeremylt       KSP ksp_coarse;
804d1d35e2fSjeremylt       PC pc_coarse;
805d1d35e2fSjeremylt       PCType pc_type;
806ccaff030SJeremy L Thompson 
807d1d35e2fSjeremylt       ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr);
808d1d35e2fSjeremylt       ierr = KSPGetType(ksp_coarse, &ksp_type); CHKERRQ(ierr);
809d1d35e2fSjeremylt       ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr);
810d1d35e2fSjeremylt       ierr = PCGetType(pc_coarse, &pc_type); CHKERRQ(ierr);
811ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
812ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
813ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
814ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
815d1d35e2fSjeremylt                          ksp_type, pc_type); CHKERRQ(ierr);
816ccaff030SJeremy L Thompson     }
817ccaff030SJeremy L Thompson   }
818ccaff030SJeremy L Thompson 
819ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
820ccaff030SJeremy L Thompson   // Compute solve time
821ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
822d1d35e2fSjeremylt   if (!app_ctx->test_mode) {
823d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm);
824ccaff030SJeremy L Thompson     CHKERRQ(ierr);
825d1d35e2fSjeremylt     ierr = MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm);
826ccaff030SJeremy L Thompson     CHKERRQ(ierr);
827ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
828ccaff030SJeremy L Thompson                        "  Performance:\n"
8297418ca88SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n"
8307418ca88SJeremy L Thompson                        "    DoFs/Sec in SNES                   : %g (%g) million\n",
831bf0f51feSjeremylt                        max_time, min_time, 1e-6*U_g_size[fine_level]*ksp_its/max_time,
832bf0f51feSjeremylt                        1e-6*U_g_size[fine_level]*ksp_its/min_time); CHKERRQ(ierr);
833ccaff030SJeremy L Thompson   }
834ccaff030SJeremy L Thompson 
835ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
836ccaff030SJeremy L Thompson   // Compute error
837ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
838d1d35e2fSjeremylt   if (app_ctx->forcing_choice == FORCE_MMS) {
839d1d35e2fSjeremylt     CeedScalar l2_error = 1., l2_U_norm = 1.;
840d1d35e2fSjeremylt     const CeedScalar *true_array;
841d1d35e2fSjeremylt     Vec error_vec, true_vec;
842ccaff030SJeremy L Thompson 
843ccaff030SJeremy L Thompson     // -- Work vectors
844d1d35e2fSjeremylt     ierr = VecDuplicate(U, &error_vec); CHKERRQ(ierr);
845d1d35e2fSjeremylt     ierr = VecSet(error_vec, 0.0); CHKERRQ(ierr);
846d1d35e2fSjeremylt     ierr = VecDuplicate(U, &true_vec); CHKERRQ(ierr);
847d1d35e2fSjeremylt     ierr = VecSet(true_vec, 0.0); CHKERRQ(ierr);
848ccaff030SJeremy L Thompson 
849ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
850d1d35e2fSjeremylt     CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln,
851d1d35e2fSjeremylt                            CEED_MEM_HOST, &true_array);
852d1d35e2fSjeremylt     ierr = VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array);
85362e9c006SJeremy L Thompson     CHKERRQ(ierr);
854d1d35e2fSjeremylt     ierr = DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec);
855ccaff030SJeremy L Thompson     CHKERRQ(ierr);
856d1d35e2fSjeremylt     ierr = VecResetArray(res_ctx->Y_loc); CHKERRQ(ierr);
857d1d35e2fSjeremylt     CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array);
858ccaff030SJeremy L Thompson 
859ccaff030SJeremy L Thompson     // -- Compute L2 error
860d1d35e2fSjeremylt     ierr = VecWAXPY(error_vec, -1.0, U, true_vec); CHKERRQ(ierr);
861d1d35e2fSjeremylt     ierr = VecNorm(error_vec, NORM_2, &l2_error); CHKERRQ(ierr);
862d1d35e2fSjeremylt     ierr = VecNorm(U, NORM_2, &l2_U_norm); CHKERRQ(ierr);
863d1d35e2fSjeremylt     l2_error /= l2_U_norm;
864ccaff030SJeremy L Thompson 
865ccaff030SJeremy L Thompson     // -- Output
866d1d35e2fSjeremylt     if (!app_ctx->test_mode || l2_error > 0.05) {
867ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
868ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
869d1d35e2fSjeremylt                          l2_error); CHKERRQ(ierr);
870ccaff030SJeremy L Thompson     }
871ccaff030SJeremy L Thompson 
872ccaff030SJeremy L Thompson     // -- Cleanup
873d1d35e2fSjeremylt     ierr = VecDestroy(&error_vec); CHKERRQ(ierr);
874d1d35e2fSjeremylt     ierr = VecDestroy(&true_vec); CHKERRQ(ierr);
8752d93065eSjeremylt   }
8762d93065eSjeremylt 
8772d93065eSjeremylt   // ---------------------------------------------------------------------------
8782d93065eSjeremylt   // Compute energy
8792d93065eSjeremylt   // ---------------------------------------------------------------------------
8808355605fSKaren (Ren) Stengel   PetscReal energy;
881d1d35e2fSjeremylt   ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy,
882a3c02c40SJeremy L Thompson                              U, &energy); CHKERRQ(ierr);
8838355605fSKaren (Ren) Stengel   if (!app_ctx->test_mode) {
8842d93065eSjeremylt     // -- Output
8852d93065eSjeremylt     ierr = PetscPrintf(comm,
88672d03b64SArash Mehraban                        "    Strain Energy                      : %.12e\n",
8872d93065eSjeremylt                        energy); CHKERRQ(ierr);
888ccaff030SJeremy L Thompson   }
8898355605fSKaren (Ren) Stengel   ierr = RegressionTests_solids(app_ctx, energy); CHKERRQ(ierr);
890ccaff030SJeremy L Thompson 
891ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
8925c25879aSJeremy L Thompson   // Output diagnostic quantities
8935c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
894d1d35e2fSjeremylt   if (app_ctx->view_soln || app_ctx->view_final_soln) {
8955c25879aSJeremy L Thompson     // -- Setup context
896d1d35e2fSjeremylt     UserMult diagnostic_ctx;
897d1d35e2fSjeremylt     ierr = PetscMalloc1(1, &diagnostic_ctx); CHKERRQ(ierr);
898d1d35e2fSjeremylt     ierr = PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)); CHKERRQ(ierr);
899d1d35e2fSjeremylt     diagnostic_ctx->dm = dm_diagnostic;
900d1d35e2fSjeremylt     diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic;
9015c25879aSJeremy L Thompson 
9025c25879aSJeremy L Thompson     // -- Compute and output
903d1d35e2fSjeremylt     ierr = ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx,
904d1d35e2fSjeremylt                                     app_ctx, U,
905d1d35e2fSjeremylt                                     ceed_data[fine_level]->elem_restr_diagnostic);
9065c25879aSJeremy L Thompson     CHKERRQ(ierr);
9075c25879aSJeremy L Thompson 
9085c25879aSJeremy L Thompson     // -- Cleanup
909d1d35e2fSjeremylt     ierr = PetscFree(diagnostic_ctx); CHKERRQ(ierr);
9105c25879aSJeremy L Thompson   }
9115c25879aSJeremy L Thompson 
9125c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
913ccaff030SJeremy L Thompson   // Free objects
914ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
915ccaff030SJeremy L Thompson   // Data in arrays per level
916d1d35e2fSjeremylt   for (PetscInt level = 0; level < num_levels; level++) {
917ccaff030SJeremy L Thompson     // Vectors
918d1d35e2fSjeremylt     ierr = VecDestroy(&U_g[level]); CHKERRQ(ierr);
919d1d35e2fSjeremylt     ierr = VecDestroy(&U_loc[level]); CHKERRQ(ierr);
920ccaff030SJeremy L Thompson 
921ccaff030SJeremy L Thompson     // Jacobian matrix and data
922d1d35e2fSjeremylt     ierr = VecDestroy(&jacob_ctx[level]->Y_loc); CHKERRQ(ierr);
923d1d35e2fSjeremylt     ierr = MatDestroy(&jacob_mat[level]); CHKERRQ(ierr);
924d1d35e2fSjeremylt     ierr = PetscFree(jacob_ctx[level]); CHKERRQ(ierr);
925ccaff030SJeremy L Thompson 
926ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
927ccaff030SJeremy L Thompson     if (level > 0) {
928d1d35e2fSjeremylt       ierr = PetscFree(prolong_restr_ctx[level]); CHKERRQ(ierr);
929d1d35e2fSjeremylt       ierr = MatDestroy(&prolong_restr_mat[level]); CHKERRQ(ierr);
930ccaff030SJeremy L Thompson     }
931ccaff030SJeremy L Thompson 
932ccaff030SJeremy L Thompson     // DM
933d1d35e2fSjeremylt     ierr = DMDestroy(&level_dms[level]); CHKERRQ(ierr);
934ccaff030SJeremy L Thompson 
935ccaff030SJeremy L Thompson     // libCEED objects
936d1d35e2fSjeremylt     ierr = CeedDataDestroy(level, ceed_data[level]); CHKERRQ(ierr);
937ccaff030SJeremy L Thompson   }
938ccaff030SJeremy L Thompson 
9398355605fSKaren (Ren) Stengel   ierr = PetscViewerDestroy(&app_ctx->energy_viewer); CHKERRQ(ierr);
9408355605fSKaren (Ren) Stengel 
941ccaff030SJeremy L Thompson   // Arrays
942d1d35e2fSjeremylt   ierr = PetscFree(U_g); CHKERRQ(ierr);
943d1d35e2fSjeremylt   ierr = PetscFree(U_loc); CHKERRQ(ierr);
944d1d35e2fSjeremylt   ierr = PetscFree(U_g_size); CHKERRQ(ierr);
945d1d35e2fSjeremylt   ierr = PetscFree(U_l_size); CHKERRQ(ierr);
946d1d35e2fSjeremylt   ierr = PetscFree(U_loc_size); CHKERRQ(ierr);
947d1d35e2fSjeremylt   ierr = PetscFree(jacob_ctx); CHKERRQ(ierr);
948d1d35e2fSjeremylt   ierr = PetscFree(jacob_mat); CHKERRQ(ierr);
949d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_ctx); CHKERRQ(ierr);
950d1d35e2fSjeremylt   ierr = PetscFree(prolong_restr_mat); CHKERRQ(ierr);
951d1d35e2fSjeremylt   ierr = PetscFree(app_ctx->level_degrees); CHKERRQ(ierr);
952d1d35e2fSjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
953ccaff030SJeremy L Thompson 
954ccaff030SJeremy L Thompson   // libCEED objects
95517f0843fSJeremy L Thompson   CeedVectorDestroy(&form_jacob_ctx->coo_values);
956d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys);
957d1d35e2fSjeremylt   CeedQFunctionContextDestroy(&ctx_phys_smoother);
958ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
959ccaff030SJeremy L Thompson 
960ccaff030SJeremy L Thompson   // PETSc objects
961ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
962ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
963d1d35e2fSjeremylt   ierr = VecDestroy(&R_loc); CHKERRQ(ierr);
964ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
965d1d35e2fSjeremylt   ierr = VecDestroy(&F_loc); CHKERRQ(ierr);
966d1d35e2fSjeremylt   ierr = VecDestroy(&neumann_bcs); CHKERRQ(ierr);
967d1d35e2fSjeremylt   ierr = VecDestroy(&bcs_loc); CHKERRQ(ierr);
968d1d35e2fSjeremylt   ierr = MatDestroy(&jacob_mat_coarse); CHKERRQ(ierr);
969ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
970d1d35e2fSjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
971d1d35e2fSjeremylt   ierr = DMDestroy(&dm_energy); CHKERRQ(ierr);
972d1d35e2fSjeremylt   ierr = DMDestroy(&dm_diagnostic); CHKERRQ(ierr);
973d1d35e2fSjeremylt   ierr = PetscFree(level_dms); CHKERRQ(ierr);
974ccaff030SJeremy L Thompson 
9755754ecacSJeremy L Thompson   // -- Function list
9765754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupPhysics);
9775754ecacSJeremy L Thompson   CHKERRQ(ierr);
9785754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel);
9795754ecacSJeremy L Thompson   CHKERRQ(ierr);
9805754ecacSJeremy L Thompson   ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedLevel);
9815754ecacSJeremy L Thompson   CHKERRQ(ierr);
9825754ecacSJeremy L Thompson 
983ccaff030SJeremy L Thompson   // Structs
984d1d35e2fSjeremylt   ierr = PetscFree(res_ctx); CHKERRQ(ierr);
985d1d35e2fSjeremylt   ierr = PetscFree(form_jacob_ctx); CHKERRQ(ierr);
986d1d35e2fSjeremylt   ierr = PetscFree(jacob_coarse_ctx); CHKERRQ(ierr);
987d1d35e2fSjeremylt   ierr = PetscFree(app_ctx); CHKERRQ(ierr);
9885754ecacSJeremy L Thompson   ierr = PetscFree(problem_functions); CHKERRQ(ierr);
989ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
990ccaff030SJeremy L Thompson 
991ccaff030SJeremy L Thompson   return PetscFinalize();
992ccaff030SJeremy L Thompson }
993