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 51f7b4142eSJeremy 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 74*13fdad4bSJeremy L Thompson PetscInt ncompe = 1, ncompd = 5; // 1 energy output, 5 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); 100f7b4142eSJeremy L Thompson if (fabs(appCtx->nuSmoother) > 1E-14) { 101f7b4142eSJeremy L Thompson ierr = PetscMalloc1(1, &physSmoother); CHKERRQ(ierr); 102f7b4142eSJeremy L Thompson ierr = PetscMemcpy(physSmoother, phys, sizeof(*phys)); CHKERRQ(ierr); 103f7b4142eSJeremy L Thompson physSmoother->nu = appCtx->nuSmoother; 104f7b4142eSJeremy 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], 1428ca58ff3SJeremy L Thompson PETSC_FALSE, ncompu + 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); 1498ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 1508ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1518ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 1528ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1538ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 1548ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1558ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "CondensedPressure"); 1568ca58ff3SJeremy L Thompson CHKERRQ(ierr); 157*13fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "TraceE"); 158*13fdad4bSJeremy L Thompson CHKERRQ(ierr); 159*13fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 160*13fdad4bSJeremy L Thompson CHKERRQ(ierr); 161*13fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 162*13fdad4bSJeremy L Thompson CHKERRQ(ierr); 163*13fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 164a3c02c40SJeremy L Thompson CHKERRQ(ierr); 165ccaff030SJeremy L Thompson } 166ccaff030SJeremy L Thompson 167ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 168ccaff030SJeremy L Thompson // Setup solution and work vectors 169ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 170ccaff030SJeremy L Thompson // Allocate arrays 171ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 172ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 173ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 174ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 175ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 176ccaff030SJeremy L Thompson 177ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 178e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 179ccaff030SJeremy L Thompson // -- Create global unknown vector U 180ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 181ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 182ccaff030SJeremy L Thompson // Note: Local size for matShell 183ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 184ccaff030SJeremy L Thompson 185ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 186ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 187ccaff030SJeremy L Thompson // Note: local size for libCEED 188ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 189ccaff030SJeremy L Thompson } 190ccaff030SJeremy L Thompson 191ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 192ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 193ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 194ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 195ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 196ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 197ccaff030SJeremy L Thompson 198ccaff030SJeremy L Thompson // Performance logging 199ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 200ccaff030SJeremy L Thompson 201ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 202ccaff030SJeremy L Thompson // Set up libCEED 203ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 204ccaff030SJeremy L Thompson // Performance logging 205ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 206ccaff030SJeremy L Thompson CHKERRQ(ierr); 207ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 208ccaff030SJeremy L Thompson 209ccaff030SJeremy L Thompson // Initalize backend 210ccaff030SJeremy L Thompson CeedInit(appCtx->ceedResource, &ceed); 211ccaff030SJeremy L Thompson if (appCtx->degree > 4 && appCtx->ceedResourceFine[0]) 212ccaff030SJeremy L Thompson CeedInit(appCtx->ceedResourceFine, &ceedFine); 213ccaff030SJeremy L Thompson 214ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 215ccaff030SJeremy L Thompson CeedVector forceCeed; 216ccaff030SJeremy L Thompson CeedScalar *f; 217ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 218ccaff030SJeremy L Thompson ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 219ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 220ccaff030SJeremy L Thompson CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f); 221ccaff030SJeremy L Thompson } 222ccaff030SJeremy L Thompson 223ccaff030SJeremy L Thompson // -- Restriction and prolongation QFunction 224ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 225ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 226ccaff030SJeremy L Thompson &qfRestrict); 227ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 228ccaff030SJeremy L Thompson &qfProlong); 229ccaff030SJeremy L Thompson } 230ccaff030SJeremy L Thompson 231ccaff030SJeremy L Thompson // -- Setup libCEED objects 232ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 233ccaff030SJeremy L Thompson // ---- Setup residual evaluator and geometric information 234ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 235ccaff030SJeremy L Thompson { 236ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine); 237ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 238a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 239a3c02c40SJeremy L Thompson levelCeed, appCtx, phys, ceedData, fineLevel, 240a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 241a3c02c40SJeremy L Thompson forceCeed, qfRestrict, qfProlong); 242ccaff030SJeremy L Thompson CHKERRQ(ierr); 243ccaff030SJeremy L Thompson } 244ccaff030SJeremy L Thompson // ---- Setup Jacobian evaluator and prolongation/restriction 245e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 246ccaff030SJeremy L Thompson if (level != fineLevel) { 247ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 248ccaff030SJeremy L Thompson } 249ccaff030SJeremy L Thompson 2501c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 251ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine); 252ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 253ccaff030SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys, 254ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 255ccaff030SJeremy L Thompson Ulocsz[level], forceCeed, qfRestrict, 256ccaff030SJeremy L Thompson qfProlong); CHKERRQ(ierr); 257ccaff030SJeremy L Thompson } 258ccaff030SJeremy L Thompson 259ccaff030SJeremy L Thompson // Performance logging 260ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 261ccaff030SJeremy L Thompson 262ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 263ccaff030SJeremy L Thompson // Setup global forcing vector 264ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 265ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 266ccaff030SJeremy L Thompson 267ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 268ccaff030SJeremy L Thompson ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 269ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 270ccaff030SJeremy L Thompson CHKERRQ(ierr); 271ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 272ccaff030SJeremy L Thompson } 273ccaff030SJeremy L Thompson 274ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 275ccaff030SJeremy L Thompson // Print problem summary 276ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 277ccaff030SJeremy L Thompson if (!appCtx->testMode) { 278ccaff030SJeremy L Thompson const char *usedresource; 279ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 280ccaff030SJeremy L Thompson 281ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 28217fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 283ccaff030SJeremy L Thompson " libCEED:\n" 284ccaff030SJeremy L Thompson " libCEED Backend : %s\n", 285ccaff030SJeremy L Thompson usedresource); CHKERRQ(ierr); 286ccaff030SJeremy L Thompson 287ccaff030SJeremy L Thompson if (ceedFine) { 288ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 289ccaff030SJeremy L Thompson " libCEED Backend - high order : %s\n", 290ccaff030SJeremy L Thompson appCtx->ceedResourceFine); CHKERRQ(ierr); 291ccaff030SJeremy L Thompson } 292ccaff030SJeremy L Thompson 293ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 294ccaff030SJeremy L Thompson " Problem:\n" 295ccaff030SJeremy L Thompson " Problem Name : %s\n" 296ccaff030SJeremy L Thompson " Forcing Function : %s\n" 297ccaff030SJeremy L Thompson " Mesh:\n" 298ccaff030SJeremy L Thompson " File : %s\n" 299ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 300ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 301ccaff030SJeremy L Thompson " Global nodes : %D\n" 302ccaff030SJeremy L Thompson " Owned nodes : %D\n" 303ccaff030SJeremy L Thompson " DoF per node : %D\n" 304ccaff030SJeremy L Thompson " Multigrid:\n" 305ccaff030SJeremy L Thompson " Type : %s\n" 306ccaff030SJeremy L Thompson " Number of Levels : %d\n", 307ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 308ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 309ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 310ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 311ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 312f892e6d1Sjeremylt (appCtx->degree == 1 && 313f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 314f892e6d1Sjeremylt "Algebraic multigrid" : 315ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 316f892e6d1Sjeremylt (appCtx->degree == 1 || 317f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 318f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 319ccaff030SJeremy L Thompson 320ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 321e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 322ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 323ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 324ccaff030SJeremy L Thompson " Level %D (%s):\n" 325ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 326ccaff030SJeremy L Thompson " Global Nodes : %D\n" 327ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 328ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 329ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 330ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 331ccaff030SJeremy L Thompson CHKERRQ(ierr); 332ccaff030SJeremy L Thompson } 333ccaff030SJeremy L Thompson } 334ccaff030SJeremy L Thompson } 335ccaff030SJeremy L Thompson 336ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 337ccaff030SJeremy L Thompson // Setup SNES 338ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 339ccaff030SJeremy L Thompson // Performance logging 340ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 341ccaff030SJeremy L Thompson CHKERRQ(ierr); 342ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 343ccaff030SJeremy L Thompson 344ccaff030SJeremy L Thompson // Create SNES 345ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 346ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 347ccaff030SJeremy L Thompson 348ccaff030SJeremy L Thompson // -- Jacobian evaluators 349ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 350ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 351e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 352ccaff030SJeremy L Thompson // -- Jacobian context for level 353ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 354ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 355f7b4142eSJeremy L Thompson Uloc[level], ceedData[level], ceed, phys, 356f7b4142eSJeremy L Thompson physSmoother, jacobCtx[level]); CHKERRQ(ierr); 357ccaff030SJeremy L Thompson 358ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 359ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 360ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 361ccaff030SJeremy L Thompson CHKERRQ(ierr); 362ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 363ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 364ccaff030SJeremy L Thompson CHKERRQ(ierr); 365ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 366ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 367ccaff030SJeremy L Thompson 368ccaff030SJeremy L Thompson } 369e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 370ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 371ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 372ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 373ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 374ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 375ccaff030SJeremy L Thompson 376ccaff030SJeremy L Thompson // -- Residual evaluation function 377ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr); 378ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 379ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 380ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 381f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 382ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 383ccaff030SJeremy L Thompson 384ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 385ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 386ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 387e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 388ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 389ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 390ccaff030SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level], 391ccaff030SJeremy L Thompson Ug[level], Uloc[level-1], Uloc[level], 392ccaff030SJeremy L Thompson ceedData[level-1], ceedData[level], ceed, 393ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 394ccaff030SJeremy L Thompson 395ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 396ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 397ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 398ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 399ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 400ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 401ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 402ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 403ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 404ccaff030SJeremy L Thompson CHKERRQ(ierr); 405ccaff030SJeremy L Thompson } 406ccaff030SJeremy L Thompson 407ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 408ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 409ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 410e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 411e3e3df41Sjeremylt // -- Jacobian Matrix 412e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 413e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 414e3e3df41Sjeremylt 415e3e3df41Sjeremylt if (appCtx->degree > 1) { 416ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 417ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 418ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 419ccaff030SJeremy L Thompson 420e3e3df41Sjeremylt // -- Jacobian function 421ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 422ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 423ccaff030SJeremy L Thompson 424ccaff030SJeremy L Thompson // -- Residual evaluation function 425ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 426ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 427ccaff030SJeremy L Thompson CHKERRQ(ierr); 428ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 429ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 430ccaff030SJeremy L Thompson 431ccaff030SJeremy L Thompson // -- Update formJacobCtx 432ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 433ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 434ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 435e3e3df41Sjeremylt } 436e3e3df41Sjeremylt } 437e3e3df41Sjeremylt 438e3e3df41Sjeremylt // Set Jacobian function 439e3e3df41Sjeremylt if (appCtx->degree > 1) { 440e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 441e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 442e3e3df41Sjeremylt } else { 443e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 444e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 445e3e3df41Sjeremylt CHKERRQ(ierr); 446e3e3df41Sjeremylt } 447ccaff030SJeremy L Thompson 448ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 449ccaff030SJeremy L Thompson // Setup KSP 450ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 451ccaff030SJeremy L Thompson { 452ccaff030SJeremy L Thompson PC pc; 453ccaff030SJeremy L Thompson KSP ksp; 454ccaff030SJeremy L Thompson 455ccaff030SJeremy L Thompson // -- KSP 456ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 457ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 458ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 459ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 460ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 461ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 462ccaff030SJeremy L Thompson 463ccaff030SJeremy L Thompson // -- Preconditioning 464ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 465ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 466ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 467ccaff030SJeremy L Thompson 468ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 469ccaff030SJeremy L Thompson // ---- No Multigrid 470ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 471ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 472f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 473f892e6d1Sjeremylt // ---- AMG for degree 1 474f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 475ccaff030SJeremy L Thompson } else { 476ccaff030SJeremy L Thompson // ---- PCMG 477ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 478ccaff030SJeremy L Thompson 479ccaff030SJeremy L Thompson // ------ PCMG levels 480ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 481e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 482ccaff030SJeremy L Thompson // -------- Smoother 483ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 484ccaff030SJeremy L Thompson PC pcSmoother; 485ccaff030SJeremy L Thompson 486ccaff030SJeremy L Thompson // ---------- Smoother KSP 487ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 488ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 489ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 490ccaff030SJeremy L Thompson 491ccaff030SJeremy L Thompson // ---------- Chebyshev options 492ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 493ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 494ccaff030SJeremy L Thompson CHKERRQ(ierr); 495ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 496ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 497ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 498ccaff030SJeremy L Thompson CHKERRQ(ierr); 499ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 500ccaff030SJeremy L Thompson CHKERRQ(ierr); 501ccaff030SJeremy L Thompson 502ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 503ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 504ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 505ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 506ccaff030SJeremy L Thompson 507ccaff030SJeremy L Thompson // -------- Work vector 508ccaff030SJeremy L Thompson if (level != fineLevel) { 509ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 510ccaff030SJeremy L Thompson } 511ccaff030SJeremy L Thompson 512ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 513ccaff030SJeremy L Thompson if (level > 0) { 514ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 515ccaff030SJeremy L Thompson CHKERRQ(ierr); 516ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 517ccaff030SJeremy L Thompson CHKERRQ(ierr); 518ccaff030SJeremy L Thompson } 519ccaff030SJeremy L Thompson } 520ccaff030SJeremy L Thompson 521ccaff030SJeremy L Thompson // ------ PCMG coarse solve 522ccaff030SJeremy L Thompson KSP kspCoarse; 523ccaff030SJeremy L Thompson PC pcCoarse; 524ccaff030SJeremy L Thompson 525ccaff030SJeremy L Thompson // -------- Coarse KSP 526ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 527ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 528ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 529ccaff030SJeremy L Thompson CHKERRQ(ierr); 530ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 531ccaff030SJeremy L Thompson 532ccaff030SJeremy L Thompson // -------- Coarse preconditioner 533ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 534ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 535ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 536ccaff030SJeremy L Thompson 537ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 538ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 539ccaff030SJeremy L Thompson 540ccaff030SJeremy L Thompson // ------ PCMG options 541ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 542ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 543ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 544ccaff030SJeremy L Thompson } 545ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 546ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 547ccaff030SJeremy L Thompson } 548038d0cd7Sjeremylt { 549038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 550481a4840SJed Brown SNESLineSearch linesearch; 551481a4840SJed Brown 552481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 553481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 554481a4840SJed Brown } 555481a4840SJed Brown 556ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 557ccaff030SJeremy L Thompson 558ccaff030SJeremy L Thompson // Performance logging 559ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 560ccaff030SJeremy L Thompson 561ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 562ccaff030SJeremy L Thompson // Set initial guess 563ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 564ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 565ccaff030SJeremy L Thompson 566ccaff030SJeremy L Thompson // View solution 567ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 568ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 569ccaff030SJeremy L Thompson } 570ccaff030SJeremy L Thompson 571ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 572ccaff030SJeremy L Thompson // Solve SNES 573ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 5745f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 5755f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 5765f24f028Sjeremylt CHKERRQ(ierr); 5775f24f028Sjeremylt 578ccaff030SJeremy L Thompson // Performance logging 579ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 580ccaff030SJeremy L Thompson CHKERRQ(ierr); 581ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 582ccaff030SJeremy L Thompson 583ccaff030SJeremy L Thompson // Timing 584ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 585ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 586ccaff030SJeremy L Thompson 587ccaff030SJeremy L Thompson // Solve for each load increment 5885f24f028Sjeremylt PetscInt increment; 5895f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 5905f24f028Sjeremylt // -- Log increment count 5915f24f028Sjeremylt if (snesMonitor) { 5925f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 5935f24f028Sjeremylt CHKERRQ(ierr); 5945f24f028Sjeremylt } 5955f24f028Sjeremylt 596ccaff030SJeremy L Thompson // -- Scale the problem 597ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 598ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 599ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 600ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 601ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 602ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson } 604ccaff030SJeremy L Thompson 605ccaff030SJeremy L Thompson // -- Solve 606ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 607ccaff030SJeremy L Thompson 608ccaff030SJeremy L Thompson // -- View solution 609560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 610ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 611ccaff030SJeremy L Thompson } 612ccaff030SJeremy L Thompson 613ccaff030SJeremy L Thompson // -- Update SNES iteration count 614ccaff030SJeremy L Thompson PetscInt its; 615ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 616ccaff030SJeremy L Thompson snesIts += its; 6173951731eSjeremylt 6183951731eSjeremylt // -- Check for divergence 6193951731eSjeremylt SNESConvergedReason reason; 6203951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6213951731eSjeremylt if (reason < 0) 6223951731eSjeremylt break; 623ccaff030SJeremy L Thompson } 624ccaff030SJeremy L Thompson 625ccaff030SJeremy L Thompson // Timing 626ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 627ccaff030SJeremy L Thompson 628ccaff030SJeremy L Thompson // Performance logging 629ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 630ccaff030SJeremy L Thompson 631ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 632ccaff030SJeremy L Thompson // Output summary 633ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 634ccaff030SJeremy L Thompson if (!appCtx->testMode) { 635ccaff030SJeremy L Thompson // -- SNES 636ccaff030SJeremy L Thompson SNESType snesType; 637ccaff030SJeremy L Thompson SNESConvergedReason reason; 638ccaff030SJeremy L Thompson PetscReal rnorm; 639ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 640ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 641ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 642ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 643ccaff030SJeremy L Thompson " SNES:\n" 644ccaff030SJeremy L Thompson " SNES Type : %s\n" 645ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 646ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 6475f24f028Sjeremylt " Completed Load Increments : %d\n" 648ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 649ccaff030SJeremy L Thompson " Final rnorm : %e\n", 650ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 6515f24f028Sjeremylt appCtx->numIncrements, increment - 1, 6525f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 653ccaff030SJeremy L Thompson 654ccaff030SJeremy L Thompson // -- KSP 655ccaff030SJeremy L Thompson KSP ksp; 656ccaff030SJeremy L Thompson KSPType kspType; 657ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 658ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 659ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 660ccaff030SJeremy L Thompson " Linear Solver:\n" 661ccaff030SJeremy L Thompson " KSP Type : %s\n", 662ccaff030SJeremy L Thompson kspType); CHKERRQ(ierr); 663ccaff030SJeremy L Thompson 664ccaff030SJeremy L Thompson // -- PC 665ccaff030SJeremy L Thompson PC pc; 666e3e3df41Sjeremylt PCType pcType; 667ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 668e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 669e3e3df41Sjeremylt ierr = PetscPrintf(comm, 670e3e3df41Sjeremylt " PC Type : %s\n", 671e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 672e3e3df41Sjeremylt 6732b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 674e3e3df41Sjeremylt PCMGType pcmgType; 675ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 676ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 677ccaff030SJeremy L Thompson " P-Multigrid:\n" 678ccaff030SJeremy L Thompson " PCMG Type : %s\n" 679ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 680ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 681ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 682ccaff030SJeremy L Thompson 683ccaff030SJeremy L Thompson // -- Coarse Solve 684ccaff030SJeremy L Thompson KSP kspCoarse; 685ccaff030SJeremy L Thompson PC pcCoarse; 686ccaff030SJeremy L Thompson PCType pcType; 687ccaff030SJeremy L Thompson 688ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 689ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 690ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 691ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 692ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 693ccaff030SJeremy L Thompson " Coarse Solve:\n" 694ccaff030SJeremy L Thompson " KSP Type : %s\n" 695ccaff030SJeremy L Thompson " PC Type : %s\n", 696ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 697ccaff030SJeremy L Thompson } 698ccaff030SJeremy L Thompson } 699ccaff030SJeremy L Thompson 700ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 701ccaff030SJeremy L Thompson // Compute solve time 702ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 703ccaff030SJeremy L Thompson if (!appCtx->testMode) { 704ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 705ccaff030SJeremy L Thompson CHKERRQ(ierr); 706ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 707ccaff030SJeremy L Thompson CHKERRQ(ierr); 708ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 709ccaff030SJeremy L Thompson " Performance:\n" 710ccaff030SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n", 711ccaff030SJeremy L Thompson maxTime, minTime); CHKERRQ(ierr); 712ccaff030SJeremy L Thompson } 713ccaff030SJeremy L Thompson 714ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 715ccaff030SJeremy L Thompson // Compute error 716ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 717ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 718ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 719ccaff030SJeremy L Thompson const CeedScalar *truearray; 720ccaff030SJeremy L Thompson Vec errorVec, trueVec; 721ccaff030SJeremy L Thompson 722ccaff030SJeremy L Thompson // -- Work vectors 723ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 724ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 725ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 726ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 727ccaff030SJeremy L Thompson 728ccaff030SJeremy L Thompson // -- Assemble global true solution vector 729ccaff030SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST, 730ccaff030SJeremy L Thompson &truearray); 731ccaff030SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr); 732ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 733ccaff030SJeremy L Thompson CHKERRQ(ierr); 734ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 735ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 736ccaff030SJeremy L Thompson 737ccaff030SJeremy L Thompson // -- Compute L2 error 738ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 739ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 740ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 741ccaff030SJeremy L Thompson l2Error /= l2Unorm; 742ccaff030SJeremy L Thompson 743ccaff030SJeremy L Thompson // -- Output 744ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 745ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 746ccaff030SJeremy L Thompson " L2 Error : %e\n", 747ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 748ccaff030SJeremy L Thompson } 749ccaff030SJeremy L Thompson 750ccaff030SJeremy L Thompson // -- Cleanup 751ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 752ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 7532d93065eSjeremylt } 7542d93065eSjeremylt 7552d93065eSjeremylt // --------------------------------------------------------------------------- 7562d93065eSjeremylt // Compute energy 7572d93065eSjeremylt // --------------------------------------------------------------------------- 7582d93065eSjeremylt if (!appCtx->testMode) { 7592d93065eSjeremylt // -- Compute L2 error 7602d93065eSjeremylt CeedScalar energy; 761a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 762a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 7632d93065eSjeremylt 7642d93065eSjeremylt // -- Output 7652d93065eSjeremylt ierr = PetscPrintf(comm, 7662d93065eSjeremylt " Strain Energy : %e\n", 7672d93065eSjeremylt energy); CHKERRQ(ierr); 768ccaff030SJeremy L Thompson } 769ccaff030SJeremy L Thompson 770ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 7715c25879aSJeremy L Thompson // Output diagnostic quantities 7725c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 7735c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 7745c25879aSJeremy L Thompson // -- Setup context 7755c25879aSJeremy L Thompson UserMult diagnosticCtx; 7765c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 7775c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 7785c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 7795c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 7805c25879aSJeremy L Thompson 7815c25879aSJeremy L Thompson // -- Compute and output 7825c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 7835c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 7845c25879aSJeremy L Thompson CHKERRQ(ierr); 7855c25879aSJeremy L Thompson 7865c25879aSJeremy L Thompson // -- Cleanup 7875c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 7885c25879aSJeremy L Thompson } 7895c25879aSJeremy L Thompson 7905c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 791ccaff030SJeremy L Thompson // Free objects 792ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 793ccaff030SJeremy L Thompson // Data in arrays per level 794e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 795ccaff030SJeremy L Thompson // Vectors 796ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 797ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 798ccaff030SJeremy L Thompson 799ccaff030SJeremy L Thompson // Jacobian matrix and data 800ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 801ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 802ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 803ccaff030SJeremy L Thompson 804ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 805ccaff030SJeremy L Thompson if (level > 0) { 806ccaff030SJeremy L Thompson ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr); 807ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 808ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 809ccaff030SJeremy L Thompson } 810ccaff030SJeremy L Thompson 811ccaff030SJeremy L Thompson // DM 812ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 813ccaff030SJeremy L Thompson 814ccaff030SJeremy L Thompson // libCEED objects 815ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson } 817ccaff030SJeremy L Thompson 818ccaff030SJeremy L Thompson // Arrays 819ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 821ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 822ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 823ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 824ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 825ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 826ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 827ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 828ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 829ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 830ccaff030SJeremy L Thompson 831ccaff030SJeremy L Thompson // libCEED objects 832ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfRestrict); 833ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfProlong); 834ccaff030SJeremy L Thompson CeedDestroy(&ceed); 835ccaff030SJeremy L Thompson CeedDestroy(&ceedFine); 836ccaff030SJeremy L Thompson 837ccaff030SJeremy L Thompson // PETSc objects 838ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 839ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 840ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 841ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 842ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 843ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 844ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 845ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 846ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 847a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 848a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 849ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 850ccaff030SJeremy L Thompson 851ccaff030SJeremy L Thompson // Structs 852ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 853ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 854ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 855ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 856ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 857f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 858ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 859ccaff030SJeremy L Thompson 860ccaff030SJeremy L Thompson return PetscFinalize(); 861ccaff030SJeremy L Thompson } 862