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 51*f7b4142eSJeremy L Thompson Physics physSmoother = NULL; // Separate context if nuSmoother set 52ccaff030SJeremy L Thompson Units units; // Contains units scaling 53ccaff030SJeremy L Thompson // PETSc objects 54ccaff030SJeremy L Thompson PetscLogStage stageDMSetup, stageLibceedSetup, 55ccaff030SJeremy L Thompson stageSnesSetup, stageSnesSolve; 56ccaff030SJeremy L Thompson DM dmOrig; // Distributed DM to clone 57a3c02c40SJeremy L Thompson DM dmEnergy, dmDiagnostic; // DMs for postprocessing 58ccaff030SJeremy L Thompson DM *levelDMs; 59ccaff030SJeremy L Thompson Vec U, *Ug, *Uloc; // U: solution, R: residual, F: forcing 60ccaff030SJeremy L Thompson Vec R, Rloc, F, Floc; // g: global, loc: local 61e3e3df41Sjeremylt SNES snes, snesCoarse = NULL; 62ccaff030SJeremy L Thompson Mat *jacobMat, jacobMatCoarse, *prolongRestrMat; 63ccaff030SJeremy L Thompson // PETSc data 64e3e3df41Sjeremylt UserMult resCtx, jacobCoarseCtx = NULL, *jacobCtx; 65ccaff030SJeremy L Thompson FormJacobCtx formJacobCtx; 66ccaff030SJeremy L Thompson UserMultProlongRestr *prolongRestrCtx; 67ccaff030SJeremy L Thompson PCMGCycleType pcmgCycleType = PC_MG_CYCLE_V; 68ccaff030SJeremy L Thompson // libCEED objects 69ccaff030SJeremy L Thompson Ceed ceed, ceedFine = NULL; 70ccaff030SJeremy L Thompson CeedData *ceedData; 71f892e6d1Sjeremylt CeedQFunction qfRestrict = NULL, qfProlong = NULL; 72ccaff030SJeremy L Thompson // Parameters 73ccaff030SJeremy L Thompson PetscInt ncompu = 3; // 3 DoFs in 3D 745c25879aSJeremy L Thompson PetscInt ncompe = 1, ncompd = 2; // 1 energy output, 2 diagnostic 75ccaff030SJeremy L Thompson PetscInt numLevels = 1, fineLevel = 0; 76ccaff030SJeremy L Thompson PetscInt *Ugsz, *Ulsz, *Ulocsz; // sz: size 77ccaff030SJeremy L Thompson PetscInt snesIts = 0; 78ccaff030SJeremy L Thompson // Timing 79ccaff030SJeremy L Thompson double startTime, elapsedTime, minTime, maxTime; 80ccaff030SJeremy L Thompson 81ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 82ccaff030SJeremy L Thompson if (ierr) 83ccaff030SJeremy L Thompson return ierr; 84ccaff030SJeremy L Thompson 85ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 86ccaff030SJeremy L Thompson // Process command line options 87ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 88ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 89ccaff030SJeremy L Thompson 90ccaff030SJeremy L Thompson // -- Set mesh file, polynomial degree, problem type 91ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr); 92ccaff030SJeremy L Thompson ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr); 93ccaff030SJeremy L Thompson numLevels = appCtx->numLevels; 94ccaff030SJeremy L Thompson fineLevel = numLevels - 1; 95ccaff030SJeremy L Thompson 96ccaff030SJeremy L Thompson // -- Set Poison's ratio, Young's Modulus 97ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr); 98ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 99ccaff030SJeremy L Thompson ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr); 100*f7b4142eSJeremy L Thompson if (fabs(appCtx->nuSmoother) > 1E-14) { 101*f7b4142eSJeremy L Thompson ierr = PetscMalloc1(1, &physSmoother); CHKERRQ(ierr); 102*f7b4142eSJeremy L Thompson ierr = PetscMemcpy(physSmoother, phys, sizeof(*phys)); CHKERRQ(ierr); 103*f7b4142eSJeremy L Thompson physSmoother->nu = appCtx->nuSmoother; 104*f7b4142eSJeremy L Thompson } 105ccaff030SJeremy L Thompson 106ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 107ccaff030SJeremy L Thompson // Setup DM 108ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 109ccaff030SJeremy L Thompson // Performance logging 110ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 111ccaff030SJeremy L Thompson CHKERRQ(ierr); 112ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 113ccaff030SJeremy L Thompson 114ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 115ccaff030SJeremy L Thompson ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 116ccaff030SJeremy L Thompson 117ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 118ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 119e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 120ccaff030SJeremy L Thompson ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 121ccaff030SJeremy L Thompson ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 122a3c02c40SJeremy L Thompson PETSC_TRUE, ncompu); CHKERRQ(ierr); 123a3c02c40SJeremy L Thompson // -- Label field components for viewing 124a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 125a3c02c40SJeremy L Thompson PetscSection section; 126a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(levelDMs[level], §ion); CHKERRQ(ierr); 127a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 128a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 129a3c02c40SJeremy L Thompson CHKERRQ(ierr); 130a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 131a3c02c40SJeremy L Thompson CHKERRQ(ierr); 132a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 133a3c02c40SJeremy L Thompson CHKERRQ(ierr); 134a3c02c40SJeremy L Thompson } 135a3c02c40SJeremy L Thompson 136a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 137a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr); 138a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel], 1395c25879aSJeremy L Thompson PETSC_FALSE, ncompe); CHKERRQ(ierr); 140a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr); 141a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel], 1425c25879aSJeremy L Thompson PETSC_FALSE, ncompd); CHKERRQ(ierr); 143a3c02c40SJeremy L Thompson { 144a3c02c40SJeremy L Thompson // -- Label field components for viewing 145a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 146a3c02c40SJeremy L Thompson PetscSection section; 147a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(dmDiagnostic, §ion); CHKERRQ(ierr); 148a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 1495c25879aSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "CondensedPressure"); 1505c25879aSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "StrainEnergyDensity"); 151a3c02c40SJeremy L Thompson CHKERRQ(ierr); 152ccaff030SJeremy L Thompson } 153ccaff030SJeremy L Thompson 154ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 155ccaff030SJeremy L Thompson // Setup solution and work vectors 156ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 157ccaff030SJeremy L Thompson // Allocate arrays 158ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 159ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 160ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 161ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 162ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 163ccaff030SJeremy L Thompson 164ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 165e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 166ccaff030SJeremy L Thompson // -- Create global unknown vector U 167ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 168ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 169ccaff030SJeremy L Thompson // Note: Local size for matShell 170ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 171ccaff030SJeremy L Thompson 172ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 173ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 174ccaff030SJeremy L Thompson // Note: local size for libCEED 175ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 176ccaff030SJeremy L Thompson } 177ccaff030SJeremy L Thompson 178ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 179ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 180ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 181ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 182ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 183ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 184ccaff030SJeremy L Thompson 185ccaff030SJeremy L Thompson // Performance logging 186ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 187ccaff030SJeremy L Thompson 188ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 189ccaff030SJeremy L Thompson // Set up libCEED 190ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 191ccaff030SJeremy L Thompson // Performance logging 192ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 193ccaff030SJeremy L Thompson CHKERRQ(ierr); 194ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 195ccaff030SJeremy L Thompson 196ccaff030SJeremy L Thompson // Initalize backend 197ccaff030SJeremy L Thompson CeedInit(appCtx->ceedResource, &ceed); 198ccaff030SJeremy L Thompson if (appCtx->degree > 4 && appCtx->ceedResourceFine[0]) 199ccaff030SJeremy L Thompson CeedInit(appCtx->ceedResourceFine, &ceedFine); 200ccaff030SJeremy L Thompson 201ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 202ccaff030SJeremy L Thompson CeedVector forceCeed; 203ccaff030SJeremy L Thompson CeedScalar *f; 204ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 205ccaff030SJeremy L Thompson ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 206ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 207ccaff030SJeremy L Thompson CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f); 208ccaff030SJeremy L Thompson } 209ccaff030SJeremy L Thompson 210ccaff030SJeremy L Thompson // -- Restriction and prolongation QFunction 211ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 212ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 213ccaff030SJeremy L Thompson &qfRestrict); 214ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 215ccaff030SJeremy L Thompson &qfProlong); 216ccaff030SJeremy L Thompson } 217ccaff030SJeremy L Thompson 218ccaff030SJeremy L Thompson // -- Setup libCEED objects 219ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 220ccaff030SJeremy L Thompson // ---- Setup residual evaluator and geometric information 221ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 222ccaff030SJeremy L Thompson { 223ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine); 224ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 225a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 226a3c02c40SJeremy L Thompson levelCeed, appCtx, phys, ceedData, fineLevel, 227a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 228a3c02c40SJeremy L Thompson forceCeed, qfRestrict, qfProlong); 229ccaff030SJeremy L Thompson CHKERRQ(ierr); 230ccaff030SJeremy L Thompson } 231ccaff030SJeremy L Thompson // ---- Setup Jacobian evaluator and prolongation/restriction 232e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 233ccaff030SJeremy L Thompson if (level != fineLevel) { 234ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 235ccaff030SJeremy L Thompson } 236ccaff030SJeremy L Thompson 2371c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 238ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine); 239ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 240ccaff030SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys, 241ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 242ccaff030SJeremy L Thompson Ulocsz[level], forceCeed, qfRestrict, 243ccaff030SJeremy L Thompson qfProlong); CHKERRQ(ierr); 244ccaff030SJeremy L Thompson } 245ccaff030SJeremy L Thompson 246ccaff030SJeremy L Thompson // Performance logging 247ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 248ccaff030SJeremy L Thompson 249ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 250ccaff030SJeremy L Thompson // Setup global forcing vector 251ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 252ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 253ccaff030SJeremy L Thompson 254ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 255ccaff030SJeremy L Thompson ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 256ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 257ccaff030SJeremy L Thompson CHKERRQ(ierr); 258ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 259ccaff030SJeremy L Thompson } 260ccaff030SJeremy L Thompson 261ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 262ccaff030SJeremy L Thompson // Print problem summary 263ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 264ccaff030SJeremy L Thompson if (!appCtx->testMode) { 265ccaff030SJeremy L Thompson const char *usedresource; 266ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 267ccaff030SJeremy L Thompson 268ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 26917fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 270ccaff030SJeremy L Thompson " libCEED:\n" 271ccaff030SJeremy L Thompson " libCEED Backend : %s\n", 272ccaff030SJeremy L Thompson usedresource); CHKERRQ(ierr); 273ccaff030SJeremy L Thompson 274ccaff030SJeremy L Thompson if (ceedFine) { 275ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 276ccaff030SJeremy L Thompson " libCEED Backend - high order : %s\n", 277ccaff030SJeremy L Thompson appCtx->ceedResourceFine); CHKERRQ(ierr); 278ccaff030SJeremy L Thompson } 279ccaff030SJeremy L Thompson 280ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 281ccaff030SJeremy L Thompson " Problem:\n" 282ccaff030SJeremy L Thompson " Problem Name : %s\n" 283ccaff030SJeremy L Thompson " Forcing Function : %s\n" 284ccaff030SJeremy L Thompson " Mesh:\n" 285ccaff030SJeremy L Thompson " File : %s\n" 286ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 287ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 288ccaff030SJeremy L Thompson " Global nodes : %D\n" 289ccaff030SJeremy L Thompson " Owned nodes : %D\n" 290ccaff030SJeremy L Thompson " DoF per node : %D\n" 291ccaff030SJeremy L Thompson " Multigrid:\n" 292ccaff030SJeremy L Thompson " Type : %s\n" 293ccaff030SJeremy L Thompson " Number of Levels : %d\n", 294ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 295ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 296ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 297ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 298ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 299f892e6d1Sjeremylt (appCtx->degree == 1 && 300f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 301f892e6d1Sjeremylt "Algebraic multigrid" : 302ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 303f892e6d1Sjeremylt (appCtx->degree == 1 || 304f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 305f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 306ccaff030SJeremy L Thompson 307ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 308e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 309ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 310ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 311ccaff030SJeremy L Thompson " Level %D (%s):\n" 312ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 313ccaff030SJeremy L Thompson " Global Nodes : %D\n" 314ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 315ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 316ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 317ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 318ccaff030SJeremy L Thompson CHKERRQ(ierr); 319ccaff030SJeremy L Thompson } 320ccaff030SJeremy L Thompson } 321ccaff030SJeremy L Thompson } 322ccaff030SJeremy L Thompson 323ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 324ccaff030SJeremy L Thompson // Setup SNES 325ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 326ccaff030SJeremy L Thompson // Performance logging 327ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 328ccaff030SJeremy L Thompson CHKERRQ(ierr); 329ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 330ccaff030SJeremy L Thompson 331ccaff030SJeremy L Thompson // Create SNES 332ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 333ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 334ccaff030SJeremy L Thompson 335ccaff030SJeremy L Thompson // -- Jacobian evaluators 336ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 337ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 338e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 339ccaff030SJeremy L Thompson // -- Jacobian context for level 340ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 341ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 342*f7b4142eSJeremy L Thompson Uloc[level], ceedData[level], ceed, phys, 343*f7b4142eSJeremy L Thompson physSmoother, jacobCtx[level]); CHKERRQ(ierr); 344ccaff030SJeremy L Thompson 345ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 346ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 347ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 348ccaff030SJeremy L Thompson CHKERRQ(ierr); 349ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 350ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 351ccaff030SJeremy L Thompson CHKERRQ(ierr); 352ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 353ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 354ccaff030SJeremy L Thompson 355ccaff030SJeremy L Thompson } 356e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 357ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 358ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 359ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 360ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 361ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 362ccaff030SJeremy L Thompson 363ccaff030SJeremy L Thompson // -- Residual evaluation function 364ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr); 365ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 366ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 367ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 368*f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 369ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 370ccaff030SJeremy L Thompson 371ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 372ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 373ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 374e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 375ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 376ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 377ccaff030SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level], 378ccaff030SJeremy L Thompson Ug[level], Uloc[level-1], Uloc[level], 379ccaff030SJeremy L Thompson ceedData[level-1], ceedData[level], ceed, 380ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 381ccaff030SJeremy L Thompson 382ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 383ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 384ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 385ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 386ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 387ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 388ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 389ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 390ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 391ccaff030SJeremy L Thompson CHKERRQ(ierr); 392ccaff030SJeremy L Thompson } 393ccaff030SJeremy L Thompson 394ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 395ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 396ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 397e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 398e3e3df41Sjeremylt // -- Jacobian Matrix 399e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 400e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 401e3e3df41Sjeremylt 402e3e3df41Sjeremylt if (appCtx->degree > 1) { 403ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 404ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 405ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 406ccaff030SJeremy L Thompson 407e3e3df41Sjeremylt // -- Jacobian function 408ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 409ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 410ccaff030SJeremy L Thompson 411ccaff030SJeremy L Thompson // -- Residual evaluation function 412ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 413ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 414ccaff030SJeremy L Thompson CHKERRQ(ierr); 415ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 416ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 417ccaff030SJeremy L Thompson 418ccaff030SJeremy L Thompson // -- Update formJacobCtx 419ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 420ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 421ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 422e3e3df41Sjeremylt } 423e3e3df41Sjeremylt } 424e3e3df41Sjeremylt 425e3e3df41Sjeremylt // Set Jacobian function 426e3e3df41Sjeremylt if (appCtx->degree > 1) { 427e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 428e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 429e3e3df41Sjeremylt } else { 430e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 431e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 432e3e3df41Sjeremylt CHKERRQ(ierr); 433e3e3df41Sjeremylt } 434ccaff030SJeremy L Thompson 435ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 436ccaff030SJeremy L Thompson // Setup KSP 437ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 438ccaff030SJeremy L Thompson { 439ccaff030SJeremy L Thompson PC pc; 440ccaff030SJeremy L Thompson KSP ksp; 441ccaff030SJeremy L Thompson 442ccaff030SJeremy L Thompson // -- KSP 443ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 444ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 445ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 446ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 447ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 448ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 449ccaff030SJeremy L Thompson 450ccaff030SJeremy L Thompson // -- Preconditioning 451ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 452ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 453ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 454ccaff030SJeremy L Thompson 455ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 456ccaff030SJeremy L Thompson // ---- No Multigrid 457ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 458ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 459f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 460f892e6d1Sjeremylt // ---- AMG for degree 1 461f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 462ccaff030SJeremy L Thompson } else { 463ccaff030SJeremy L Thompson // ---- PCMG 464ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 465ccaff030SJeremy L Thompson 466ccaff030SJeremy L Thompson // ------ PCMG levels 467ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 468e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 469ccaff030SJeremy L Thompson // -------- Smoother 470ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 471ccaff030SJeremy L Thompson PC pcSmoother; 472ccaff030SJeremy L Thompson 473ccaff030SJeremy L Thompson // ---------- Smoother KSP 474ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 475ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 476ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 477ccaff030SJeremy L Thompson 478ccaff030SJeremy L Thompson // ---------- Chebyshev options 479ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 480ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 481ccaff030SJeremy L Thompson CHKERRQ(ierr); 482ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 483ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 484ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 485ccaff030SJeremy L Thompson CHKERRQ(ierr); 486ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 487ccaff030SJeremy L Thompson CHKERRQ(ierr); 488ccaff030SJeremy L Thompson 489ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 490ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 491ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 492ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 493ccaff030SJeremy L Thompson 494ccaff030SJeremy L Thompson // -------- Work vector 495ccaff030SJeremy L Thompson if (level != fineLevel) { 496ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 497ccaff030SJeremy L Thompson } 498ccaff030SJeremy L Thompson 499ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 500ccaff030SJeremy L Thompson if (level > 0) { 501ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 502ccaff030SJeremy L Thompson CHKERRQ(ierr); 503ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 504ccaff030SJeremy L Thompson CHKERRQ(ierr); 505ccaff030SJeremy L Thompson } 506ccaff030SJeremy L Thompson } 507ccaff030SJeremy L Thompson 508ccaff030SJeremy L Thompson // ------ PCMG coarse solve 509ccaff030SJeremy L Thompson KSP kspCoarse; 510ccaff030SJeremy L Thompson PC pcCoarse; 511ccaff030SJeremy L Thompson 512ccaff030SJeremy L Thompson // -------- Coarse KSP 513ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 514ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 515ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 516ccaff030SJeremy L Thompson CHKERRQ(ierr); 517ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 518ccaff030SJeremy L Thompson 519ccaff030SJeremy L Thompson // -------- Coarse preconditioner 520ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 521ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 522ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 523ccaff030SJeremy L Thompson 524ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 525ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 526ccaff030SJeremy L Thompson 527ccaff030SJeremy L Thompson // ------ PCMG options 528ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 529ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 530ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 531ccaff030SJeremy L Thompson } 532ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 533ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 534ccaff030SJeremy L Thompson } 535038d0cd7Sjeremylt { 536038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 537481a4840SJed Brown SNESLineSearch linesearch; 538481a4840SJed Brown 539481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 540481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 541481a4840SJed Brown } 542481a4840SJed Brown 543ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 544ccaff030SJeremy L Thompson 545ccaff030SJeremy L Thompson // Performance logging 546ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 547ccaff030SJeremy L Thompson 548ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 549ccaff030SJeremy L Thompson // Set initial guess 550ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 551ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 552ccaff030SJeremy L Thompson 553ccaff030SJeremy L Thompson // View solution 554ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 555ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 556ccaff030SJeremy L Thompson } 557ccaff030SJeremy L Thompson 558ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 559ccaff030SJeremy L Thompson // Solve SNES 560ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 5615f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 5625f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 5635f24f028Sjeremylt CHKERRQ(ierr); 5645f24f028Sjeremylt 565ccaff030SJeremy L Thompson // Performance logging 566ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 567ccaff030SJeremy L Thompson CHKERRQ(ierr); 568ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 569ccaff030SJeremy L Thompson 570ccaff030SJeremy L Thompson // Timing 571ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 572ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 573ccaff030SJeremy L Thompson 574ccaff030SJeremy L Thompson // Solve for each load increment 5755f24f028Sjeremylt PetscInt increment; 5765f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 5775f24f028Sjeremylt // -- Log increment count 5785f24f028Sjeremylt if (snesMonitor) { 5795f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 5805f24f028Sjeremylt CHKERRQ(ierr); 5815f24f028Sjeremylt } 5825f24f028Sjeremylt 583ccaff030SJeremy L Thompson // -- Scale the problem 584ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 585ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 586ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 587ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 588ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 589ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 590ccaff030SJeremy L Thompson } 591ccaff030SJeremy L Thompson 592ccaff030SJeremy L Thompson // -- Solve 593ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 594ccaff030SJeremy L Thompson 595ccaff030SJeremy L Thompson // -- View solution 596560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 597ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 598ccaff030SJeremy L Thompson } 599ccaff030SJeremy L Thompson 600ccaff030SJeremy L Thompson // -- Update SNES iteration count 601ccaff030SJeremy L Thompson PetscInt its; 602ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson snesIts += its; 6043951731eSjeremylt 6053951731eSjeremylt // -- Check for divergence 6063951731eSjeremylt SNESConvergedReason reason; 6073951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6083951731eSjeremylt if (reason < 0) 6093951731eSjeremylt break; 610ccaff030SJeremy L Thompson } 611ccaff030SJeremy L Thompson 612ccaff030SJeremy L Thompson // Timing 613ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 614ccaff030SJeremy L Thompson 615ccaff030SJeremy L Thompson // Performance logging 616ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 617ccaff030SJeremy L Thompson 618560e6f00SJeremy L Thompson // View solution 619560e6f00SJeremy L Thompson if (appCtx->viewFinalSoln && !appCtx->viewSoln) { 620560e6f00SJeremy L Thompson ierr = ViewSolution(comm, U, increment - 1, 1.0); CHKERRQ(ierr); 621560e6f00SJeremy L Thompson } 622560e6f00SJeremy L Thompson 623ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 624ccaff030SJeremy L Thompson // Output summary 625ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 626ccaff030SJeremy L Thompson if (!appCtx->testMode) { 627ccaff030SJeremy L Thompson // -- SNES 628ccaff030SJeremy L Thompson SNESType snesType; 629ccaff030SJeremy L Thompson SNESConvergedReason reason; 630ccaff030SJeremy L Thompson PetscReal rnorm; 631ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 632ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 633ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 634ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 635ccaff030SJeremy L Thompson " SNES:\n" 636ccaff030SJeremy L Thompson " SNES Type : %s\n" 637ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 638ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 6395f24f028Sjeremylt " Completed Load Increments : %d\n" 640ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 641ccaff030SJeremy L Thompson " Final rnorm : %e\n", 642ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 6435f24f028Sjeremylt appCtx->numIncrements, increment - 1, 6445f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 645ccaff030SJeremy L Thompson 646ccaff030SJeremy L Thompson // -- KSP 647ccaff030SJeremy L Thompson KSP ksp; 648ccaff030SJeremy L Thompson KSPType kspType; 649ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 650ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 651ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 652ccaff030SJeremy L Thompson " Linear Solver:\n" 653ccaff030SJeremy L Thompson " KSP Type : %s\n", 654ccaff030SJeremy L Thompson kspType); CHKERRQ(ierr); 655ccaff030SJeremy L Thompson 656ccaff030SJeremy L Thompson // -- PC 657ccaff030SJeremy L Thompson PC pc; 658e3e3df41Sjeremylt PCType pcType; 659ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 660e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 661e3e3df41Sjeremylt ierr = PetscPrintf(comm, 662e3e3df41Sjeremylt " PC Type : %s\n", 663e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 664e3e3df41Sjeremylt 6652b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 666e3e3df41Sjeremylt PCMGType pcmgType; 667ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 668ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 669ccaff030SJeremy L Thompson " P-Multigrid:\n" 670ccaff030SJeremy L Thompson " PCMG Type : %s\n" 671ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 672ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 673ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 674ccaff030SJeremy L Thompson 675ccaff030SJeremy L Thompson // -- Coarse Solve 676ccaff030SJeremy L Thompson KSP kspCoarse; 677ccaff030SJeremy L Thompson PC pcCoarse; 678ccaff030SJeremy L Thompson PCType pcType; 679ccaff030SJeremy L Thompson 680ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 681ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 682ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 683ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 684ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 685ccaff030SJeremy L Thompson " Coarse Solve:\n" 686ccaff030SJeremy L Thompson " KSP Type : %s\n" 687ccaff030SJeremy L Thompson " PC Type : %s\n", 688ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 689ccaff030SJeremy L Thompson } 690ccaff030SJeremy L Thompson } 691ccaff030SJeremy L Thompson 692ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 693ccaff030SJeremy L Thompson // Compute solve time 694ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 695ccaff030SJeremy L Thompson if (!appCtx->testMode) { 696ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 697ccaff030SJeremy L Thompson CHKERRQ(ierr); 698ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 699ccaff030SJeremy L Thompson CHKERRQ(ierr); 700ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 701ccaff030SJeremy L Thompson " Performance:\n" 702ccaff030SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n", 703ccaff030SJeremy L Thompson maxTime, minTime); CHKERRQ(ierr); 704ccaff030SJeremy L Thompson } 705ccaff030SJeremy L Thompson 706ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 707ccaff030SJeremy L Thompson // Compute error 708ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 709ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 710ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 711ccaff030SJeremy L Thompson const CeedScalar *truearray; 712ccaff030SJeremy L Thompson Vec errorVec, trueVec; 713ccaff030SJeremy L Thompson 714ccaff030SJeremy L Thompson // -- Work vectors 715ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 716ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 717ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 718ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 719ccaff030SJeremy L Thompson 720ccaff030SJeremy L Thompson // -- Assemble global true solution vector 721ccaff030SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST, 722ccaff030SJeremy L Thompson &truearray); 723ccaff030SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr); 724ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 725ccaff030SJeremy L Thompson CHKERRQ(ierr); 726ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 727ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 728ccaff030SJeremy L Thompson 729ccaff030SJeremy L Thompson // -- Compute L2 error 730ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 731ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 732ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 733ccaff030SJeremy L Thompson l2Error /= l2Unorm; 734ccaff030SJeremy L Thompson 735ccaff030SJeremy L Thompson // -- Output 736ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 737ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 738ccaff030SJeremy L Thompson " L2 Error : %e\n", 739ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 740ccaff030SJeremy L Thompson } 741ccaff030SJeremy L Thompson 742ccaff030SJeremy L Thompson // -- Cleanup 743ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 744ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 7452d93065eSjeremylt } 7462d93065eSjeremylt 7472d93065eSjeremylt // --------------------------------------------------------------------------- 7482d93065eSjeremylt // Compute energy 7492d93065eSjeremylt // --------------------------------------------------------------------------- 7502d93065eSjeremylt if (!appCtx->testMode) { 7512d93065eSjeremylt // -- Compute L2 error 7522d93065eSjeremylt CeedScalar energy; 753a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 754a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 7552d93065eSjeremylt 7562d93065eSjeremylt // -- Output 7572d93065eSjeremylt ierr = PetscPrintf(comm, 7582d93065eSjeremylt " Strain Energy : %e\n", 7592d93065eSjeremylt energy); CHKERRQ(ierr); 760ccaff030SJeremy L Thompson } 761ccaff030SJeremy L Thompson 762ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 7635c25879aSJeremy L Thompson // Output diagnostic quantities 7645c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 7655c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 7665c25879aSJeremy L Thompson // -- Setup context 7675c25879aSJeremy L Thompson UserMult diagnosticCtx; 7685c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 7695c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 7705c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 7715c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 7725c25879aSJeremy L Thompson 7735c25879aSJeremy L Thompson // -- Compute and output 7745c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 7755c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 7765c25879aSJeremy L Thompson CHKERRQ(ierr); 7775c25879aSJeremy L Thompson 7785c25879aSJeremy L Thompson // -- Cleanup 7795c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 7805c25879aSJeremy L Thompson } 7815c25879aSJeremy L Thompson 7825c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 783ccaff030SJeremy L Thompson // Free objects 784ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 785ccaff030SJeremy L Thompson // Data in arrays per level 786e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 787ccaff030SJeremy L Thompson // Vectors 788ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 789ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 790ccaff030SJeremy L Thompson 791ccaff030SJeremy L Thompson // Jacobian matrix and data 792ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 793ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 794ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 795ccaff030SJeremy L Thompson 796ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 797ccaff030SJeremy L Thompson if (level > 0) { 798ccaff030SJeremy L Thompson ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr); 799ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 800ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 801ccaff030SJeremy L Thompson } 802ccaff030SJeremy L Thompson 803ccaff030SJeremy L Thompson // DM 804ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 805ccaff030SJeremy L Thompson 806ccaff030SJeremy L Thompson // libCEED objects 807ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 808ccaff030SJeremy L Thompson } 809ccaff030SJeremy L Thompson 810ccaff030SJeremy L Thompson // Arrays 811ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 812ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 813ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 814ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 815ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 817ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 818ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 819ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 821ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 822ccaff030SJeremy L Thompson 823ccaff030SJeremy L Thompson // libCEED objects 824ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfRestrict); 825ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfProlong); 826ccaff030SJeremy L Thompson CeedDestroy(&ceed); 827ccaff030SJeremy L Thompson CeedDestroy(&ceedFine); 828ccaff030SJeremy L Thompson 829ccaff030SJeremy L Thompson // PETSc objects 830ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 831ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 832ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 833ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 834ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 835ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 836ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 837ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 838ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 839a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 840a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 841ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 842ccaff030SJeremy L Thompson 843ccaff030SJeremy L Thompson // Structs 844ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 845ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 846ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 847ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 848ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 849*f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 850ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 851ccaff030SJeremy L Thompson 852ccaff030SJeremy L Thompson return PetscFinalize(); 853ccaff030SJeremy L Thompson } 854