xref: /libCEED/examples/solids/elasticity.c (revision e3e3df415f43d87f4bb6b086c337ffc22456aaae)
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
31ccaff030SJeremy 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
32ccaff030SJeremy 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
33ccaff030SJeremy L Thompson //
34ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes
35ccaff030SJeremy L Thompson //
367a3aec2eSjeremylt //TESTARGS -ceed {ceed_resource} -test -degree 2 -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
51ccaff030SJeremy L Thompson   Units          units;                  // Contains units scaling
52ccaff030SJeremy L Thompson   // PETSc objects
53ccaff030SJeremy L Thompson   PetscLogStage  stageDMSetup, stageLibceedSetup,
54ccaff030SJeremy L Thompson                  stageSnesSetup, stageSnesSolve;
55ccaff030SJeremy L Thompson   DM             dmOrig;                 // Distributed DM to clone
56ccaff030SJeremy L Thompson   DM             *levelDMs;
57ccaff030SJeremy L Thompson   Vec            U, *Ug, *Uloc;          // U: solution, R: residual, F: forcing
58ccaff030SJeremy L Thompson   Vec            R, Rloc, F, Floc;       // g: global, loc: local
59*e3e3df41Sjeremylt   SNES           snes, snesCoarse = NULL;
60ccaff030SJeremy L Thompson   Mat            *jacobMat, jacobMatCoarse, *prolongRestrMat;
61ccaff030SJeremy L Thompson   // PETSc data
62*e3e3df41Sjeremylt   UserMult       resCtx, jacobCoarseCtx = NULL, *jacobCtx;
63ccaff030SJeremy L Thompson   FormJacobCtx   formJacobCtx;
64ccaff030SJeremy L Thompson   UserMultProlongRestr *prolongRestrCtx;
65ccaff030SJeremy L Thompson   PCMGCycleType  pcmgCycleType = PC_MG_CYCLE_V;
66ccaff030SJeremy L Thompson   // libCEED objects
67ccaff030SJeremy L Thompson   Ceed           ceed, ceedFine = NULL;
68ccaff030SJeremy L Thompson   CeedData       *ceedData;
69f892e6d1Sjeremylt   CeedQFunction  qfRestrict = NULL, qfProlong = NULL;
70ccaff030SJeremy L Thompson   // Parameters
71ccaff030SJeremy L Thompson   PetscInt       ncompu = 3;             // 3 DoFs in 3D
72ccaff030SJeremy L Thompson   PetscInt       numLevels = 1, fineLevel = 0;
73ccaff030SJeremy L Thompson   PetscInt       *Ugsz, *Ulsz, *Ulocsz;  // sz: size
74ccaff030SJeremy L Thompson   PetscInt       snesIts = 0;
75ccaff030SJeremy L Thompson   // Timing
76ccaff030SJeremy L Thompson   double         startTime, elapsedTime, minTime, maxTime;
77ccaff030SJeremy L Thompson 
78ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
79ccaff030SJeremy L Thompson   if (ierr)
80ccaff030SJeremy L Thompson     return ierr;
81ccaff030SJeremy L Thompson 
82ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
83ccaff030SJeremy L Thompson   // Process command line options
84ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
85ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
86ccaff030SJeremy L Thompson 
87ccaff030SJeremy L Thompson   // -- Set mesh file, polynomial degree, problem type
88ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr);
89ccaff030SJeremy L Thompson   ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr);
90ccaff030SJeremy L Thompson   numLevels = appCtx->numLevels;
91ccaff030SJeremy L Thompson   fineLevel = numLevels - 1;
92ccaff030SJeremy L Thompson 
93ccaff030SJeremy L Thompson   // -- Set Poison's ratio, Young's Modulus
94ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr);
95ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
96ccaff030SJeremy L Thompson   ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr);
97ccaff030SJeremy L Thompson 
98ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
99ccaff030SJeremy L Thompson   // Setup DM
100ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
101ccaff030SJeremy L Thompson   // Performance logging
102ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup);
103ccaff030SJeremy L Thompson   CHKERRQ(ierr);
104ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr);
105ccaff030SJeremy L Thompson 
106ccaff030SJeremy L Thompson   // -- Create distributed DM from mesh file
107ccaff030SJeremy L Thompson   ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr);
108ccaff030SJeremy L Thompson 
109ccaff030SJeremy L Thompson   // -- Setup DM by polynomial degree
110ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr);
111*e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
112ccaff030SJeremy L Thompson     ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr);
113ccaff030SJeremy L Thompson     ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level],
114ccaff030SJeremy L Thompson                            ncompu); CHKERRQ(ierr);
115ccaff030SJeremy L Thompson   }
116ccaff030SJeremy L Thompson 
117ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
118ccaff030SJeremy L Thompson   // Setup solution and work vectors
119ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
120ccaff030SJeremy L Thompson   // Allocate arrays
121ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr);
122ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr);
123ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr);
124ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr);
125ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr);
126ccaff030SJeremy L Thompson 
127ccaff030SJeremy L Thompson   // -- Setup solution vectors for each level
128*e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
129ccaff030SJeremy L Thompson     // -- Create global unknown vector U
130ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr);
131ccaff030SJeremy L Thompson     ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr);
132ccaff030SJeremy L Thompson     // Note: Local size for matShell
133ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr);
134ccaff030SJeremy L Thompson 
135ccaff030SJeremy L Thompson     // -- Create local unknown vector Uloc
136ccaff030SJeremy L Thompson     ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr);
137ccaff030SJeremy L Thompson     // Note: local size for libCEED
138ccaff030SJeremy L Thompson     ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr);
139ccaff030SJeremy L Thompson   }
140ccaff030SJeremy L Thompson 
141ccaff030SJeremy L Thompson   // -- Create residual and forcing vectors
142ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr);
143ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr);
144ccaff030SJeremy L Thompson   ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr);
145ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr);
146ccaff030SJeremy L Thompson   ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr);
147ccaff030SJeremy L Thompson 
148ccaff030SJeremy L Thompson   // Performance logging
149ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
150ccaff030SJeremy L Thompson 
151ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
152ccaff030SJeremy L Thompson   // Set up libCEED
153ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
154ccaff030SJeremy L Thompson   // Performance logging
155ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup);
156ccaff030SJeremy L Thompson   CHKERRQ(ierr);
157ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr);
158ccaff030SJeremy L Thompson 
159ccaff030SJeremy L Thompson   // Initalize backend
160ccaff030SJeremy L Thompson   CeedInit(appCtx->ceedResource, &ceed);
161ccaff030SJeremy L Thompson   if (appCtx->degree > 4 && appCtx->ceedResourceFine[0])
162ccaff030SJeremy L Thompson     CeedInit(appCtx->ceedResourceFine, &ceedFine);
163ccaff030SJeremy L Thompson 
164ccaff030SJeremy L Thompson   // -- Create libCEED local forcing vector
165ccaff030SJeremy L Thompson   CeedVector forceCeed;
166ccaff030SJeremy L Thompson   CeedScalar *f;
167ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
168ccaff030SJeremy L Thompson     ierr = VecGetArray(Floc, &f); CHKERRQ(ierr);
169ccaff030SJeremy L Thompson     CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed);
170ccaff030SJeremy L Thompson     CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f);
171ccaff030SJeremy L Thompson   }
172ccaff030SJeremy L Thompson 
173ccaff030SJeremy L Thompson   // -- Restriction and prolongation QFunction
174ccaff030SJeremy L Thompson   if (appCtx->multigridChoice != MULTIGRID_NONE) {
175ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP,
176ccaff030SJeremy L Thompson                                 &qfRestrict);
177ccaff030SJeremy L Thompson     CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE,
178ccaff030SJeremy L Thompson                                 &qfProlong);
179ccaff030SJeremy L Thompson   }
180ccaff030SJeremy L Thompson 
181ccaff030SJeremy L Thompson   // -- Setup libCEED objects
182ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr);
183ccaff030SJeremy L Thompson   // ---- Setup residual evaluator and geometric information
184ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr);
185ccaff030SJeremy L Thompson   {
186ccaff030SJeremy L Thompson     bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine);
187ccaff030SJeremy L Thompson     Ceed levelCeed = highOrder ? ceedFine : ceed;
188ccaff030SJeremy L Thompson     ierr = SetupLibceedFineLevel(levelDMs[fineLevel], levelCeed, appCtx,
189ccaff030SJeremy L Thompson                                  phys, ceedData, fineLevel, ncompu,
190ccaff030SJeremy L Thompson                                  Ugsz[fineLevel], Ulocsz[fineLevel], forceCeed,
191ccaff030SJeremy L Thompson                                  qfRestrict, qfProlong);
192ccaff030SJeremy L Thompson     CHKERRQ(ierr);
193ccaff030SJeremy L Thompson   }
194ccaff030SJeremy L Thompson   // ---- Setup Jacobian evaluator and prolongation/restriction
195*e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
196ccaff030SJeremy L Thompson     if (level != fineLevel) {
197ccaff030SJeremy L Thompson       ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr);
198ccaff030SJeremy L Thompson     }
199ccaff030SJeremy L Thompson 
2001c8142b3Sjeremylt     // Note: use high order ceed, if specified and degree > 4
201ccaff030SJeremy L Thompson     bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine);
202ccaff030SJeremy L Thompson     Ceed levelCeed = highOrder ? ceedFine : ceed;
203ccaff030SJeremy L Thompson     ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys,
204ccaff030SJeremy L Thompson                              ceedData,  level, ncompu, Ugsz[level],
205ccaff030SJeremy L Thompson                              Ulocsz[level], forceCeed, qfRestrict,
206ccaff030SJeremy L Thompson                              qfProlong); CHKERRQ(ierr);
207ccaff030SJeremy L Thompson   }
208ccaff030SJeremy L Thompson 
209ccaff030SJeremy L Thompson   // Performance logging
210ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
211ccaff030SJeremy L Thompson 
212ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
213ccaff030SJeremy L Thompson   // Setup global forcing vector
214ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
215ccaff030SJeremy L Thompson   ierr = VecZeroEntries(F); CHKERRQ(ierr);
216ccaff030SJeremy L Thompson 
217ccaff030SJeremy L Thompson   if (appCtx->forcingChoice != FORCE_NONE) {
218ccaff030SJeremy L Thompson     ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr);
219ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F);
220ccaff030SJeremy L Thompson     CHKERRQ(ierr);
221ccaff030SJeremy L Thompson     CeedVectorDestroy(&forceCeed);
222ccaff030SJeremy L Thompson   }
223ccaff030SJeremy L Thompson 
224ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
225ccaff030SJeremy L Thompson   // Print problem summary
226ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
227ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
228ccaff030SJeremy L Thompson     const char *usedresource;
229ccaff030SJeremy L Thompson     CeedGetResource(ceed, &usedresource);
230ccaff030SJeremy L Thompson 
231ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
232ccaff030SJeremy L Thompson                        "\n-- Elastisticy Example - libCEED + PETSc --\n"
233ccaff030SJeremy L Thompson                        "  libCEED:\n"
234ccaff030SJeremy L Thompson                        "    libCEED Backend                    : %s\n",
235ccaff030SJeremy L Thompson                        usedresource); CHKERRQ(ierr);
236ccaff030SJeremy L Thompson 
237ccaff030SJeremy L Thompson     if (ceedFine) {
238ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
239ccaff030SJeremy L Thompson                          "    libCEED Backend - high order       : %s\n",
240ccaff030SJeremy L Thompson                          appCtx->ceedResourceFine); CHKERRQ(ierr);
241ccaff030SJeremy L Thompson     }
242ccaff030SJeremy L Thompson 
243ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
244ccaff030SJeremy L Thompson                        "  Problem:\n"
245ccaff030SJeremy L Thompson                        "    Problem Name                       : %s\n"
246ccaff030SJeremy L Thompson                        "    Forcing Function                   : %s\n"
247ccaff030SJeremy L Thompson                        "  Mesh:\n"
248ccaff030SJeremy L Thompson                        "    File                               : %s\n"
249ccaff030SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
250ccaff030SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
251ccaff030SJeremy L Thompson                        "    Global nodes                       : %D\n"
252ccaff030SJeremy L Thompson                        "    Owned nodes                        : %D\n"
253ccaff030SJeremy L Thompson                        "    DoF per node                       : %D\n"
254ccaff030SJeremy L Thompson                        "  Multigrid:\n"
255ccaff030SJeremy L Thompson                        "    Type                               : %s\n"
256ccaff030SJeremy L Thompson                        "    Number of Levels                   : %d\n",
257ccaff030SJeremy L Thompson                        problemTypesForDisp[appCtx->problemChoice],
258ccaff030SJeremy L Thompson                        forcingTypesForDisp[appCtx->forcingChoice],
259ccaff030SJeremy L Thompson                        appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh",
260ccaff030SJeremy L Thompson                        appCtx->degree + 1, appCtx->degree + 1,
261ccaff030SJeremy L Thompson                        Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu,
262f892e6d1Sjeremylt                        (appCtx->degree == 1 &&
263f892e6d1Sjeremylt                         appCtx->multigridChoice != MULTIGRID_NONE) ?
264f892e6d1Sjeremylt                        "Algebraic multigrid" :
265ccaff030SJeremy L Thompson                        multigridTypesForDisp[appCtx->multigridChoice],
266f892e6d1Sjeremylt                        (appCtx->degree == 1 ||
267f892e6d1Sjeremylt                         appCtx->multigridChoice == MULTIGRID_NONE) ?
268f892e6d1Sjeremylt                        0 : numLevels); CHKERRQ(ierr);
269ccaff030SJeremy L Thompson 
270ccaff030SJeremy L Thompson     if (appCtx->multigridChoice != MULTIGRID_NONE) {
271*e3e3df41Sjeremylt       for (PetscInt i = 0; i < 2; i++) {
272ccaff030SJeremy L Thompson         CeedInt level = i ? fineLevel : 0;
273ccaff030SJeremy L Thompson         ierr = PetscPrintf(comm,
274ccaff030SJeremy L Thompson                            "    Level %D (%s):\n"
275ccaff030SJeremy L Thompson                            "      Number of 1D Basis Nodes (p)     : %d\n"
276ccaff030SJeremy L Thompson                            "      Global Nodes                     : %D\n"
277ccaff030SJeremy L Thompson                            "      Owned Nodes                      : %D\n",
278ccaff030SJeremy L Thompson                            level, i ? "fine" : "coarse",
279ccaff030SJeremy L Thompson                            appCtx->levelDegrees[level] + 1,
280ccaff030SJeremy L Thompson                            Ugsz[level]/ncompu, Ulsz[level]/ncompu);
281ccaff030SJeremy L Thompson         CHKERRQ(ierr);
282ccaff030SJeremy L Thompson       }
283ccaff030SJeremy L Thompson     }
284ccaff030SJeremy L Thompson   }
285ccaff030SJeremy L Thompson 
286ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
287ccaff030SJeremy L Thompson   // Setup SNES
288ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
289ccaff030SJeremy L Thompson   // Performance logging
290ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup);
291ccaff030SJeremy L Thompson   CHKERRQ(ierr);
292ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr);
293ccaff030SJeremy L Thompson 
294ccaff030SJeremy L Thompson   // Create SNES
295ccaff030SJeremy L Thompson   ierr = SNESCreate(comm, &snes); CHKERRQ(ierr);
296ccaff030SJeremy L Thompson   ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr);
297ccaff030SJeremy L Thompson 
298ccaff030SJeremy L Thompson   // -- Jacobian evaluators
299ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr);
300ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr);
301*e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
302ccaff030SJeremy L Thompson     // -- Jacobian context for level
303ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr);
304ccaff030SJeremy L Thompson     ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level],
305ccaff030SJeremy L Thompson                             Uloc[level], ceedData[level], ceed,
306ccaff030SJeremy L Thompson                             jacobCtx[level]); CHKERRQ(ierr);
307ccaff030SJeremy L Thompson 
308ccaff030SJeremy L Thompson     // -- Form Action of Jacobian on delta_u
309ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level],
310ccaff030SJeremy L Thompson                           Ugsz[level], jacobCtx[level], &jacobMat[level]);
311ccaff030SJeremy L Thompson     CHKERRQ(ierr);
312ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT,
313ccaff030SJeremy L Thompson                                 (void (*)(void))ApplyJacobian_Ceed);
314ccaff030SJeremy L Thompson     CHKERRQ(ierr);
315ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL,
316ccaff030SJeremy L Thompson                                 (void(*)(void))GetDiag_Ceed);
317ccaff030SJeremy L Thompson 
318ccaff030SJeremy L Thompson   }
319*e3e3df41Sjeremylt   // Note: FormJacobian updates Jacobian matrices on each level
320ccaff030SJeremy L Thompson   //   and assembles the Jpre matrix, if needed
321ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr);
322ccaff030SJeremy L Thompson   formJacobCtx->jacobCtx = jacobCtx;
323ccaff030SJeremy L Thompson   formJacobCtx->numLevels = numLevels;
324ccaff030SJeremy L Thompson   formJacobCtx->jacobMat = jacobMat;
325ccaff030SJeremy L Thompson 
326ccaff030SJeremy L Thompson   // -- Residual evaluation function
327ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr);
328ccaff030SJeremy L Thompson   ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel],
329ccaff030SJeremy L Thompson                      sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr);
330ccaff030SJeremy L Thompson   resCtx->op = ceedData[fineLevel]->opApply;
331ccaff030SJeremy L Thompson   ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr);
332ccaff030SJeremy L Thompson 
333ccaff030SJeremy L Thompson   // -- Prolongation/Restriction evaluation
334ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr);
335ccaff030SJeremy L Thompson   ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr);
336*e3e3df41Sjeremylt   for (PetscInt level = 1; level < numLevels; level++) {
337ccaff030SJeremy L Thompson     // ---- Prolongation/restriction context for level
338ccaff030SJeremy L Thompson     ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr);
339ccaff030SJeremy L Thompson     ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level],
340ccaff030SJeremy L Thompson                                    Ug[level], Uloc[level-1], Uloc[level],
341ccaff030SJeremy L Thompson                                    ceedData[level-1], ceedData[level], ceed,
342ccaff030SJeremy L Thompson                                    prolongRestrCtx[level]); CHKERRQ(ierr);
343ccaff030SJeremy L Thompson 
344ccaff030SJeremy L Thompson     // ---- Form Action of Jacobian on delta_u
345ccaff030SJeremy L Thompson     ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level],
346ccaff030SJeremy L Thompson                           Ugsz[level-1], prolongRestrCtx[level],
347ccaff030SJeremy L Thompson                           &prolongRestrMat[level]); CHKERRQ(ierr);
348ccaff030SJeremy L Thompson     // Note: In PCMG, restriction is the transpose of prolongation
349ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT,
350ccaff030SJeremy L Thompson                                 (void (*)(void))Prolong_Ceed);
351ccaff030SJeremy L Thompson     ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE,
352ccaff030SJeremy L Thompson                                 (void (*)(void))Restrict_Ceed);
353ccaff030SJeremy L Thompson     CHKERRQ(ierr);
354ccaff030SJeremy L Thompson   }
355ccaff030SJeremy L Thompson 
356ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
357ccaff030SJeremy L Thompson   // Setup dummy SNES for AMG coarse solve
358ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
359*e3e3df41Sjeremylt   if (appCtx->multigridChoice != MULTIGRID_NONE) {
360*e3e3df41Sjeremylt     // -- Jacobian Matrix
361*e3e3df41Sjeremylt     ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr);
362*e3e3df41Sjeremylt     ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr);
363*e3e3df41Sjeremylt 
364*e3e3df41Sjeremylt     if (appCtx->degree > 1) {
365ccaff030SJeremy L Thompson       ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr);
366ccaff030SJeremy L Thompson       ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr);
367ccaff030SJeremy L Thompson       ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr);
368ccaff030SJeremy L Thompson 
369*e3e3df41Sjeremylt       // -- Jacobian function
370ccaff030SJeremy L Thompson       ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL,
371ccaff030SJeremy L Thompson                              NULL); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson 
373ccaff030SJeremy L Thompson       // -- Residual evaluation function
374ccaff030SJeremy L Thompson       ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr);
375ccaff030SJeremy L Thompson       ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0]));
376ccaff030SJeremy L Thompson       CHKERRQ(ierr);
377ccaff030SJeremy L Thompson       ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed,
378ccaff030SJeremy L Thompson                              jacobCoarseCtx); CHKERRQ(ierr);
379ccaff030SJeremy L Thompson 
380ccaff030SJeremy L Thompson       // -- Update formJacobCtx
381ccaff030SJeremy L Thompson       formJacobCtx->Ucoarse = Ug[0];
382ccaff030SJeremy L Thompson       formJacobCtx->snesCoarse = snesCoarse;
383ccaff030SJeremy L Thompson       formJacobCtx->jacobMatCoarse = jacobMatCoarse;
384*e3e3df41Sjeremylt     }
385*e3e3df41Sjeremylt   }
386*e3e3df41Sjeremylt 
387*e3e3df41Sjeremylt   // Set Jacobian function
388*e3e3df41Sjeremylt   if (appCtx->degree > 1) {
389*e3e3df41Sjeremylt       ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel],
390*e3e3df41Sjeremylt                              FormJacobian, formJacobCtx); CHKERRQ(ierr);
391*e3e3df41Sjeremylt   } else {
392*e3e3df41Sjeremylt     ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse,
393*e3e3df41Sjeremylt                            SNESComputeJacobianDefaultColor, NULL);
394*e3e3df41Sjeremylt     CHKERRQ(ierr);
395*e3e3df41Sjeremylt   }
396ccaff030SJeremy L Thompson 
397ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
398ccaff030SJeremy L Thompson   // Setup KSP
399ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
400ccaff030SJeremy L Thompson   {
401ccaff030SJeremy L Thompson     PC pc;
402ccaff030SJeremy L Thompson     KSP ksp;
403ccaff030SJeremy L Thompson 
404ccaff030SJeremy L Thompson     // -- KSP
405ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
406ccaff030SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
408ccaff030SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
409ccaff030SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson     ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr);
411ccaff030SJeremy L Thompson 
412ccaff030SJeremy L Thompson     // -- Preconditioning
413ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
414ccaff030SJeremy L Thompson     ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr);
415ccaff030SJeremy L Thompson     ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr);
416ccaff030SJeremy L Thompson 
417ccaff030SJeremy L Thompson     if (appCtx->multigridChoice == MULTIGRID_NONE) {
418ccaff030SJeremy L Thompson       // ---- No Multigrid
419ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
420ccaff030SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
421f892e6d1Sjeremylt     } else if (appCtx->degree == 1) {
422f892e6d1Sjeremylt       // ---- AMG for degree 1
423f892e6d1Sjeremylt       ierr = KSPSetType(ksp, KSPPREONLY); CHKERRQ(ierr);
424f892e6d1Sjeremylt       ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr);
425ccaff030SJeremy L Thompson     } else {
426ccaff030SJeremy L Thompson       // ---- PCMG
427ccaff030SJeremy L Thompson       ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
428ccaff030SJeremy L Thompson 
429ccaff030SJeremy L Thompson       // ------ PCMG levels
430ccaff030SJeremy L Thompson       ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr);
431*e3e3df41Sjeremylt       for (PetscInt level = 0; level < numLevels; level++) {
432ccaff030SJeremy L Thompson         // -------- Smoother
433ccaff030SJeremy L Thompson         KSP kspSmoother, kspEst;
434ccaff030SJeremy L Thompson         PC pcSmoother;
435ccaff030SJeremy L Thompson 
436ccaff030SJeremy L Thompson         // ---------- Smoother KSP
437ccaff030SJeremy L Thompson         ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr);
438ccaff030SJeremy L Thompson         ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr);
439ccaff030SJeremy L Thompson         ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr);
440ccaff030SJeremy L Thompson 
441ccaff030SJeremy L Thompson         // ---------- Chebyshev options
442ccaff030SJeremy L Thompson         ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr);
443ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1);
444ccaff030SJeremy L Thompson         CHKERRQ(ierr);
445ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr);
446ccaff030SJeremy L Thompson         ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr);
447ccaff030SJeremy L Thompson         ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE);
448ccaff030SJeremy L Thompson         CHKERRQ(ierr);
449ccaff030SJeremy L Thompson         ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]);
450ccaff030SJeremy L Thompson         CHKERRQ(ierr);
451ccaff030SJeremy L Thompson 
452ccaff030SJeremy L Thompson         // ---------- Smoother preconditioner
453ccaff030SJeremy L Thompson         ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr);
454ccaff030SJeremy L Thompson         ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr);
455ccaff030SJeremy L Thompson         ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
456ccaff030SJeremy L Thompson 
457ccaff030SJeremy L Thompson         // -------- Work vector
458ccaff030SJeremy L Thompson         if (level != fineLevel) {
459ccaff030SJeremy L Thompson           ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr);
460ccaff030SJeremy L Thompson         }
461ccaff030SJeremy L Thompson 
462ccaff030SJeremy L Thompson         // -------- Level prolongation/restriction operator
463ccaff030SJeremy L Thompson         if (level > 0) {
464ccaff030SJeremy L Thompson           ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]);
465ccaff030SJeremy L Thompson           CHKERRQ(ierr);
466ccaff030SJeremy L Thompson           ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]);
467ccaff030SJeremy L Thompson           CHKERRQ(ierr);
468ccaff030SJeremy L Thompson         }
469ccaff030SJeremy L Thompson       }
470ccaff030SJeremy L Thompson 
471ccaff030SJeremy L Thompson       // ------ PCMG coarse solve
472ccaff030SJeremy L Thompson       KSP kspCoarse;
473ccaff030SJeremy L Thompson       PC pcCoarse;
474ccaff030SJeremy L Thompson 
475ccaff030SJeremy L Thompson       // -------- Coarse KSP
476ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson       ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson       ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse);
479ccaff030SJeremy L Thompson       CHKERRQ(ierr);
480ccaff030SJeremy L Thompson       ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson 
482ccaff030SJeremy L Thompson       // -------- Coarse preconditioner
483ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
484ccaff030SJeremy L Thompson       ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr);
485ccaff030SJeremy L Thompson       ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr);
486ccaff030SJeremy L Thompson 
487ccaff030SJeremy L Thompson       ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr);
488ccaff030SJeremy L Thompson       ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr);
489ccaff030SJeremy L Thompson 
490ccaff030SJeremy L Thompson       // ------ PCMG options
491ccaff030SJeremy L Thompson       ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson       ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson       ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr);
494ccaff030SJeremy L Thompson     }
495ccaff030SJeremy L Thompson     ierr = KSPSetFromOptions(ksp);
496ccaff030SJeremy L Thompson     ierr = PCSetFromOptions(pc);
497ccaff030SJeremy L Thompson   }
498038d0cd7Sjeremylt   {
499038d0cd7Sjeremylt     // Default to critical-point (CP) line search (related to Wolfe's curvature condition)
500481a4840SJed Brown     SNESLineSearch linesearch;
501481a4840SJed Brown 
502481a4840SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr);
503481a4840SJed Brown     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr);
504481a4840SJed Brown   }
505481a4840SJed Brown 
506ccaff030SJeremy L Thompson   ierr = SNESSetFromOptions(snes); CHKERRQ(ierr);
507ccaff030SJeremy L Thompson 
508ccaff030SJeremy L Thompson   // Performance logging
509ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
510ccaff030SJeremy L Thompson 
511ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
512ccaff030SJeremy L Thompson   // Set initial guess
513ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
514ccaff030SJeremy L Thompson   ierr = VecSet(U, 0.0); CHKERRQ(ierr);
515ccaff030SJeremy L Thompson 
516ccaff030SJeremy L Thompson   // View solution
517ccaff030SJeremy L Thompson   if (appCtx->viewSoln) {
518ccaff030SJeremy L Thompson     ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson   }
520ccaff030SJeremy L Thompson 
521ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
522ccaff030SJeremy L Thompson   // Solve SNES
523ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
5245f24f028Sjeremylt   PetscBool snesMonitor = PETSC_FALSE;
5255f24f028Sjeremylt   ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor);
5265f24f028Sjeremylt   CHKERRQ(ierr);
5275f24f028Sjeremylt 
528ccaff030SJeremy L Thompson   // Performance logging
529ccaff030SJeremy L Thompson   ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve);
530ccaff030SJeremy L Thompson   CHKERRQ(ierr);
531ccaff030SJeremy L Thompson   ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr);
532ccaff030SJeremy L Thompson 
533ccaff030SJeremy L Thompson   // Timing
534ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   startTime = MPI_Wtime();
536ccaff030SJeremy L Thompson 
537ccaff030SJeremy L Thompson   // Solve for each load increment
5385f24f028Sjeremylt   PetscInt increment;
5395f24f028Sjeremylt   for (increment = 1; increment <= appCtx->numIncrements; increment++) {
5405f24f028Sjeremylt     // -- Log increment count
5415f24f028Sjeremylt     if (snesMonitor) {
5425f24f028Sjeremylt       ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1);
5435f24f028Sjeremylt       CHKERRQ(ierr);
5445f24f028Sjeremylt     }
5455f24f028Sjeremylt 
546ccaff030SJeremy L Thompson     // -- Scale the problem
547ccaff030SJeremy L Thompson     PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements,
548ccaff030SJeremy L Thompson                 scalingFactor = loadIncrement /
549ccaff030SJeremy L Thompson                                 (increment == 1 ? 1 : resCtx->loadIncrement);
550ccaff030SJeremy L Thompson     resCtx->loadIncrement = loadIncrement;
551ccaff030SJeremy L Thompson     if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) {
552ccaff030SJeremy L Thompson       ierr = VecScale(F, scalingFactor); CHKERRQ(ierr);
553ccaff030SJeremy L Thompson     }
554ccaff030SJeremy L Thompson 
555ccaff030SJeremy L Thompson     // -- Solve
556ccaff030SJeremy L Thompson     ierr = SNESSolve(snes, F, U); CHKERRQ(ierr);
557ccaff030SJeremy L Thompson 
558ccaff030SJeremy L Thompson     // -- View solution
559ccaff030SJeremy L Thompson     if (appCtx->viewSoln) {
560ccaff030SJeremy L Thompson       ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson     }
562ccaff030SJeremy L Thompson 
563ccaff030SJeremy L Thompson     // -- Update SNES iteration count
564ccaff030SJeremy L Thompson     PetscInt its;
565ccaff030SJeremy L Thompson     ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson     snesIts += its;
5673951731eSjeremylt 
5683951731eSjeremylt     // -- Check for divergence
5693951731eSjeremylt     SNESConvergedReason reason;
5703951731eSjeremylt     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
5713951731eSjeremylt     if (reason < 0)
5723951731eSjeremylt       break;
573ccaff030SJeremy L Thompson   }
574ccaff030SJeremy L Thompson 
575ccaff030SJeremy L Thompson   // Timing
576ccaff030SJeremy L Thompson   elapsedTime = MPI_Wtime() - startTime;
577ccaff030SJeremy L Thompson 
578ccaff030SJeremy L Thompson   // Performance logging
579ccaff030SJeremy L Thompson   ierr = PetscLogStagePop();
580ccaff030SJeremy L Thompson 
581ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
582ccaff030SJeremy L Thompson   // Output summary
583ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
584ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
585ccaff030SJeremy L Thompson     // -- SNES
586ccaff030SJeremy L Thompson     SNESType snesType;
587ccaff030SJeremy L Thompson     SNESConvergedReason reason;
588ccaff030SJeremy L Thompson     PetscReal rnorm;
589ccaff030SJeremy L Thompson     ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson     ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr);
591ccaff030SJeremy L Thompson     ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
593ccaff030SJeremy L Thompson                        "  SNES:\n"
594ccaff030SJeremy L Thompson                        "    SNES Type                          : %s\n"
595ccaff030SJeremy L Thompson                        "    SNES Convergence                   : %s\n"
596ccaff030SJeremy L Thompson                        "    Number of Load Increments          : %d\n"
5975f24f028Sjeremylt                        "    Completed Load Increments          : %d\n"
598ccaff030SJeremy L Thompson                        "    Total SNES Iterations              : %D\n"
599ccaff030SJeremy L Thompson                        "    Final rnorm                        : %e\n",
600ccaff030SJeremy L Thompson                        snesType, SNESConvergedReasons[reason],
6015f24f028Sjeremylt                        appCtx->numIncrements, increment - 1,
6025f24f028Sjeremylt                        snesIts, (double)rnorm); CHKERRQ(ierr);
603ccaff030SJeremy L Thompson 
604ccaff030SJeremy L Thompson     // -- KSP
605ccaff030SJeremy L Thompson     KSP ksp;
606ccaff030SJeremy L Thompson     KSPType kspType;
607ccaff030SJeremy L Thompson     ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson     ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr);
609ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
610ccaff030SJeremy L Thompson                        "  Linear Solver:\n"
611ccaff030SJeremy L Thompson                        "    KSP Type                           : %s\n",
612ccaff030SJeremy L Thompson                        kspType); CHKERRQ(ierr);
613ccaff030SJeremy L Thompson 
614ccaff030SJeremy L Thompson     // -- PC
615ccaff030SJeremy L Thompson     PC pc;
616*e3e3df41Sjeremylt     PCType pcType;
617ccaff030SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
618*e3e3df41Sjeremylt     ierr = PCGetType(pc, &pcType); CHKERRQ(ierr);
619*e3e3df41Sjeremylt     ierr = PetscPrintf(comm,
620*e3e3df41Sjeremylt                        "    PC Type                            : %s\n",
621*e3e3df41Sjeremylt                        pcType); CHKERRQ(ierr);
622*e3e3df41Sjeremylt 
623*e3e3df41Sjeremylt     if (appCtx->multigridChoice != MULTIGRID_NONE && appCtx->degree > 1) {
624*e3e3df41Sjeremylt       PCMGType pcmgType;
625ccaff030SJeremy L Thompson       ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr);
626ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
627ccaff030SJeremy L Thompson                          "  P-Multigrid:\n"
628ccaff030SJeremy L Thompson                          "    PCMG Type                          : %s\n"
629ccaff030SJeremy L Thompson                          "    PCMG Cycle Type                    : %s\n",
630ccaff030SJeremy L Thompson                          PCMGTypes[pcmgType],
631ccaff030SJeremy L Thompson                          PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr);
632ccaff030SJeremy L Thompson 
633ccaff030SJeremy L Thompson       // -- Coarse Solve
634ccaff030SJeremy L Thompson       KSP kspCoarse;
635ccaff030SJeremy L Thompson       PC pcCoarse;
636ccaff030SJeremy L Thompson       PCType pcType;
637ccaff030SJeremy L Thompson 
638ccaff030SJeremy L Thompson       ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr);
639ccaff030SJeremy L Thompson       ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr);
640ccaff030SJeremy L Thompson       ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr);
641ccaff030SJeremy L Thompson       ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
643ccaff030SJeremy L Thompson                          "    Coarse Solve:\n"
644ccaff030SJeremy L Thompson                          "      KSP Type                         : %s\n"
645ccaff030SJeremy L Thompson                          "      PC Type                          : %s\n",
646ccaff030SJeremy L Thompson                          kspType, pcType); CHKERRQ(ierr);
647ccaff030SJeremy L Thompson     }
648ccaff030SJeremy L Thompson   }
649ccaff030SJeremy L Thompson 
650ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
651ccaff030SJeremy L Thompson   // Compute solve time
652ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
653ccaff030SJeremy L Thompson   if (!appCtx->testMode) {
654ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm);
655ccaff030SJeremy L Thompson     CHKERRQ(ierr);
656ccaff030SJeremy L Thompson     ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm);
657ccaff030SJeremy L Thompson     CHKERRQ(ierr);
658ccaff030SJeremy L Thompson     ierr = PetscPrintf(comm,
659ccaff030SJeremy L Thompson                        "  Performance:\n"
660ccaff030SJeremy L Thompson                        "    SNES Solve Time                    : %g (%g) sec\n",
661ccaff030SJeremy L Thompson                        maxTime, minTime); CHKERRQ(ierr);
662ccaff030SJeremy L Thompson   }
663ccaff030SJeremy L Thompson 
664ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
665ccaff030SJeremy L Thompson   // Compute error
666ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
667ccaff030SJeremy L Thompson   if (appCtx->forcingChoice == FORCE_MMS) {
668ccaff030SJeremy L Thompson     CeedScalar l2Error = 1., l2Unorm = 1.;
669ccaff030SJeremy L Thompson     const CeedScalar *truearray;
670ccaff030SJeremy L Thompson     Vec errorVec, trueVec;
671ccaff030SJeremy L Thompson 
672ccaff030SJeremy L Thompson     // -- Work vectors
673ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr);
674ccaff030SJeremy L Thompson     ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr);
675ccaff030SJeremy L Thompson     ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr);
676ccaff030SJeremy L Thompson     ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr);
677ccaff030SJeremy L Thompson 
678ccaff030SJeremy L Thompson     // -- Assemble global true solution vector
679ccaff030SJeremy L Thompson     CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST,
680ccaff030SJeremy L Thompson                            &truearray);
681ccaff030SJeremy L Thompson     ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr);
682ccaff030SJeremy L Thompson     ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec);
683ccaff030SJeremy L Thompson     CHKERRQ(ierr);
684ccaff030SJeremy L Thompson     ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr);
685ccaff030SJeremy L Thompson     CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray);
686ccaff030SJeremy L Thompson 
687ccaff030SJeremy L Thompson     // -- Compute L2 error
688ccaff030SJeremy L Thompson     ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr);
689ccaff030SJeremy L Thompson     ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr);
690ccaff030SJeremy L Thompson     ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr);
691ccaff030SJeremy L Thompson     l2Error /= l2Unorm;
692ccaff030SJeremy L Thompson 
693ccaff030SJeremy L Thompson     // -- Output
694ccaff030SJeremy L Thompson     if (!appCtx->testMode || l2Error > 0.05) {
695ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm,
696ccaff030SJeremy L Thompson                          "    L2 Error                           : %e\n",
697ccaff030SJeremy L Thompson                          l2Error); CHKERRQ(ierr);
698ccaff030SJeremy L Thompson     }
699ccaff030SJeremy L Thompson 
700ccaff030SJeremy L Thompson     // -- Cleanup
701ccaff030SJeremy L Thompson     ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
702ccaff030SJeremy L Thompson     ierr = VecDestroy(&trueVec); CHKERRQ(ierr);
703ccaff030SJeremy L Thompson   }
704ccaff030SJeremy L Thompson 
705ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
706ccaff030SJeremy L Thompson   // Free objects
707ccaff030SJeremy L Thompson   // ---------------------------------------------------------------------------
708ccaff030SJeremy L Thompson   // Data in arrays per level
709*e3e3df41Sjeremylt   for (PetscInt level = 0; level < numLevels; level++) {
710ccaff030SJeremy L Thompson     // Vectors
711ccaff030SJeremy L Thompson     ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr);
712ccaff030SJeremy L Thompson     ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr);
713ccaff030SJeremy L Thompson 
714ccaff030SJeremy L Thompson     // Jacobian matrix and data
715ccaff030SJeremy L Thompson     ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr);
716ccaff030SJeremy L Thompson     ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr);
717ccaff030SJeremy L Thompson     ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr);
718ccaff030SJeremy L Thompson 
719ccaff030SJeremy L Thompson     // Prolongation/Restriction matrix and data
720ccaff030SJeremy L Thompson     if (level > 0) {
721ccaff030SJeremy L Thompson       ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr);
722ccaff030SJeremy L Thompson       ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr);
723ccaff030SJeremy L Thompson       ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr);
724ccaff030SJeremy L Thompson     }
725ccaff030SJeremy L Thompson 
726ccaff030SJeremy L Thompson     // DM
727ccaff030SJeremy L Thompson     ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr);
728ccaff030SJeremy L Thompson 
729ccaff030SJeremy L Thompson     // libCEED objects
730ccaff030SJeremy L Thompson     ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr);
731ccaff030SJeremy L Thompson   }
732ccaff030SJeremy L Thompson 
733ccaff030SJeremy L Thompson   // Arrays
734ccaff030SJeremy L Thompson   ierr = PetscFree(Ug); CHKERRQ(ierr);
735ccaff030SJeremy L Thompson   ierr = PetscFree(Uloc); CHKERRQ(ierr);
736ccaff030SJeremy L Thompson   ierr = PetscFree(Ugsz); CHKERRQ(ierr);
737ccaff030SJeremy L Thompson   ierr = PetscFree(Ulsz); CHKERRQ(ierr);
738ccaff030SJeremy L Thompson   ierr = PetscFree(Ulocsz); CHKERRQ(ierr);
739ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCtx); CHKERRQ(ierr);
740ccaff030SJeremy L Thompson   ierr = PetscFree(jacobMat); CHKERRQ(ierr);
741ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr);
742ccaff030SJeremy L Thompson   ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr);
743ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr);
744ccaff030SJeremy L Thompson   ierr = PetscFree(ceedData); CHKERRQ(ierr);
745ccaff030SJeremy L Thompson 
746ccaff030SJeremy L Thompson   // libCEED objects
747ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfRestrict);
748ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qfProlong);
749ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
750ccaff030SJeremy L Thompson   CeedDestroy(&ceedFine);
751ccaff030SJeremy L Thompson 
752ccaff030SJeremy L Thompson   // PETSc objects
753ccaff030SJeremy L Thompson   ierr = VecDestroy(&U); CHKERRQ(ierr);
754ccaff030SJeremy L Thompson   ierr = VecDestroy(&R); CHKERRQ(ierr);
755ccaff030SJeremy L Thompson   ierr = VecDestroy(&Rloc); CHKERRQ(ierr);
756ccaff030SJeremy L Thompson   ierr = VecDestroy(&F); CHKERRQ(ierr);
757ccaff030SJeremy L Thompson   ierr = VecDestroy(&Floc); CHKERRQ(ierr);
758ccaff030SJeremy L Thompson   ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr);
759ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snes); CHKERRQ(ierr);
760ccaff030SJeremy L Thompson   ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr);
761ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmOrig); CHKERRQ(ierr);
762ccaff030SJeremy L Thompson   ierr = PetscFree(levelDMs); CHKERRQ(ierr);
763ccaff030SJeremy L Thompson 
764ccaff030SJeremy L Thompson   // Structs
765ccaff030SJeremy L Thompson   ierr = PetscFree(resCtx); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson   ierr = PetscFree(formJacobCtx); CHKERRQ(ierr);
767ccaff030SJeremy L Thompson   ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr);
768ccaff030SJeremy L Thompson   ierr = PetscFree(appCtx); CHKERRQ(ierr);
769ccaff030SJeremy L Thompson   ierr = PetscFree(phys); CHKERRQ(ierr);
770ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
771ccaff030SJeremy L Thompson 
772ccaff030SJeremy L Thompson   return PetscFinalize();
773ccaff030SJeremy L Thompson }
774