xref: /libCEED/examples/solids/elasticity.c (revision ccaff0309dc399f656ea11018b919b30feb8b0fa)
1*ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2*ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3*ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4*ccaff030SJeremy L Thompson //
5*ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6*ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7*ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8*ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9*ccaff030SJeremy L Thompson //
10*ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12*ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13*ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14*ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15*ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16*ccaff030SJeremy L Thompson 
17*ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Elasticity
18*ccaff030SJeremy L Thompson //
19*ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve
20*ccaff030SJeremy L Thompson //   elasticity problems.
21*ccaff030SJeremy L Thompson //
22*ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex.
23*ccaff030SJeremy L Thompson //
24*ccaff030SJeremy L Thompson // Build with:
25*ccaff030SJeremy L Thompson //
26*ccaff030SJeremy L Thompson //     make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
27*ccaff030SJeremy L Thompson //
28*ccaff030SJeremy L Thompson // Sample runs:
29*ccaff030SJeremy L Thompson //
30*ccaff030SJeremy L Thompson //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem linElas -forcing mms
31*ccaff030SJeremy L Thompson //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperSS -forcing none -ceed /cpu/self
32*ccaff030SJeremy L Thompson //     ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperFS -forcing none -ceed /gpu/occa
33*ccaff030SJeremy L Thompson //
34*ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
35*ccaff030SJeremy L Thompson //
36*ccaff030SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -degree 2 -nu 0.3 -E 1
37*ccaff030SJeremy L Thompson 
38*ccaff030SJeremy L Thompson /// @file
39*ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex
40*ccaff030SJeremy L Thompson 
41*ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n";
42*ccaff030SJeremy L Thompson 
43*ccaff030SJeremy L Thompson #include "elasticity.h"
44*ccaff030SJeremy L Thompson 
45*ccaff030SJeremy L Thompson int main(int argc, char **argv) {
46*ccaff030SJeremy L Thompson   PetscInt       ierr;
47*ccaff030SJeremy L Thompson   MPI_Comm       comm;
48*ccaff030SJeremy L Thompson   // Context structs
49*ccaff030SJeremy L Thompson   AppCtx         appCtx;                 // Contains problem options
50*ccaff030SJeremy L Thompson   Physics        phys;                   // Contains physical constants
51*ccaff030SJeremy L Thompson   Units          units;                  // Contains units scaling
52*ccaff030SJeremy L Thompson   // PETSc objects
53*ccaff030SJeremy L Thompson   PetscLogStage  stageDMSetup, stageLibceedSetup,
54*ccaff030SJeremy L Thompson                  stageSnesSetup, stageSnesSolve;
55*ccaff030SJeremy L Thompson   DM             dmOrig;                 // Distributed DM to clone
56*ccaff030SJeremy L Thompson   DM             *levelDMs;
57*ccaff030SJeremy L Thompson   Vec            U, *Ug, *Uloc;          // U: solution, R: residual, F: forcing
58*ccaff030SJeremy L Thompson   Vec            R, Rloc, F, Floc;       // g: global, loc: local
59*ccaff030SJeremy L Thompson   SNES           snes, snesCoarse;
60*ccaff030SJeremy L Thompson   Mat            *jacobMat, jacobMatCoarse, *prolongRestrMat;
61*ccaff030SJeremy L Thompson   // PETSc data
62*ccaff030SJeremy L Thompson   UserMult       resCtx, jacobCoarseCtx, *jacobCtx;
63*ccaff030SJeremy L Thompson   FormJacobCtx   formJacobCtx;
64*ccaff030SJeremy L Thompson   UserMultProlongRestr *prolongRestrCtx;
65*ccaff030SJeremy L Thompson   PCMGCycleType  pcmgCycleType = PC_MG_CYCLE_V;
66*ccaff030SJeremy L Thompson   // libCEED objects
67*ccaff030SJeremy L Thompson   Ceed           ceed, ceedFine = NULL;
68*ccaff030SJeremy L Thompson   CeedData       *ceedData;
69*ccaff030SJeremy L Thompson   CeedQFunction  qfRestrict, qfProlong;
70*ccaff030SJeremy L Thompson   // Parameters
71*ccaff030SJeremy L Thompson   PetscInt       ncompu = 3;             // 3 DoFs in 3D
72*ccaff030SJeremy L Thompson   PetscInt       numLevels = 1, fineLevel = 0;
73*ccaff030SJeremy L Thompson   PetscInt       *Ugsz, *Ulsz, *Ulocsz;  // sz: size
74*ccaff030SJeremy L Thompson   PetscInt       snesIts = 0;
75*ccaff030SJeremy L Thompson   // Timing
76*ccaff030SJeremy L Thompson   double         startTime, elapsedTime, minTime, maxTime;
77*ccaff030SJeremy L Thompson 
78*ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
79*ccaff030SJeremy L Thompson   if (ierr)
80*ccaff030SJeremy L Thompson     return ierr;
81*ccaff030SJeremy L Thompson 
82*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
83*ccaff030SJeremy L Thompson   // Process command line options
84*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
85*ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
86*ccaff030SJeremy L Thompson 
87*ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
88*ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr);
89*ccaff030SJeremy L Thompson   ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr);
90*ccaff030SJeremy L Thompson   numLevels = appCtx->numLevels;
91*ccaff030SJeremy L Thompson   fineLevel = numLevels - 1;
92*ccaff030SJeremy L Thompson 
93*ccaff030SJeremy L Thompson   // -- Set Poison's ratio, Young's Modulus
94*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
95*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
96*ccaff030SJeremy L Thompson   ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr);
97*ccaff030SJeremy L Thompson 
98*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
99*ccaff030SJeremy L Thompson   // Setup DM
100*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
101*ccaff030SJeremy L Thompson   // Performance logging
102*ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup);
103*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
104*ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr);
105*ccaff030SJeremy L Thompson 
106*ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
107*ccaff030SJeremy L Thompson   ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr);
108*ccaff030SJeremy L Thompson 
109*ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
110*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr);
111*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
112*ccaff030SJeremy L Thompson     ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr);
113*ccaff030SJeremy L Thompson     ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level],
114*ccaff030SJeremy L Thompson                            ncompu); CHKERRQ(ierr);
115*ccaff030SJeremy L Thompson   }
116*ccaff030SJeremy L Thompson 
117*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
118*ccaff030SJeremy L Thompson   // Setup solution and work vectors
119*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
120*ccaff030SJeremy L Thompson   // Allocate arrays
121*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr);
122*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr);
123*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr);
124*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr);
125*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr);
126*ccaff030SJeremy L Thompson 
127*ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
128*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
129*ccaff030SJeremy L Thompson     // -- Create global unknown vector U
130*ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr);
131*ccaff030SJeremy L Thompson     ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr);
132*ccaff030SJeremy L Thompson     // Note: Local size for matShell
133*ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr);
134*ccaff030SJeremy L Thompson 
135*ccaff030SJeremy L Thompson     // -- Create local unknown vector Uloc
136*ccaff030SJeremy L Thompson     ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr);
137*ccaff030SJeremy L Thompson     // Note: local size for libCEED
138*ccaff030SJeremy L Thompson     ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr);
139*ccaff030SJeremy L Thompson   }
140*ccaff030SJeremy L Thompson 
141*ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
142*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr);
143*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr);
144*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr);
145*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr);
146*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr);
147*ccaff030SJeremy L Thompson 
148*ccaff030SJeremy L Thompson   // Performance logging
149*ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
150*ccaff030SJeremy L Thompson 
151*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
152*ccaff030SJeremy L Thompson   // Set up libCEED
153*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
154*ccaff030SJeremy L Thompson   // Performance logging
155*ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup);
156*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
157*ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr);
158*ccaff030SJeremy L Thompson 
159*ccaff030SJeremy L Thompson   // Initalize backend
160*ccaff030SJeremy L Thompson   CeedInit(appCtx->ceedResource, &ceed);
161*ccaff030SJeremy L Thompson   if (appCtx->degree > 4 && appCtx->ceedResourceFine[0])
162*ccaff030SJeremy L Thompson     CeedInit(appCtx->ceedResourceFine, &ceedFine);
163*ccaff030SJeremy L Thompson 
164*ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
165*ccaff030SJeremy L Thompson   CeedVector forceCeed;
166*ccaff030SJeremy L Thompson   CeedScalar *f;
167*ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
168*ccaff030SJeremy L Thompson     ierr = VecGetArray(Floc, &f); CHKERRQ(ierr);
169*ccaff030SJeremy L Thompson     CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed);
170*ccaff030SJeremy L Thompson     CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f);
171*ccaff030SJeremy L Thompson   }
172*ccaff030SJeremy L Thompson 
173*ccaff030SJeremy L Thompson   // -- Restriction and prolongation QFunction
174*ccaff030SJeremy L Thompson   if (appCtx->multigridChoice != MULTIGRID_NONE) {
175*ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
176*ccaff030SJeremy L Thompson                                 &qfRestrict);
177*ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
178*ccaff030SJeremy L Thompson                                 &qfProlong);
179*ccaff030SJeremy L Thompson   }
180*ccaff030SJeremy L Thompson 
181*ccaff030SJeremy L Thompson   // -- Setup libCEED objects
182*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr);
183*ccaff030SJeremy L Thompson   // ---- Setup residual evaluator and geometric information
184*ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr);
185*ccaff030SJeremy L Thompson   {
186*ccaff030SJeremy L Thompson     bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine);
187*ccaff030SJeremy L Thompson     Ceed levelCeed = highOrder ? ceedFine : ceed;
188*ccaff030SJeremy L Thompson     ierr = SetupLibceedFineLevel(levelDMs[fineLevel], levelCeed, appCtx,
189*ccaff030SJeremy L Thompson                                  phys, ceedData, fineLevel, ncompu,
190*ccaff030SJeremy L Thompson                                  Ugsz[fineLevel], Ulocsz[fineLevel], forceCeed,
191*ccaff030SJeremy L Thompson                                  qfRestrict, qfProlong);
192*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
193*ccaff030SJeremy L Thompson   }
194*ccaff030SJeremy L Thompson   // ---- Setup Jacobian evaluator and prolongation/restriction
195*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
196*ccaff030SJeremy L Thompson     if (level != fineLevel) {
197*ccaff030SJeremy L Thompson       ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr);
198*ccaff030SJeremy L Thompson     }
199*ccaff030SJeremy L Thompson 
200*ccaff030SJeremy L Thompson     // Note: use high order ceed, if specified and degree > 3
201*ccaff030SJeremy L Thompson     bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine);
202*ccaff030SJeremy L Thompson     Ceed levelCeed = highOrder ? ceedFine : ceed;
203*ccaff030SJeremy L Thompson     ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys,
204*ccaff030SJeremy L Thompson                              ceedData,  level, ncompu, Ugsz[level],
205*ccaff030SJeremy L Thompson                              Ulocsz[level], forceCeed, qfRestrict,
206*ccaff030SJeremy L Thompson                              qfProlong); CHKERRQ(ierr);
207*ccaff030SJeremy L Thompson   }
208*ccaff030SJeremy L Thompson 
209*ccaff030SJeremy L Thompson   // Performance logging
210*ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
211*ccaff030SJeremy L Thompson 
212*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
213*ccaff030SJeremy L Thompson   // Setup global forcing vector
214*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
215*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
216*ccaff030SJeremy L Thompson 
217*ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
218*ccaff030SJeremy L Thompson     ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr);
219*ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F);
220*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
221*ccaff030SJeremy L Thompson     CeedVectorDestroy(&forceCeed);
222*ccaff030SJeremy L Thompson   }
223*ccaff030SJeremy L Thompson 
224*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
225*ccaff030SJeremy L Thompson   // Print problem summary
226*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
227*ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
228*ccaff030SJeremy L Thompson     const char *usedresource;
229*ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
230*ccaff030SJeremy L Thompson 
231*ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
232*ccaff030SJeremy L Thompson                        "\n-- Elastisticy Example - libCEED + PETSc --\n"
233*ccaff030SJeremy L Thompson                        "  libCEED:\n"
234*ccaff030SJeremy L Thompson                        "    libCEED Backend                    : %s\n",
235*ccaff030SJeremy L Thompson                        usedresource); CHKERRQ(ierr);
236*ccaff030SJeremy L Thompson 
237*ccaff030SJeremy L Thompson     if (ceedFine) {
238*ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
239*ccaff030SJeremy L Thompson                          "    libCEED Backend - high order       : %s\n",
240*ccaff030SJeremy L Thompson                          appCtx->ceedResourceFine); CHKERRQ(ierr);
241*ccaff030SJeremy L Thompson     }
242*ccaff030SJeremy L Thompson 
243*ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
244*ccaff030SJeremy L Thompson                        "  Problem:\n"
245*ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
246*ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
247*ccaff030SJeremy L Thompson                        "  Mesh:\n"
248*ccaff030SJeremy L Thompson                        "    File                               : %s\n"
249*ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
250*ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
251*ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
252*ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
253*ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
254*ccaff030SJeremy L Thompson                        "  Multigrid:\n"
255*ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
256*ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
257*ccaff030SJeremy L Thompson                        problemTypesForDisp[appCtx->problemChoice],
258*ccaff030SJeremy L Thompson                        forcingTypesForDisp[appCtx->forcingChoice],
259*ccaff030SJeremy L Thompson                        appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh",
260*ccaff030SJeremy L Thompson                        appCtx->degree + 1, appCtx->degree + 1,
261*ccaff030SJeremy L Thompson                        Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu,
262*ccaff030SJeremy L Thompson                        multigridTypesForDisp[appCtx->multigridChoice],
263*ccaff030SJeremy L Thompson                        numLevels); CHKERRQ(ierr);
264*ccaff030SJeremy L Thompson 
265*ccaff030SJeremy L Thompson     if (appCtx->multigridChoice != MULTIGRID_NONE) {
266*ccaff030SJeremy L Thompson       for (int i = 0; i < 2; i++) {
267*ccaff030SJeremy L Thompson         CeedInt level = i ? fineLevel : 0;
268*ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
269*ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
270*ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
271*ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
272*ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
273*ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
274*ccaff030SJeremy L Thompson                            appCtx->levelDegrees[level] + 1,
275*ccaff030SJeremy L Thompson                            Ugsz[level]/ncompu, Ulsz[level]/ncompu);
276*ccaff030SJeremy L Thompson         CHKERRQ(ierr);
277*ccaff030SJeremy L Thompson       }
278*ccaff030SJeremy L Thompson     }
279*ccaff030SJeremy L Thompson   }
280*ccaff030SJeremy L Thompson 
281*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
282*ccaff030SJeremy L Thompson   // Setup SNES
283*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
284*ccaff030SJeremy L Thompson   // Performance logging
285*ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup);
286*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
287*ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr);
288*ccaff030SJeremy L Thompson 
289*ccaff030SJeremy L Thompson   // Create SNES
290*ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
291*ccaff030SJeremy L Thompson   ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr);
292*ccaff030SJeremy L Thompson 
293*ccaff030SJeremy L Thompson   // -- Jacobian evaluators
294*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr);
295*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr);
296*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
297*ccaff030SJeremy L Thompson     // -- Jacobian context for level
298*ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr);
299*ccaff030SJeremy L Thompson     ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level],
300*ccaff030SJeremy L Thompson                             Uloc[level], ceedData[level], ceed,
301*ccaff030SJeremy L Thompson                             jacobCtx[level]); CHKERRQ(ierr);
302*ccaff030SJeremy L Thompson 
303*ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
304*ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level],
305*ccaff030SJeremy L Thompson                           Ugsz[level], jacobCtx[level], &jacobMat[level]);
306*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
307*ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT,
308*ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
309*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
310*ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL,
311*ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
312*ccaff030SJeremy L Thompson 
313*ccaff030SJeremy L Thompson   }
314*ccaff030SJeremy L Thompson   // Note: FormJacobian updates the state count of the Jacobian diagonals
315*ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
316*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr);
317*ccaff030SJeremy L Thompson   formJacobCtx->jacobCtx = jacobCtx;
318*ccaff030SJeremy L Thompson   formJacobCtx->numLevels = numLevels;
319*ccaff030SJeremy L Thompson   formJacobCtx->jacobMat = jacobMat;
320*ccaff030SJeremy L Thompson   ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel],
321*ccaff030SJeremy L Thompson                          FormJacobian, formJacobCtx); CHKERRQ(ierr);
322*ccaff030SJeremy L Thompson 
323*ccaff030SJeremy L Thompson   // -- Residual evaluation function
324*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr);
325*ccaff030SJeremy L Thompson   ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel],
326*ccaff030SJeremy L Thompson                      sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr);
327*ccaff030SJeremy L Thompson   resCtx->op = ceedData[fineLevel]->opApply;
328*ccaff030SJeremy L Thompson   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr);
329*ccaff030SJeremy L Thompson 
330*ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
331*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr);
332*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr);
333*ccaff030SJeremy L Thompson   for (int level = 1; level < numLevels; level++) {
334*ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
335*ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr);
336*ccaff030SJeremy L Thompson     ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level],
337*ccaff030SJeremy L Thompson                                    Ug[level], Uloc[level-1], Uloc[level],
338*ccaff030SJeremy L Thompson                                    ceedData[level-1], ceedData[level], ceed,
339*ccaff030SJeremy L Thompson                                    prolongRestrCtx[level]); CHKERRQ(ierr);
340*ccaff030SJeremy L Thompson 
341*ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
342*ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level],
343*ccaff030SJeremy L Thompson                           Ugsz[level-1], prolongRestrCtx[level],
344*ccaff030SJeremy L Thompson                           &prolongRestrMat[level]); CHKERRQ(ierr);
345*ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
346*ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT,
347*ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
348*ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE,
349*ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
350*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
351*ccaff030SJeremy L Thompson   }
352*ccaff030SJeremy L Thompson 
353*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
354*ccaff030SJeremy L Thompson   // Setup dummy SNES for AMG coarse solve
355*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
356*ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr);
357*ccaff030SJeremy L Thompson   ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr);
358*ccaff030SJeremy L Thompson   ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr);
359*ccaff030SJeremy L Thompson 
360*ccaff030SJeremy L Thompson   // -- Jacobian matrix
361*ccaff030SJeremy L Thompson   ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr);
362*ccaff030SJeremy L Thompson   ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr);
363*ccaff030SJeremy L Thompson   ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL,
364*ccaff030SJeremy L Thompson                          NULL); CHKERRQ(ierr);
365*ccaff030SJeremy L Thompson 
366*ccaff030SJeremy L Thompson   // -- Residual evaluation function
367*ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr);
368*ccaff030SJeremy L Thompson   ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0]));
369*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
370*ccaff030SJeremy L Thompson   ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed,
371*ccaff030SJeremy L Thompson                          jacobCoarseCtx); CHKERRQ(ierr);
372*ccaff030SJeremy L Thompson 
373*ccaff030SJeremy L Thompson   // -- Update formJacobCtx
374*ccaff030SJeremy L Thompson   formJacobCtx->Ucoarse = Ug[0];
375*ccaff030SJeremy L Thompson   formJacobCtx->snesCoarse = snesCoarse;
376*ccaff030SJeremy L Thompson   formJacobCtx->jacobMatCoarse = jacobMatCoarse;
377*ccaff030SJeremy L Thompson 
378*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
379*ccaff030SJeremy L Thompson   // Setup KSP
380*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
381*ccaff030SJeremy L Thompson   {
382*ccaff030SJeremy L Thompson     PC pc;
383*ccaff030SJeremy L Thompson     KSP ksp;
384*ccaff030SJeremy L Thompson 
385*ccaff030SJeremy L Thompson     // -- KSP
386*ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
387*ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
388*ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
389*ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
390*ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
391*ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
392*ccaff030SJeremy L Thompson 
393*ccaff030SJeremy L Thompson     // -- Preconditioning
394*ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
395*ccaff030SJeremy L Thompson     ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr);
396*ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
397*ccaff030SJeremy L Thompson 
398*ccaff030SJeremy L Thompson     if (appCtx->multigridChoice == MULTIGRID_NONE) {
399*ccaff030SJeremy L Thompson       // ---- No Multigrid
400*ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
401*ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
402*ccaff030SJeremy L Thompson     } else {
403*ccaff030SJeremy L Thompson       // ---- PCMG
404*ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
405*ccaff030SJeremy L Thompson 
406*ccaff030SJeremy L Thompson       // ------ PCMG levels
407*ccaff030SJeremy L Thompson       ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr);
408*ccaff030SJeremy L Thompson       for (int level = 0; level < numLevels; level++) {
409*ccaff030SJeremy L Thompson         // -------- Smoother
410*ccaff030SJeremy L Thompson         KSP kspSmoother, kspEst;
411*ccaff030SJeremy L Thompson         PC pcSmoother;
412*ccaff030SJeremy L Thompson 
413*ccaff030SJeremy L Thompson         // ---------- Smoother KSP
414*ccaff030SJeremy L Thompson         ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr);
415*ccaff030SJeremy L Thompson         ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr);
416*ccaff030SJeremy L Thompson         ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr);
417*ccaff030SJeremy L Thompson 
418*ccaff030SJeremy L Thompson         // ---------- Chebyshev options
419*ccaff030SJeremy L Thompson         ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr);
420*ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1);
421*ccaff030SJeremy L Thompson         CHKERRQ(ierr);
422*ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr);
423*ccaff030SJeremy L Thompson         ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr);
424*ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE);
425*ccaff030SJeremy L Thompson         CHKERRQ(ierr);
426*ccaff030SJeremy L Thompson         ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]);
427*ccaff030SJeremy L Thompson         CHKERRQ(ierr);
428*ccaff030SJeremy L Thompson 
429*ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
430*ccaff030SJeremy L Thompson         ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr);
431*ccaff030SJeremy L Thompson         ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr);
432*ccaff030SJeremy L Thompson         ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
433*ccaff030SJeremy L Thompson 
434*ccaff030SJeremy L Thompson         // -------- Work vector
435*ccaff030SJeremy L Thompson         if (level != fineLevel) {
436*ccaff030SJeremy L Thompson           ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr);
437*ccaff030SJeremy L Thompson         }
438*ccaff030SJeremy L Thompson 
439*ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
440*ccaff030SJeremy L Thompson         if (level > 0) {
441*ccaff030SJeremy L Thompson           ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]);
442*ccaff030SJeremy L Thompson           CHKERRQ(ierr);
443*ccaff030SJeremy L Thompson           ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]);
444*ccaff030SJeremy L Thompson           CHKERRQ(ierr);
445*ccaff030SJeremy L Thompson         }
446*ccaff030SJeremy L Thompson       }
447*ccaff030SJeremy L Thompson 
448*ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
449*ccaff030SJeremy L Thompson       KSP kspCoarse;
450*ccaff030SJeremy L Thompson       PC pcCoarse;
451*ccaff030SJeremy L Thompson 
452*ccaff030SJeremy L Thompson       // -------- Coarse KSP
453*ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
454*ccaff030SJeremy L Thompson       ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr);
455*ccaff030SJeremy L Thompson       ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse);
456*ccaff030SJeremy L Thompson       CHKERRQ(ierr);
457*ccaff030SJeremy L Thompson       ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr);
458*ccaff030SJeremy L Thompson 
459*ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
460*ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
461*ccaff030SJeremy L Thompson       ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr);
462*ccaff030SJeremy L Thompson       ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr);
463*ccaff030SJeremy L Thompson 
464*ccaff030SJeremy L Thompson       ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr);
465*ccaff030SJeremy L Thompson       ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr);
466*ccaff030SJeremy L Thompson 
467*ccaff030SJeremy L Thompson       // ------ PCMG options
468*ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
469*ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
470*ccaff030SJeremy L Thompson       ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr);
471*ccaff030SJeremy L Thompson     }
472*ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
473*ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
474*ccaff030SJeremy L Thompson   }
475*ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
476*ccaff030SJeremy L Thompson 
477*ccaff030SJeremy L Thompson   // Performance logging
478*ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
479*ccaff030SJeremy L Thompson 
480*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
481*ccaff030SJeremy L Thompson   // Set initial guess
482*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
483*ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
484*ccaff030SJeremy L Thompson 
485*ccaff030SJeremy L Thompson   // View solution
486*ccaff030SJeremy L Thompson   if (appCtx->viewSoln) {
487*ccaff030SJeremy L Thompson     ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr);
488*ccaff030SJeremy L Thompson   }
489*ccaff030SJeremy L Thompson 
490*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
491*ccaff030SJeremy L Thompson   // Solve SNES
492*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
493*ccaff030SJeremy L Thompson   // Performance logging
494*ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve);
495*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
496*ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr);
497*ccaff030SJeremy L Thompson 
498*ccaff030SJeremy L Thompson   // Timing
499*ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
500*ccaff030SJeremy L Thompson   startTime = MPI_Wtime();
501*ccaff030SJeremy L Thompson 
502*ccaff030SJeremy L Thompson   // Solve for each load increment
503*ccaff030SJeremy L Thompson   for (PetscInt increment = 1; increment <= appCtx->numIncrements; increment++) {
504*ccaff030SJeremy L Thompson     // -- Scale the problem
505*ccaff030SJeremy L Thompson     PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements,
506*ccaff030SJeremy L Thompson                 scalingFactor = loadIncrement /
507*ccaff030SJeremy L Thompson                                 (increment == 1 ? 1 : resCtx->loadIncrement);
508*ccaff030SJeremy L Thompson     resCtx->loadIncrement = loadIncrement;
509*ccaff030SJeremy L Thompson     if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) {
510*ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
511*ccaff030SJeremy L Thompson     }
512*ccaff030SJeremy L Thompson 
513*ccaff030SJeremy L Thompson     // -- Solve
514*ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
515*ccaff030SJeremy L Thompson 
516*ccaff030SJeremy L Thompson     // -- View solution
517*ccaff030SJeremy L Thompson     if (appCtx->viewSoln) {
518*ccaff030SJeremy L Thompson       ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr);
519*ccaff030SJeremy L Thompson     }
520*ccaff030SJeremy L Thompson 
521*ccaff030SJeremy L Thompson     // -- Update SNES iteration count
522*ccaff030SJeremy L Thompson     PetscInt its;
523*ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
524*ccaff030SJeremy L Thompson     snesIts += its;
525*ccaff030SJeremy L Thompson   }
526*ccaff030SJeremy L Thompson 
527*ccaff030SJeremy L Thompson   // Timing
528*ccaff030SJeremy L Thompson   elapsedTime = MPI_Wtime() - startTime;
529*ccaff030SJeremy L Thompson 
530*ccaff030SJeremy L Thompson   // Performance logging
531*ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
532*ccaff030SJeremy L Thompson 
533*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
534*ccaff030SJeremy L Thompson   // Output summary
535*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
536*ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
537*ccaff030SJeremy L Thompson     // -- SNES
538*ccaff030SJeremy L Thompson     SNESType snesType;
539*ccaff030SJeremy L Thompson     SNESConvergedReason reason;
540*ccaff030SJeremy L Thompson     PetscReal rnorm;
541*ccaff030SJeremy L Thompson     ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr);
542*ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
543*ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
544*ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
545*ccaff030SJeremy L Thompson                        "  SNES:\n"
546*ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
547*ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
548*ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
549*ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
550*ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
551*ccaff030SJeremy L Thompson                        snesType, SNESConvergedReasons[reason],
552*ccaff030SJeremy L Thompson                        appCtx->numIncrements, snesIts, (double)rnorm);
553*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
554*ccaff030SJeremy L Thompson 
555*ccaff030SJeremy L Thompson     // -- KSP
556*ccaff030SJeremy L Thompson     KSP ksp;
557*ccaff030SJeremy L Thompson     KSPType kspType;
558*ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
559*ccaff030SJeremy L Thompson     ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr);
560*ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
561*ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
562*ccaff030SJeremy L Thompson                        "    KSP Type                           : %s\n",
563*ccaff030SJeremy L Thompson                        kspType); CHKERRQ(ierr);
564*ccaff030SJeremy L Thompson 
565*ccaff030SJeremy L Thompson     // -- PC
566*ccaff030SJeremy L Thompson     if (appCtx->multigridChoice != MULTIGRID_NONE) {
567*ccaff030SJeremy L Thompson       PC pc;
568*ccaff030SJeremy L Thompson       PCMGType pcmgType;
569*ccaff030SJeremy L Thompson       ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
570*ccaff030SJeremy L Thompson       ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr);
571*ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
572*ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
573*ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
574*ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
575*ccaff030SJeremy L Thompson                          PCMGTypes[pcmgType],
576*ccaff030SJeremy L Thompson                          PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr);
577*ccaff030SJeremy L Thompson 
578*ccaff030SJeremy L Thompson       // -- Coarse Solve
579*ccaff030SJeremy L Thompson       KSP kspCoarse;
580*ccaff030SJeremy L Thompson       PC pcCoarse;
581*ccaff030SJeremy L Thompson       PCType pcType;
582*ccaff030SJeremy L Thompson 
583*ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
584*ccaff030SJeremy L Thompson       ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr);
585*ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
586*ccaff030SJeremy L Thompson       ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr);
587*ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
588*ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
589*ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
590*ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
591*ccaff030SJeremy L Thompson                          kspType, pcType); CHKERRQ(ierr);
592*ccaff030SJeremy L Thompson     }
593*ccaff030SJeremy L Thompson   }
594*ccaff030SJeremy L Thompson 
595*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
596*ccaff030SJeremy L Thompson   // Compute solve time
597*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
598*ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
599*ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm);
600*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
601*ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm);
602*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
603*ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
604*ccaff030SJeremy L Thompson                        "  Performance:\n"
605*ccaff030SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n",
606*ccaff030SJeremy L Thompson                        maxTime, minTime); CHKERRQ(ierr);
607*ccaff030SJeremy L Thompson   }
608*ccaff030SJeremy L Thompson 
609*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
610*ccaff030SJeremy L Thompson   // Compute error
611*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
612*ccaff030SJeremy L Thompson   if (appCtx->forcingChoice == FORCE_MMS) {
613*ccaff030SJeremy L Thompson     CeedScalar l2Error = 1., l2Unorm = 1.;
614*ccaff030SJeremy L Thompson     const CeedScalar *truearray;
615*ccaff030SJeremy L Thompson     Vec errorVec, trueVec;
616*ccaff030SJeremy L Thompson 
617*ccaff030SJeremy L Thompson     // -- Work vectors
618*ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr);
619*ccaff030SJeremy L Thompson     ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr);
620*ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr);
621*ccaff030SJeremy L Thompson     ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr);
622*ccaff030SJeremy L Thompson 
623*ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
624*ccaff030SJeremy L Thompson     CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST,
625*ccaff030SJeremy L Thompson                            &truearray);
626*ccaff030SJeremy L Thompson     ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr);
627*ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec);
628*ccaff030SJeremy L Thompson     CHKERRQ(ierr);
629*ccaff030SJeremy L Thompson     ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr);
630*ccaff030SJeremy L Thompson     CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray);
631*ccaff030SJeremy L Thompson 
632*ccaff030SJeremy L Thompson     // -- Compute L2 error
633*ccaff030SJeremy L Thompson     ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr);
634*ccaff030SJeremy L Thompson     ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr);
635*ccaff030SJeremy L Thompson     ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr);
636*ccaff030SJeremy L Thompson     l2Error /= l2Unorm;
637*ccaff030SJeremy L Thompson 
638*ccaff030SJeremy L Thompson     // -- Output
639*ccaff030SJeremy L Thompson     if (!appCtx->testMode || l2Error > 0.05) {
640*ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
641*ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
642*ccaff030SJeremy L Thompson                          l2Error); CHKERRQ(ierr);
643*ccaff030SJeremy L Thompson     }
644*ccaff030SJeremy L Thompson 
645*ccaff030SJeremy L Thompson     // -- Cleanup
646*ccaff030SJeremy L Thompson     ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
647*ccaff030SJeremy L Thompson     ierr = VecDestroy(&trueVec); CHKERRQ(ierr);
648*ccaff030SJeremy L Thompson   }
649*ccaff030SJeremy L Thompson 
650*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
651*ccaff030SJeremy L Thompson   // Free objects
652*ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
653*ccaff030SJeremy L Thompson   // Data in arrays per level
654*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
655*ccaff030SJeremy L Thompson     // Vectors
656*ccaff030SJeremy L Thompson     ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr);
657*ccaff030SJeremy L Thompson     ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr);
658*ccaff030SJeremy L Thompson 
659*ccaff030SJeremy L Thompson     // Jacobian matrix and data
660*ccaff030SJeremy L Thompson     ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr);
661*ccaff030SJeremy L Thompson     ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr);
662*ccaff030SJeremy L Thompson     ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr);
663*ccaff030SJeremy L Thompson 
664*ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
665*ccaff030SJeremy L Thompson     if (level > 0) {
666*ccaff030SJeremy L Thompson       ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr);
667*ccaff030SJeremy L Thompson       ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr);
668*ccaff030SJeremy L Thompson       ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr);
669*ccaff030SJeremy L Thompson     }
670*ccaff030SJeremy L Thompson 
671*ccaff030SJeremy L Thompson     // DM
672*ccaff030SJeremy L Thompson     ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr);
673*ccaff030SJeremy L Thompson 
674*ccaff030SJeremy L Thompson     // libCEED objects
675*ccaff030SJeremy L Thompson     ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr);
676*ccaff030SJeremy L Thompson   }
677*ccaff030SJeremy L Thompson 
678*ccaff030SJeremy L Thompson   // Arrays
679*ccaff030SJeremy L Thompson   ierr = PetscFree(Ug); CHKERRQ(ierr);
680*ccaff030SJeremy L Thompson   ierr = PetscFree(Uloc); CHKERRQ(ierr);
681*ccaff030SJeremy L Thompson   ierr = PetscFree(Ugsz); CHKERRQ(ierr);
682*ccaff030SJeremy L Thompson   ierr = PetscFree(Ulsz); CHKERRQ(ierr);
683*ccaff030SJeremy L Thompson   ierr = PetscFree(Ulocsz); CHKERRQ(ierr);
684*ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCtx); CHKERRQ(ierr);
685*ccaff030SJeremy L Thompson   ierr = PetscFree(jacobMat); CHKERRQ(ierr);
686*ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr);
687*ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr);
688*ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr);
689*ccaff030SJeremy L Thompson   ierr = PetscFree(ceedData); CHKERRQ(ierr);
690*ccaff030SJeremy L Thompson 
691*ccaff030SJeremy L Thompson   // libCEED objects
692*ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfRestrict);
693*ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfProlong);
694*ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
695*ccaff030SJeremy L Thompson   CeedDestroy(&ceedFine);
696*ccaff030SJeremy L Thompson 
697*ccaff030SJeremy L Thompson   // PETSc objects
698*ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
699*ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
700*ccaff030SJeremy L Thompson   ierr = VecDestroy(&Rloc); CHKERRQ(ierr);
701*ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
702*ccaff030SJeremy L Thompson   ierr = VecDestroy(&Floc); CHKERRQ(ierr);
703*ccaff030SJeremy L Thompson   ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr);
704*ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
705*ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr);
706*ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
707*ccaff030SJeremy L Thompson   ierr = PetscFree(levelDMs); CHKERRQ(ierr);
708*ccaff030SJeremy L Thompson 
709*ccaff030SJeremy L Thompson   // Structs
710*ccaff030SJeremy L Thompson   ierr = PetscFree(resCtx); CHKERRQ(ierr);
711*ccaff030SJeremy L Thompson   ierr = PetscFree(formJacobCtx); CHKERRQ(ierr);
712*ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr);
713*ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx); CHKERRQ(ierr);
714*ccaff030SJeremy L Thompson   ierr = PetscFree(phys); CHKERRQ(ierr);
715*ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
716*ccaff030SJeremy L Thompson 
717*ccaff030SJeremy L Thompson   return PetscFinalize();
718*ccaff030SJeremy L Thompson }
719