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