xref: /libCEED/examples/solids/elasticity.c (revision d99fa3c5cd91a1690aedf0679cbf290d44fec74c)
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 //
30ccaff030SJeremy L Thompson //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem linElas -forcing mms
31aee2786aSjeremylt //     ./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 hyperSS -forcing none -ceed /cpu/self
32b04a8a52SLeila Ghaffari //     ./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 hyperFS -forcing none -ceed /gpu/occa
33ccaff030SJeremy L Thompson //
34ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
35ccaff030SJeremy L Thompson //
36da5d3694SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -degree 3 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson /// @file
39ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson #include "elasticity.h"
44ccaff030SJeremy L Thompson 
45ccaff030SJeremy L Thompson int main(int argc, char **argv) {
46ccaff030SJeremy L Thompson   PetscInt       ierr;
47ccaff030SJeremy L Thompson   MPI_Comm       comm;
48ccaff030SJeremy L Thompson   // Context structs
49ccaff030SJeremy L Thompson   AppCtx         appCtx;                 // Contains problem options
50ccaff030SJeremy L Thompson   Physics        phys;                   // Contains physical constants
51f7b4142eSJeremy L Thompson   Physics        physSmoother = NULL;    // Separate context if nuSmoother set
52ccaff030SJeremy L Thompson   Units          units;                  // Contains units scaling
53ccaff030SJeremy L Thompson   // PETSc objects
54ccaff030SJeremy L Thompson   PetscLogStage  stageDMSetup, stageLibceedSetup,
55ccaff030SJeremy L Thompson                  stageSnesSetup, stageSnesSolve;
56ccaff030SJeremy L Thompson   DM             dmOrig;                 // Distributed DM to clone
57a3c02c40SJeremy L Thompson   DM             dmEnergy, dmDiagnostic; // DMs for postprocessing
58ccaff030SJeremy L Thompson   DM             *levelDMs;
59ccaff030SJeremy L Thompson   Vec            U, *Ug, *Uloc;          // U: solution, R: residual, F: forcing
60ccaff030SJeremy L Thompson   Vec            R, Rloc, F, Floc;       // g: global, loc: local
61e3e3df41Sjeremylt   SNES           snes, snesCoarse = NULL;
62ccaff030SJeremy L Thompson   Mat            *jacobMat, jacobMatCoarse, *prolongRestrMat;
63ccaff030SJeremy L Thompson   // PETSc data
64e3e3df41Sjeremylt   UserMult       resCtx, jacobCoarseCtx = NULL, *jacobCtx;
65ccaff030SJeremy L Thompson   FormJacobCtx   formJacobCtx;
66ccaff030SJeremy L Thompson   UserMultProlongRestr *prolongRestrCtx;
67ccaff030SJeremy L Thompson   PCMGCycleType  pcmgCycleType = PC_MG_CYCLE_V;
68ccaff030SJeremy L Thompson   // libCEED objects
69*d99fa3c5SJeremy L Thompson   Ceed           ceed;
70ccaff030SJeremy L Thompson   CeedData       *ceedData;
71f892e6d1Sjeremylt   CeedQFunction  qfRestrict = NULL, qfProlong = NULL;
72ccaff030SJeremy L Thompson   // Parameters
73ccaff030SJeremy L Thompson   PetscInt       ncompu = 3;             // 3 DoFs in 3D
7413fdad4bSJeremy L Thompson   PetscInt       ncompe = 1, ncompd = 5; // 1 energy output, 5 diagnostic
75ccaff030SJeremy L Thompson   PetscInt       numLevels = 1, fineLevel = 0;
76ccaff030SJeremy L Thompson   PetscInt       *Ugsz, *Ulsz, *Ulocsz;  // sz: size
777418ca88SJeremy L Thompson   PetscInt       snesIts = 0, kspIts = 0;
78ccaff030SJeremy L Thompson   // Timing
79ccaff030SJeremy L Thompson   double         startTime, elapsedTime, minTime, maxTime;
80ccaff030SJeremy L Thompson 
81ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
82ccaff030SJeremy L Thompson   if (ierr)
83ccaff030SJeremy L Thompson     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
91ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr);
92ccaff030SJeremy L Thompson   ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr);
93ccaff030SJeremy L Thompson   numLevels = appCtx->numLevels;
94ccaff030SJeremy L Thompson   fineLevel = numLevels - 1;
95ccaff030SJeremy L Thompson 
96ccaff030SJeremy L Thompson   // -- Set Poison's ratio, Young's Modulus
97ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
98ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
99ccaff030SJeremy L Thompson   ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr);
100f7b4142eSJeremy L Thompson   if (fabs(appCtx->nuSmoother) > 1E-14) {
101f7b4142eSJeremy L Thompson     ierr = PetscMalloc1(1, &physSmoother); CHKERRQ(ierr);
102f7b4142eSJeremy L Thompson     ierr = PetscMemcpy(physSmoother, phys, sizeof(*phys)); CHKERRQ(ierr);
103f7b4142eSJeremy L Thompson     physSmoother->nu = appCtx->nuSmoother;
104f7b4142eSJeremy L Thompson   }
105ccaff030SJeremy L Thompson 
106ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
10762e9c006SJeremy L Thompson   // Initalize libCEED
10862e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
10962e9c006SJeremy L Thompson   // Initalize backend
11062e9c006SJeremy L Thompson   CeedInit(appCtx->ceedResource, &ceed);
11162e9c006SJeremy L Thompson 
11262e9c006SJeremy L Thompson   // Check preferred MemType
11362e9c006SJeremy L Thompson   CeedMemType memTypeBackend;
11462e9c006SJeremy L Thompson   CeedGetPreferredMemType(ceed, &memTypeBackend);
11562e9c006SJeremy L Thompson   if (!appCtx->setMemTypeRequest)
11662e9c006SJeremy L Thompson     appCtx->memTypeRequested = memTypeBackend;
11762e9c006SJeremy L Thompson   else if (!appCtx->petscHaveCuda && appCtx->memTypeRequested == CEED_MEM_DEVICE)
11862e9c006SJeremy L Thompson     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
11962e9c006SJeremy L Thompson              "PETSc was not built with CUDA. "
12062e9c006SJeremy L Thompson              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
12162e9c006SJeremy L Thompson 
12262e9c006SJeremy L Thompson   // ---------------------------------------------------------------------------
123ccaff030SJeremy L Thompson   // Setup DM
124ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
125ccaff030SJeremy L Thompson   // Performance logging
126ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup);
127ccaff030SJeremy L Thompson   CHKERRQ(ierr);
128ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr);
129ccaff030SJeremy L Thompson 
130ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
131ccaff030SJeremy L Thompson   ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr);
132ccaff030SJeremy L Thompson 
133ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
134ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr);
135e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
136ccaff030SJeremy L Thompson     ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr);
13762e9c006SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_DEVICE) {
13862e9c006SJeremy L Thompson       ierr = DMSetVecType(levelDMs[level], VECCUDA); CHKERRQ(ierr);
13962e9c006SJeremy L Thompson     }
140ccaff030SJeremy L Thompson     ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level],
141a3c02c40SJeremy L Thompson                            PETSC_TRUE, ncompu); CHKERRQ(ierr);
142a3c02c40SJeremy L Thompson     // -- Label field components for viewing
143a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
144a3c02c40SJeremy L Thompson     PetscSection section;
145a3c02c40SJeremy L Thompson     ierr = DMGetLocalSection(levelDMs[level], &section); CHKERRQ(ierr);
146a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr);
147a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
148a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
149a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
150a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
151a3c02c40SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
152a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
153a3c02c40SJeremy L Thompson   }
154a3c02c40SJeremy L Thompson 
155a3c02c40SJeremy L Thompson   // -- Setup postprocessing DMs
156a3c02c40SJeremy L Thompson   ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr);
157a3c02c40SJeremy L Thompson   ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel],
1585c25879aSJeremy L Thompson                          PETSC_FALSE, ncompe); CHKERRQ(ierr);
159a3c02c40SJeremy L Thompson   ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr);
160a3c02c40SJeremy L Thompson   ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel],
1618ca58ff3SJeremy L Thompson                          PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr);
16262e9c006SJeremy L Thompson   if (appCtx->memTypeRequested == CEED_MEM_DEVICE) {
16362e9c006SJeremy L Thompson     ierr = DMSetVecType(dmEnergy, VECCUDA); CHKERRQ(ierr);
16462e9c006SJeremy L Thompson     ierr = DMSetVecType(dmDiagnostic, VECCUDA); CHKERRQ(ierr);
16562e9c006SJeremy L Thompson   }
166a3c02c40SJeremy L Thompson   {
167a3c02c40SJeremy L Thompson     // -- Label field components for viewing
168a3c02c40SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
169a3c02c40SJeremy L Thompson     PetscSection section;
170a3c02c40SJeremy L Thompson     ierr = DMGetLocalSection(dmDiagnostic, &section); CHKERRQ(ierr);
171a3c02c40SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr);
1728ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX");
1738ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1748ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY");
1758ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1768ca58ff3SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ");
1778ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
178379387d4SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure");
1798ca58ff3SJeremy L Thompson     CHKERRQ(ierr);
1807ab18febSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain");
18113fdad4bSJeremy L Thompson     CHKERRQ(ierr);
18213fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2");
18313fdad4bSJeremy L Thompson     CHKERRQ(ierr);
18413fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 6, "detJ");
18513fdad4bSJeremy L Thompson     CHKERRQ(ierr);
18613fdad4bSJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity");
187a3c02c40SJeremy L Thompson     CHKERRQ(ierr);
188ccaff030SJeremy L Thompson   }
189ccaff030SJeremy L Thompson 
190ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
191ccaff030SJeremy L Thompson   // Setup solution and work vectors
192ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
193ccaff030SJeremy L Thompson   // Allocate arrays
194ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr);
195ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr);
196ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr);
197ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr);
198ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr);
199ccaff030SJeremy L Thompson 
200ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
201e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
202ccaff030SJeremy L Thompson     // -- Create global unknown vector U
203ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr);
204ccaff030SJeremy L Thompson     ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr);
205ccaff030SJeremy L Thompson     // Note: Local size for matShell
206ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr);
207ccaff030SJeremy L Thompson 
208ccaff030SJeremy L Thompson     // -- Create local unknown vector Uloc
209ccaff030SJeremy L Thompson     ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr);
210ccaff030SJeremy L Thompson     // Note: local size for libCEED
211ccaff030SJeremy L Thompson     ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr);
212ccaff030SJeremy L Thompson   }
213ccaff030SJeremy L Thompson 
214ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
215ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr);
216ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr);
217ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr);
218ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr);
219ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr);
220ccaff030SJeremy L Thompson 
221ccaff030SJeremy L Thompson   // Performance logging
222ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
223ccaff030SJeremy L Thompson 
224ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
225ccaff030SJeremy L Thompson   // Set up libCEED
226ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
227ccaff030SJeremy L Thompson   // Performance logging
228ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup);
229ccaff030SJeremy L Thompson   CHKERRQ(ierr);
230ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr);
231ccaff030SJeremy L Thompson 
232ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
233ccaff030SJeremy L Thompson   CeedVector forceCeed;
234ccaff030SJeremy L Thompson   CeedScalar *f;
235ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
23662e9c006SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
237ccaff030SJeremy L Thompson       ierr = VecGetArray(Floc, &f); CHKERRQ(ierr);
23862e9c006SJeremy L Thompson     } else {
23962e9c006SJeremy L Thompson       ierr = VecCUDAGetArray(Floc, &f); CHKERRQ(ierr);
24062e9c006SJeremy L Thompson     }
241ccaff030SJeremy L Thompson     CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed);
24262e9c006SJeremy L Thompson     CeedVectorSetArray(forceCeed, appCtx->memTypeRequested, CEED_USE_POINTER, f);
243ccaff030SJeremy L Thompson   }
244ccaff030SJeremy L Thompson 
245ccaff030SJeremy L Thompson   // -- Restriction and prolongation QFunction
246ccaff030SJeremy L Thompson   if (appCtx->multigridChoice != MULTIGRID_NONE) {
247ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
248ccaff030SJeremy L Thompson                                 &qfRestrict);
249ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
250ccaff030SJeremy L Thompson                                 &qfProlong);
251ccaff030SJeremy L Thompson   }
252ccaff030SJeremy L Thompson 
253ccaff030SJeremy L Thompson   // -- Setup libCEED objects
254ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr);
255*d99fa3c5SJeremy L Thompson   // ---- Setup residual, Jacobian evaluator and geometric information
256ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr);
257ccaff030SJeremy L Thompson   {
258a3c02c40SJeremy L Thompson     ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic,
259*d99fa3c5SJeremy L Thompson                                  ceed, appCtx, phys, ceedData, fineLevel,
260a3c02c40SJeremy L Thompson                                  ncompu, Ugsz[fineLevel], Ulocsz[fineLevel],
261a3c02c40SJeremy L Thompson                                  forceCeed, qfRestrict, qfProlong);
262ccaff030SJeremy L Thompson     CHKERRQ(ierr);
263ccaff030SJeremy L Thompson   }
264*d99fa3c5SJeremy L Thompson   // ---- Setup coarse Jacobian evaluator and prolongation/restriction
265*d99fa3c5SJeremy L Thompson   for (PetscInt level = numLevels - 2; level >= 0; level--) {
266ccaff030SJeremy L Thompson     ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr);
267*d99fa3c5SJeremy L Thompson 
268*d99fa3c5SJeremy L Thompson     // Get global communication restriction
269*d99fa3c5SJeremy L Thompson     ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr);
270*d99fa3c5SJeremy L Thompson     ierr = VecSet(Uloc[level+1], 1.0); CHKERRQ(ierr);
271*d99fa3c5SJeremy L Thompson     ierr = DMLocalToGlobal(levelDMs[level+1], Uloc[level+1], ADD_VALUES,
272*d99fa3c5SJeremy L Thompson                            Ug[level+1]); CHKERRQ(ierr);
273*d99fa3c5SJeremy L Thompson     ierr = DMGlobalToLocal(levelDMs[level+1], Ug[level+1], INSERT_VALUES,
274*d99fa3c5SJeremy L Thompson                            Uloc[level+1]); CHKERRQ(ierr);
275*d99fa3c5SJeremy L Thompson 
276*d99fa3c5SJeremy L Thompson     // Place in libCEED array
277*d99fa3c5SJeremy L Thompson     const PetscScalar *m;
278*d99fa3c5SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
279*d99fa3c5SJeremy L Thompson       ierr = VecGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr);
280*d99fa3c5SJeremy L Thompson     } else {
281*d99fa3c5SJeremy L Thompson       ierr = VecCUDAGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr);
282ccaff030SJeremy L Thompson     }
283*d99fa3c5SJeremy L Thompson     CeedVectorSetArray(ceedData[level+1]->xceed, appCtx->memTypeRequested,
284*d99fa3c5SJeremy L Thompson                        CEED_USE_POINTER, (CeedScalar *)m);
285ccaff030SJeremy L Thompson 
2861c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
287*d99fa3c5SJeremy L Thompson     ierr = SetupLibceedLevel(levelDMs[level], ceed, appCtx, phys,
288ccaff030SJeremy L Thompson                              ceedData, level, ncompu, Ugsz[level],
289*d99fa3c5SJeremy L Thompson                              Ulocsz[level], ceedData[level+1]->xceed, qfRestrict,
290ccaff030SJeremy L Thompson                              qfProlong); CHKERRQ(ierr);
291*d99fa3c5SJeremy L Thompson 
292*d99fa3c5SJeremy L Thompson     // Restore PETSc vector
293*d99fa3c5SJeremy L Thompson     CeedVectorTakeArray(ceedData[level+1]->xceed, appCtx->memTypeRequested,
294*d99fa3c5SJeremy L Thompson                         (CeedScalar **)&m);
295*d99fa3c5SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
296*d99fa3c5SJeremy L Thompson       ierr = VecRestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr);
297*d99fa3c5SJeremy L Thompson     } else {
298*d99fa3c5SJeremy L Thompson       ierr = VecCUDARestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr);
299*d99fa3c5SJeremy L Thompson     }
300*d99fa3c5SJeremy L Thompson     ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr);
301*d99fa3c5SJeremy L Thompson     ierr = VecZeroEntries(Uloc[level+1]); CHKERRQ(ierr);
302ccaff030SJeremy L Thompson   }
303ccaff030SJeremy L Thompson 
304ccaff030SJeremy L Thompson   // Performance logging
305ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
306ccaff030SJeremy L Thompson 
307ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
308ccaff030SJeremy L Thompson   // Setup global forcing vector
309ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
310ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
311ccaff030SJeremy L Thompson 
312ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
3136a6c615bSJeremy L Thompson     CeedVectorTakeArray(forceCeed, appCtx->memTypeRequested, NULL);
31462e9c006SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
315ccaff030SJeremy L Thompson       ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr);
31662e9c006SJeremy L Thompson     } else {
31762e9c006SJeremy L Thompson       ierr = VecCUDARestoreArray(Floc, &f); CHKERRQ(ierr);
31862e9c006SJeremy L Thompson     }
319ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F);
320ccaff030SJeremy L Thompson     CHKERRQ(ierr);
321ccaff030SJeremy L Thompson     CeedVectorDestroy(&forceCeed);
322ccaff030SJeremy L Thompson   }
323ccaff030SJeremy L Thompson 
324ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
325ccaff030SJeremy L Thompson   // Print problem summary
326ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
327ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
328ccaff030SJeremy L Thompson     const char *usedresource;
329ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
330ccaff030SJeremy L Thompson 
331ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
33217fd649aSJeremy L Thompson                        "\n-- Elasticity Example - libCEED + PETSc --\n"
333ccaff030SJeremy L Thompson                        "  libCEED:\n"
33462e9c006SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
33562e9c006SJeremy L Thompson                        "    libCEED Backend MemType            : %s\n"
33662e9c006SJeremy L Thompson                        "    libCEED User Requested MemType     : %s\n",
33762e9c006SJeremy L Thompson                        usedresource, CeedMemTypes[memTypeBackend],
33862e9c006SJeremy L Thompson                        (appCtx->setMemTypeRequest) ?
33962e9c006SJeremy L Thompson                        CeedMemTypes[appCtx->memTypeRequested] : "none");
34062e9c006SJeremy L Thompson     CHKERRQ(ierr);
341ccaff030SJeremy L Thompson 
34262e9c006SJeremy L Thompson     VecType vecType;
34362e9c006SJeremy L Thompson     ierr = VecGetType(U, &vecType); CHKERRQ(ierr);
34462e9c006SJeremy L Thompson     ierr = PetscPrintf(comm,
34562e9c006SJeremy L Thompson                        "  PETSc:\n"
34662e9c006SJeremy L Thompson                        "    PETSc Vec Type                     : %s\n",
34762e9c006SJeremy L Thompson                        vecType); CHKERRQ(ierr);
34862e9c006SJeremy L Thompson 
349ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
350ccaff030SJeremy L Thompson                        "  Problem:\n"
351ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
352ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
353ccaff030SJeremy L Thompson                        "  Mesh:\n"
354ccaff030SJeremy L Thompson                        "    File                               : %s\n"
355ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
356ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
357ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
358ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
359ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
360ccaff030SJeremy L Thompson                        "  Multigrid:\n"
361ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
362ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
363ccaff030SJeremy L Thompson                        problemTypesForDisp[appCtx->problemChoice],
364ccaff030SJeremy L Thompson                        forcingTypesForDisp[appCtx->forcingChoice],
365ccaff030SJeremy L Thompson                        appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh",
366ccaff030SJeremy L Thompson                        appCtx->degree + 1, appCtx->degree + 1,
367ccaff030SJeremy L Thompson                        Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu,
368f892e6d1Sjeremylt                        (appCtx->degree == 1 &&
369f892e6d1Sjeremylt                         appCtx->multigridChoice != MULTIGRID_NONE) ?
370f892e6d1Sjeremylt                        "Algebraic multigrid" :
371ccaff030SJeremy L Thompson                        multigridTypesForDisp[appCtx->multigridChoice],
372f892e6d1Sjeremylt                        (appCtx->degree == 1 ||
373f892e6d1Sjeremylt                         appCtx->multigridChoice == MULTIGRID_NONE) ?
374f892e6d1Sjeremylt                        0 : numLevels); CHKERRQ(ierr);
375ccaff030SJeremy L Thompson 
376ccaff030SJeremy L Thompson     if (appCtx->multigridChoice != MULTIGRID_NONE) {
377e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
378ccaff030SJeremy L Thompson         CeedInt level = i ? fineLevel : 0;
379ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
380ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
381ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
382ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
383ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
384ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
385ccaff030SJeremy L Thompson                            appCtx->levelDegrees[level] + 1,
386ccaff030SJeremy L Thompson                            Ugsz[level]/ncompu, Ulsz[level]/ncompu);
387ccaff030SJeremy L Thompson         CHKERRQ(ierr);
388ccaff030SJeremy L Thompson       }
389ccaff030SJeremy L Thompson     }
390ccaff030SJeremy L Thompson   }
391ccaff030SJeremy L Thompson 
392ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
393ccaff030SJeremy L Thompson   // Setup SNES
394ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
395ccaff030SJeremy L Thompson   // Performance logging
396ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup);
397ccaff030SJeremy L Thompson   CHKERRQ(ierr);
398ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr);
399ccaff030SJeremy L Thompson 
400ccaff030SJeremy L Thompson   // Create SNES
401ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
402ccaff030SJeremy L Thompson   ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr);
403ccaff030SJeremy L Thompson 
404ccaff030SJeremy L Thompson   // -- Jacobian evaluators
405ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr);
406ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr);
407e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
408ccaff030SJeremy L Thompson     // -- Jacobian context for level
409ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson     ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level],
411f7b4142eSJeremy L Thompson                             Uloc[level], ceedData[level], ceed, phys,
412f7b4142eSJeremy L Thompson                             physSmoother, jacobCtx[level]); CHKERRQ(ierr);
413ccaff030SJeremy L Thompson 
414ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
415ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level],
416ccaff030SJeremy L Thompson                           Ugsz[level], jacobCtx[level], &jacobMat[level]);
417ccaff030SJeremy L Thompson     CHKERRQ(ierr);
418ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT,
419ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
420ccaff030SJeremy L Thompson     CHKERRQ(ierr);
421ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL,
422ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
423a3658badSJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_DEVICE) {
424a3658badSJeremy L Thompson       ierr = MatShellSetVecType(jacobMat[level], VECCUDA); CHKERRQ(ierr);
425a3658badSJeremy L Thompson     }
426ccaff030SJeremy L Thompson   }
427e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
428ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
429ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr);
430ccaff030SJeremy L Thompson   formJacobCtx->jacobCtx = jacobCtx;
431ccaff030SJeremy L Thompson   formJacobCtx->numLevels = numLevels;
432ccaff030SJeremy L Thompson   formJacobCtx->jacobMat = jacobMat;
433ccaff030SJeremy L Thompson 
434ccaff030SJeremy L Thompson   // -- Residual evaluation function
435ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr);
436ccaff030SJeremy L Thompson   ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel],
437ccaff030SJeremy L Thompson                      sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr);
438ccaff030SJeremy L Thompson   resCtx->op = ceedData[fineLevel]->opApply;
439f7b4142eSJeremy L Thompson   resCtx->qf = ceedData[fineLevel]->qfApply;
440ccaff030SJeremy L Thompson   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson 
442ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
443ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr);
444ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr);
445e3e3df41Sjeremylt   for (PetscInt level = 1; level < numLevels; level++) {
446ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
447ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr);
44862e9c006SJeremy L Thompson     ierr = SetupProlongRestrictCtx(comm, appCtx, levelDMs[level-1],
44962e9c006SJeremy L Thompson                                    levelDMs[level], Ug[level], Uloc[level-1],
45062e9c006SJeremy L Thompson                                    Uloc[level], ceedData[level-1],
45162e9c006SJeremy L Thompson                                    ceedData[level], ceed,
452ccaff030SJeremy L Thompson                                    prolongRestrCtx[level]); CHKERRQ(ierr);
453ccaff030SJeremy L Thompson 
454ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
455ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level],
456ccaff030SJeremy L Thompson                           Ugsz[level-1], prolongRestrCtx[level],
457ccaff030SJeremy L Thompson                           &prolongRestrMat[level]); CHKERRQ(ierr);
458ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
459ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT,
460ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
461ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE,
462ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
463ccaff030SJeremy L Thompson     CHKERRQ(ierr);
464a3658badSJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_DEVICE) {
465a3658badSJeremy L Thompson       ierr = MatShellSetVecType(prolongRestrMat[level], VECCUDA); CHKERRQ(ierr);
466a3658badSJeremy L Thompson     }
467ccaff030SJeremy L Thompson   }
468ccaff030SJeremy L Thompson 
469ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
470ccaff030SJeremy L Thompson   // Setup dummy SNES for AMG coarse solve
471ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
472e3e3df41Sjeremylt   if (appCtx->multigridChoice != MULTIGRID_NONE) {
473e3e3df41Sjeremylt     // -- Jacobian Matrix
474e3e3df41Sjeremylt     ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr);
475e3e3df41Sjeremylt     ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr);
476e3e3df41Sjeremylt 
477e3e3df41Sjeremylt     if (appCtx->degree > 1) {
478ccaff030SJeremy L Thompson       ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr);
479ccaff030SJeremy L Thompson       ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr);
480ccaff030SJeremy L Thompson       ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson 
482e3e3df41Sjeremylt       // -- Jacobian function
483ccaff030SJeremy L Thompson       ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL,
484ccaff030SJeremy L Thompson                              NULL); CHKERRQ(ierr);
485ccaff030SJeremy L Thompson 
486ccaff030SJeremy L Thompson       // -- Residual evaluation function
487ccaff030SJeremy L Thompson       ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr);
488ccaff030SJeremy L Thompson       ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0]));
489ccaff030SJeremy L Thompson       CHKERRQ(ierr);
490ccaff030SJeremy L Thompson       ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed,
491ccaff030SJeremy L Thompson                              jacobCoarseCtx); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson 
493ccaff030SJeremy L Thompson       // -- Update formJacobCtx
494ccaff030SJeremy L Thompson       formJacobCtx->Ucoarse = Ug[0];
495ccaff030SJeremy L Thompson       formJacobCtx->snesCoarse = snesCoarse;
496ccaff030SJeremy L Thompson       formJacobCtx->jacobMatCoarse = jacobMatCoarse;
497e3e3df41Sjeremylt     }
498e3e3df41Sjeremylt   }
499e3e3df41Sjeremylt 
500e3e3df41Sjeremylt   // Set Jacobian function
501e3e3df41Sjeremylt   if (appCtx->degree > 1) {
502e3e3df41Sjeremylt     ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel],
503e3e3df41Sjeremylt                            FormJacobian, formJacobCtx); CHKERRQ(ierr);
504e3e3df41Sjeremylt   } else {
505e3e3df41Sjeremylt     ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse,
506e3e3df41Sjeremylt                            SNESComputeJacobianDefaultColor, NULL);
507e3e3df41Sjeremylt     CHKERRQ(ierr);
508e3e3df41Sjeremylt   }
509ccaff030SJeremy L Thompson 
510ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
511ccaff030SJeremy L Thompson   // Setup KSP
512ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
513ccaff030SJeremy L Thompson   {
514ccaff030SJeremy L Thompson     PC pc;
515ccaff030SJeremy L Thompson     KSP ksp;
516ccaff030SJeremy L Thompson 
517ccaff030SJeremy L Thompson     // -- KSP
518ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
520ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
521ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
522ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
523ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
524ccaff030SJeremy L Thompson 
525ccaff030SJeremy L Thompson     // -- Preconditioning
526ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
527ccaff030SJeremy L Thompson     ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr);
528ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
529ccaff030SJeremy L Thompson 
530ccaff030SJeremy L Thompson     if (appCtx->multigridChoice == MULTIGRID_NONE) {
531ccaff030SJeremy L Thompson       // ---- No Multigrid
532ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
534f892e6d1Sjeremylt     } else if (appCtx->degree == 1) {
535f892e6d1Sjeremylt       // ---- AMG for degree 1
536f892e6d1Sjeremylt       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson     } else {
538ccaff030SJeremy L Thompson       // ---- PCMG
539ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson 
541ccaff030SJeremy L Thompson       // ------ PCMG levels
542ccaff030SJeremy L Thompson       ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr);
543e3e3df41Sjeremylt       for (PetscInt level = 0; level < numLevels; level++) {
544ccaff030SJeremy L Thompson         // -------- Smoother
545ccaff030SJeremy L Thompson         KSP kspSmoother, kspEst;
546ccaff030SJeremy L Thompson         PC pcSmoother;
547ccaff030SJeremy L Thompson 
548ccaff030SJeremy L Thompson         // ---------- Smoother KSP
549ccaff030SJeremy L Thompson         ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr);
550ccaff030SJeremy L Thompson         ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson         ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson         // ---------- Chebyshev options
554ccaff030SJeremy L Thompson         ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr);
555ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1);
556ccaff030SJeremy L Thompson         CHKERRQ(ierr);
557ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson         ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE);
560ccaff030SJeremy L Thompson         CHKERRQ(ierr);
561ccaff030SJeremy L Thompson         ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]);
562ccaff030SJeremy L Thompson         CHKERRQ(ierr);
563ccaff030SJeremy L Thompson 
564ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
565ccaff030SJeremy L Thompson         ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson         ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr);
567ccaff030SJeremy L Thompson         ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
568ccaff030SJeremy L Thompson 
569ccaff030SJeremy L Thompson         // -------- Work vector
570ccaff030SJeremy L Thompson         if (level != fineLevel) {
571ccaff030SJeremy L Thompson           ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr);
572ccaff030SJeremy L Thompson         }
573ccaff030SJeremy L Thompson 
574ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
575ccaff030SJeremy L Thompson         if (level > 0) {
576ccaff030SJeremy L Thompson           ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]);
577ccaff030SJeremy L Thompson           CHKERRQ(ierr);
578ccaff030SJeremy L Thompson           ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]);
579ccaff030SJeremy L Thompson           CHKERRQ(ierr);
580ccaff030SJeremy L Thompson         }
581ccaff030SJeremy L Thompson       }
582ccaff030SJeremy L Thompson 
583ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
584ccaff030SJeremy L Thompson       KSP kspCoarse;
585ccaff030SJeremy L Thompson       PC pcCoarse;
586ccaff030SJeremy L Thompson 
587ccaff030SJeremy L Thompson       // -------- Coarse KSP
588ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
589ccaff030SJeremy L Thompson       ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson       ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse);
591ccaff030SJeremy L Thompson       CHKERRQ(ierr);
592ccaff030SJeremy L Thompson       ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr);
593ccaff030SJeremy L Thompson 
594ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
595ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
596ccaff030SJeremy L Thompson       ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr);
597ccaff030SJeremy L Thompson       ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr);
598ccaff030SJeremy L Thompson 
599ccaff030SJeremy L Thompson       ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr);
600ccaff030SJeremy L Thompson       ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr);
601ccaff030SJeremy L Thompson 
602ccaff030SJeremy L Thompson       // ------ PCMG options
603ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
604ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
605ccaff030SJeremy L Thompson       ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr);
606ccaff030SJeremy L Thompson     }
607ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
608ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
609ccaff030SJeremy L Thompson   }
610038d0cd7Sjeremylt   {
611038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
612481a4840SJed Brown     SNESLineSearch linesearch;
613481a4840SJed Brown 
614481a4840SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr);
615481a4840SJed Brown     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr);
616481a4840SJed Brown   }
617481a4840SJed Brown 
618ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
619ccaff030SJeremy L Thompson 
620ccaff030SJeremy L Thompson   // Performance logging
621ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
622ccaff030SJeremy L Thompson 
623ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
624ccaff030SJeremy L Thompson   // Set initial guess
625ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
626f81c27eaSJed Brown   ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr);
627ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
628ccaff030SJeremy L Thompson 
629ccaff030SJeremy L Thompson   // View solution
630ccaff030SJeremy L Thompson   if (appCtx->viewSoln) {
631ccaff030SJeremy L Thompson     ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr);
632ccaff030SJeremy L Thompson   }
633ccaff030SJeremy L Thompson 
634ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
635ccaff030SJeremy L Thompson   // Solve SNES
636ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
6375f24f028Sjeremylt   PetscBool snesMonitor = PETSC_FALSE;
6385f24f028Sjeremylt   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor);
6395f24f028Sjeremylt   CHKERRQ(ierr);
6405f24f028Sjeremylt 
641ccaff030SJeremy L Thompson   // Performance logging
642ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve);
643ccaff030SJeremy L Thompson   CHKERRQ(ierr);
644ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson 
646ccaff030SJeremy L Thompson   // Timing
647ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
648ccaff030SJeremy L Thompson   startTime = MPI_Wtime();
649ccaff030SJeremy L Thompson 
650ccaff030SJeremy L Thompson   // Solve for each load increment
6515f24f028Sjeremylt   PetscInt increment;
6525f24f028Sjeremylt   for (increment = 1; increment <= appCtx->numIncrements; increment++) {
6535f24f028Sjeremylt     // -- Log increment count
6545f24f028Sjeremylt     if (snesMonitor) {
6555f24f028Sjeremylt       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
6565f24f028Sjeremylt       CHKERRQ(ierr);
6575f24f028Sjeremylt     }
6585f24f028Sjeremylt 
659ccaff030SJeremy L Thompson     // -- Scale the problem
660ccaff030SJeremy L Thompson     PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements,
661ccaff030SJeremy L Thompson                 scalingFactor = loadIncrement /
662ccaff030SJeremy L Thompson                                 (increment == 1 ? 1 : resCtx->loadIncrement);
663ccaff030SJeremy L Thompson     resCtx->loadIncrement = loadIncrement;
664ccaff030SJeremy L Thompson     if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) {
665ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
666ccaff030SJeremy L Thompson     }
667ccaff030SJeremy L Thompson 
668ccaff030SJeremy L Thompson     // -- Solve
669ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
670ccaff030SJeremy L Thompson 
671ccaff030SJeremy L Thompson     // -- View solution
672560e6f00SJeremy L Thompson     if (appCtx->viewSoln) {
673ccaff030SJeremy L Thompson       ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr);
674ccaff030SJeremy L Thompson     }
675ccaff030SJeremy L Thompson 
676ccaff030SJeremy L Thompson     // -- Update SNES iteration count
677ccaff030SJeremy L Thompson     PetscInt its;
678ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
679ccaff030SJeremy L Thompson     snesIts += its;
6807418ca88SJeremy L Thompson     ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr);
6817418ca88SJeremy L Thompson     kspIts += its;
6823951731eSjeremylt 
6833951731eSjeremylt     // -- Check for divergence
6843951731eSjeremylt     SNESConvergedReason reason;
6853951731eSjeremylt     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
6863951731eSjeremylt     if (reason < 0)
6873951731eSjeremylt       break;
688ccaff030SJeremy L Thompson   }
689ccaff030SJeremy L Thompson 
690ccaff030SJeremy L Thompson   // Timing
691ccaff030SJeremy L Thompson   elapsedTime = MPI_Wtime() - startTime;
692ccaff030SJeremy L Thompson 
693ccaff030SJeremy L Thompson   // Performance logging
694ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
695ccaff030SJeremy L Thompson 
696ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
697ccaff030SJeremy L Thompson   // Output summary
698ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
699ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
700ccaff030SJeremy L Thompson     // -- SNES
701ccaff030SJeremy L Thompson     SNESType snesType;
702ccaff030SJeremy L Thompson     SNESConvergedReason reason;
703ccaff030SJeremy L Thompson     PetscReal rnorm;
704ccaff030SJeremy L Thompson     ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr);
705ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
706ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
707ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
708ccaff030SJeremy L Thompson                        "  SNES:\n"
709ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
710ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
711ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
7125f24f028Sjeremylt                        "    Completed Load Increments          : %d\n"
713ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
714ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
715ccaff030SJeremy L Thompson                        snesType, SNESConvergedReasons[reason],
7165f24f028Sjeremylt                        appCtx->numIncrements, increment - 1,
7175f24f028Sjeremylt                        snesIts, (double)rnorm); CHKERRQ(ierr);
718ccaff030SJeremy L Thompson 
719ccaff030SJeremy L Thompson     // -- KSP
720ccaff030SJeremy L Thompson     KSP ksp;
721ccaff030SJeremy L Thompson     KSPType kspType;
722ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
723ccaff030SJeremy L Thompson     ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr);
724ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
725ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
7267418ca88SJeremy L Thompson                        "    KSP Type                           : %s\n"
7277418ca88SJeremy L Thompson                        "    Total KSP Iterations               : %D\n",
7287418ca88SJeremy L Thompson                        kspType, kspIts); CHKERRQ(ierr);
729ccaff030SJeremy L Thompson 
730ccaff030SJeremy L Thompson     // -- PC
731ccaff030SJeremy L Thompson     PC pc;
732e3e3df41Sjeremylt     PCType pcType;
733ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
734e3e3df41Sjeremylt     ierr = PCGetType(pc, &pcType); CHKERRQ(ierr);
735e3e3df41Sjeremylt     ierr = PetscPrintf(comm,
736e3e3df41Sjeremylt                        "    PC Type                            : %s\n",
737e3e3df41Sjeremylt                        pcType); CHKERRQ(ierr);
738e3e3df41Sjeremylt 
7392b901a5eSjeremylt     if (!strcmp(pcType, PCMG)) {
740e3e3df41Sjeremylt       PCMGType pcmgType;
741ccaff030SJeremy L Thompson       ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr);
742ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
743ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
744ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
745ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
746ccaff030SJeremy L Thompson                          PCMGTypes[pcmgType],
747ccaff030SJeremy L Thompson                          PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr);
748ccaff030SJeremy L Thompson 
749ccaff030SJeremy L Thompson       // -- Coarse Solve
750ccaff030SJeremy L Thompson       KSP kspCoarse;
751ccaff030SJeremy L Thompson       PC pcCoarse;
752ccaff030SJeremy L Thompson       PCType pcType;
753ccaff030SJeremy L Thompson 
754ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
755ccaff030SJeremy L Thompson       ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr);
756ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
757ccaff030SJeremy L Thompson       ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr);
758ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
759ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
760ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
761ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
762ccaff030SJeremy L Thompson                          kspType, pcType); CHKERRQ(ierr);
763ccaff030SJeremy L Thompson     }
764ccaff030SJeremy L Thompson   }
765ccaff030SJeremy L Thompson 
766ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
767ccaff030SJeremy L Thompson   // Compute solve time
768ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
769ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
770ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm);
771ccaff030SJeremy L Thompson     CHKERRQ(ierr);
772ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm);
773ccaff030SJeremy L Thompson     CHKERRQ(ierr);
774ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
775ccaff030SJeremy L Thompson                        "  Performance:\n"
7767418ca88SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n"
7777418ca88SJeremy L Thompson                        "    DoFs/Sec in SNES                   : %g (%g) million\n",
7787418ca88SJeremy L Thompson                        maxTime, minTime, 1e-6*Ugsz[fineLevel]*kspIts/maxTime,
7797418ca88SJeremy L Thompson                        1e-6*Ugsz[fineLevel]*kspIts/minTime); CHKERRQ(ierr);
780ccaff030SJeremy L Thompson   }
781ccaff030SJeremy L Thompson 
782ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
783ccaff030SJeremy L Thompson   // Compute error
784ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
785ccaff030SJeremy L Thompson   if (appCtx->forcingChoice == FORCE_MMS) {
786ccaff030SJeremy L Thompson     CeedScalar l2Error = 1., l2Unorm = 1.;
787ccaff030SJeremy L Thompson     const CeedScalar *truearray;
788ccaff030SJeremy L Thompson     Vec errorVec, trueVec;
789ccaff030SJeremy L Thompson 
790ccaff030SJeremy L Thompson     // -- Work vectors
791ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr);
792ccaff030SJeremy L Thompson     ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr);
793ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr);
794ccaff030SJeremy L Thompson     ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr);
795ccaff030SJeremy L Thompson 
796ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
79762e9c006SJeremy L Thompson     CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln,
79862e9c006SJeremy L Thompson                            appCtx->memTypeRequested, &truearray);
79962e9c006SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
80062e9c006SJeremy L Thompson       ierr = VecPlaceArray(resCtx->Yloc, (PetscScalar *)truearray);
80162e9c006SJeremy L Thompson       CHKERRQ(ierr);
80262e9c006SJeremy L Thompson     } else {
80362e9c006SJeremy L Thompson       ierr = VecCUDAPlaceArray(resCtx->Yloc, (PetscScalar *)truearray);
80462e9c006SJeremy L Thompson       CHKERRQ(ierr);
80562e9c006SJeremy L Thompson     }
806ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec);
807ccaff030SJeremy L Thompson     CHKERRQ(ierr);
80862e9c006SJeremy L Thompson     if (appCtx->memTypeRequested == CEED_MEM_HOST) {
809ccaff030SJeremy L Thompson       ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr);
81062e9c006SJeremy L Thompson     } else {
81162e9c006SJeremy L Thompson       ierr = VecCUDAResetArray(resCtx->Yloc); CHKERRQ(ierr);
81262e9c006SJeremy L Thompson     }
813ccaff030SJeremy L Thompson     CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray);
814ccaff030SJeremy L Thompson 
815ccaff030SJeremy L Thompson     // -- Compute L2 error
816ccaff030SJeremy L Thompson     ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr);
817ccaff030SJeremy L Thompson     ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr);
818ccaff030SJeremy L Thompson     ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr);
819ccaff030SJeremy L Thompson     l2Error /= l2Unorm;
820ccaff030SJeremy L Thompson 
821ccaff030SJeremy L Thompson     // -- Output
822ccaff030SJeremy L Thompson     if (!appCtx->testMode || l2Error > 0.05) {
823ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
824ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
825ccaff030SJeremy L Thompson                          l2Error); CHKERRQ(ierr);
826ccaff030SJeremy L Thompson     }
827ccaff030SJeremy L Thompson 
828ccaff030SJeremy L Thompson     // -- Cleanup
829ccaff030SJeremy L Thompson     ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
830ccaff030SJeremy L Thompson     ierr = VecDestroy(&trueVec); CHKERRQ(ierr);
8312d93065eSjeremylt   }
8322d93065eSjeremylt 
8332d93065eSjeremylt   // ---------------------------------------------------------------------------
8342d93065eSjeremylt   // Compute energy
8352d93065eSjeremylt   // ---------------------------------------------------------------------------
8362d93065eSjeremylt   if (!appCtx->testMode) {
8372d93065eSjeremylt     // -- Compute L2 error
8382d93065eSjeremylt     CeedScalar energy;
839a3c02c40SJeremy L Thompson     ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy,
840a3c02c40SJeremy L Thompson                                U, &energy); CHKERRQ(ierr);
8412d93065eSjeremylt 
8422d93065eSjeremylt     // -- Output
8432d93065eSjeremylt     ierr = PetscPrintf(comm,
8442d93065eSjeremylt                        "    Strain Energy                      : %e\n",
8452d93065eSjeremylt                        energy); CHKERRQ(ierr);
846ccaff030SJeremy L Thompson   }
847ccaff030SJeremy L Thompson 
848ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
8495c25879aSJeremy L Thompson   // Output diagnostic quantities
8505c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
8515c25879aSJeremy L Thompson   if (appCtx->viewSoln || appCtx->viewFinalSoln) {
8525c25879aSJeremy L Thompson     // -- Setup context
8535c25879aSJeremy L Thompson     UserMult diagnosticCtx;
8545c25879aSJeremy L Thompson     ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr);
8555c25879aSJeremy L Thompson     ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr);
8565c25879aSJeremy L Thompson     diagnosticCtx->dm = dmDiagnostic;
8575c25879aSJeremy L Thompson     diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic;
8585c25879aSJeremy L Thompson 
8595c25879aSJeremy L Thompson     // -- Compute and output
8605c25879aSJeremy L Thompson     ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U,
8615c25879aSJeremy L Thompson                                     ceedData[fineLevel]->ErestrictDiagnostic);
8625c25879aSJeremy L Thompson     CHKERRQ(ierr);
8635c25879aSJeremy L Thompson 
8645c25879aSJeremy L Thompson     // -- Cleanup
8655c25879aSJeremy L Thompson     ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr);
8665c25879aSJeremy L Thompson   }
8675c25879aSJeremy L Thompson 
8685c25879aSJeremy L Thompson   // ---------------------------------------------------------------------------
869ccaff030SJeremy L Thompson   // Free objects
870ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
871ccaff030SJeremy L Thompson   // Data in arrays per level
872e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
873ccaff030SJeremy L Thompson     // Vectors
874ccaff030SJeremy L Thompson     ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr);
875ccaff030SJeremy L Thompson     ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr);
876ccaff030SJeremy L Thompson 
877ccaff030SJeremy L Thompson     // Jacobian matrix and data
878ccaff030SJeremy L Thompson     ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr);
879ccaff030SJeremy L Thompson     ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr);
880ccaff030SJeremy L Thompson     ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr);
881ccaff030SJeremy L Thompson 
882ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
883ccaff030SJeremy L Thompson     if (level > 0) {
884ccaff030SJeremy L Thompson       ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr);
885ccaff030SJeremy L Thompson       ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr);
886ccaff030SJeremy L Thompson     }
887ccaff030SJeremy L Thompson 
888ccaff030SJeremy L Thompson     // DM
889ccaff030SJeremy L Thompson     ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr);
890ccaff030SJeremy L Thompson 
891ccaff030SJeremy L Thompson     // libCEED objects
892ccaff030SJeremy L Thompson     ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr);
893ccaff030SJeremy L Thompson   }
894ccaff030SJeremy L Thompson 
895ccaff030SJeremy L Thompson   // Arrays
896ccaff030SJeremy L Thompson   ierr = PetscFree(Ug); CHKERRQ(ierr);
897ccaff030SJeremy L Thompson   ierr = PetscFree(Uloc); CHKERRQ(ierr);
898ccaff030SJeremy L Thompson   ierr = PetscFree(Ugsz); CHKERRQ(ierr);
899ccaff030SJeremy L Thompson   ierr = PetscFree(Ulsz); CHKERRQ(ierr);
900ccaff030SJeremy L Thompson   ierr = PetscFree(Ulocsz); CHKERRQ(ierr);
901ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCtx); CHKERRQ(ierr);
902ccaff030SJeremy L Thompson   ierr = PetscFree(jacobMat); CHKERRQ(ierr);
903ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr);
904ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr);
905ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr);
906ccaff030SJeremy L Thompson   ierr = PetscFree(ceedData); CHKERRQ(ierr);
907ccaff030SJeremy L Thompson 
908ccaff030SJeremy L Thompson   // libCEED objects
909ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfRestrict);
910ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfProlong);
911ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
912ccaff030SJeremy L Thompson 
913ccaff030SJeremy L Thompson   // PETSc objects
914ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
915ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
916ccaff030SJeremy L Thompson   ierr = VecDestroy(&Rloc); CHKERRQ(ierr);
917ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
918ccaff030SJeremy L Thompson   ierr = VecDestroy(&Floc); CHKERRQ(ierr);
919ccaff030SJeremy L Thompson   ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr);
920ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
921ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr);
922ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
923a3c02c40SJeremy L Thompson   ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr);
924a3c02c40SJeremy L Thompson   ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr);
925ccaff030SJeremy L Thompson   ierr = PetscFree(levelDMs); CHKERRQ(ierr);
926ccaff030SJeremy L Thompson 
927ccaff030SJeremy L Thompson   // Structs
928ccaff030SJeremy L Thompson   ierr = PetscFree(resCtx); CHKERRQ(ierr);
929ccaff030SJeremy L Thompson   ierr = PetscFree(formJacobCtx); CHKERRQ(ierr);
930ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr);
931ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx); CHKERRQ(ierr);
932ccaff030SJeremy L Thompson   ierr = PetscFree(phys); CHKERRQ(ierr);
933f7b4142eSJeremy L Thompson   ierr = PetscFree(physSmoother); CHKERRQ(ierr);
934ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
935ccaff030SJeremy L Thompson 
936ccaff030SJeremy L Thompson   return PetscFinalize();
937ccaff030SJeremy L Thompson }
938