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