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 32b04a8a52SLeila Ghaffari // ./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 // 36*da5d3694SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -degree 3 -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 7413fdad4bSJeremy 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 // --------------------------------------------------------------------------- 10762e9c006SJeremy L Thompson // Initalize libCEED 10862e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 10962e9c006SJeremy L Thompson // Initalize backend 11062e9c006SJeremy L Thompson CeedInit(appCtx->ceedResource, &ceed); 11162e9c006SJeremy L Thompson if (appCtx->degree > 4 && appCtx->ceedResourceFine[0]) 11262e9c006SJeremy L Thompson CeedInit(appCtx->ceedResourceFine, &ceedFine); 11362e9c006SJeremy L Thompson 11462e9c006SJeremy L Thompson // Check preferred MemType 11562e9c006SJeremy L Thompson CeedMemType memTypeBackend; 11662e9c006SJeremy L Thompson CeedGetPreferredMemType(ceed, &memTypeBackend); 11762e9c006SJeremy L Thompson if (!appCtx->setMemTypeRequest) 11862e9c006SJeremy L Thompson appCtx->memTypeRequested = memTypeBackend; 11962e9c006SJeremy L Thompson else if (!appCtx->petscHaveCuda && appCtx->memTypeRequested == CEED_MEM_DEVICE) 12062e9c006SJeremy L Thompson SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 12162e9c006SJeremy L Thompson "PETSc was not built with CUDA. " 12262e9c006SJeremy L Thompson "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 12362e9c006SJeremy L Thompson 12462e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 125ccaff030SJeremy L Thompson // Setup DM 126ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 127ccaff030SJeremy L Thompson // Performance logging 128ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 129ccaff030SJeremy L Thompson CHKERRQ(ierr); 130ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 131ccaff030SJeremy L Thompson 132ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 133ccaff030SJeremy L Thompson ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 134ccaff030SJeremy L Thompson 135ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 136ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 137e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 138ccaff030SJeremy L Thompson ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 13962e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 14062e9c006SJeremy L Thompson ierr = DMSetVecType(levelDMs[level], VECCUDA); CHKERRQ(ierr); 14162e9c006SJeremy L Thompson } 142ccaff030SJeremy L Thompson ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 143a3c02c40SJeremy L Thompson PETSC_TRUE, ncompu); CHKERRQ(ierr); 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(levelDMs[level], §ion); CHKERRQ(ierr); 148a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 149a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 150a3c02c40SJeremy L Thompson CHKERRQ(ierr); 151a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 152a3c02c40SJeremy L Thompson CHKERRQ(ierr); 153a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 154a3c02c40SJeremy L Thompson CHKERRQ(ierr); 155a3c02c40SJeremy L Thompson } 156a3c02c40SJeremy L Thompson 157a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 158a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr); 159a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel], 1605c25879aSJeremy L Thompson PETSC_FALSE, ncompe); CHKERRQ(ierr); 161a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr); 162a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel], 1638ca58ff3SJeremy L Thompson PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr); 16462e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 16562e9c006SJeremy L Thompson ierr = DMSetVecType(dmEnergy, VECCUDA); CHKERRQ(ierr); 16662e9c006SJeremy L Thompson ierr = DMSetVecType(dmDiagnostic, VECCUDA); CHKERRQ(ierr); 16762e9c006SJeremy L Thompson } 168a3c02c40SJeremy L Thompson { 169a3c02c40SJeremy L Thompson // -- Label field components for viewing 170a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 171a3c02c40SJeremy L Thompson PetscSection section; 172a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(dmDiagnostic, §ion); CHKERRQ(ierr); 173a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 1748ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 1758ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1768ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 1778ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1788ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 1798ca58ff3SJeremy L Thompson CHKERRQ(ierr); 180379387d4SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure"); 1818ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1827ab18febSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"); 18313fdad4bSJeremy L Thompson CHKERRQ(ierr); 18413fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 18513fdad4bSJeremy L Thompson CHKERRQ(ierr); 18613fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 18713fdad4bSJeremy L Thompson CHKERRQ(ierr); 18813fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 189a3c02c40SJeremy L Thompson CHKERRQ(ierr); 190ccaff030SJeremy L Thompson } 191ccaff030SJeremy L Thompson 192ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 193ccaff030SJeremy L Thompson // Setup solution and work vectors 194ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 195ccaff030SJeremy L Thompson // Allocate arrays 196ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 197ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 198ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 199ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 200ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 201ccaff030SJeremy L Thompson 202ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 203e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 204ccaff030SJeremy L Thompson // -- Create global unknown vector U 205ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 206ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 207ccaff030SJeremy L Thompson // Note: Local size for matShell 208ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 209ccaff030SJeremy L Thompson 210ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 211ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 212ccaff030SJeremy L Thompson // Note: local size for libCEED 213ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 214ccaff030SJeremy L Thompson } 215ccaff030SJeremy L Thompson 216ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 217ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 218ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 219ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 220ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 221ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 222ccaff030SJeremy L Thompson 223ccaff030SJeremy L Thompson // Performance logging 224ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 225ccaff030SJeremy L Thompson 226ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 227ccaff030SJeremy L Thompson // Set up libCEED 228ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 229ccaff030SJeremy L Thompson // Performance logging 230ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 231ccaff030SJeremy L Thompson CHKERRQ(ierr); 232ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 233ccaff030SJeremy L Thompson 234ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 235ccaff030SJeremy L Thompson CeedVector forceCeed; 236ccaff030SJeremy L Thompson CeedScalar *f; 237ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 23862e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 239ccaff030SJeremy L Thompson ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 24062e9c006SJeremy L Thompson } else { 24162e9c006SJeremy L Thompson ierr = VecCUDAGetArray(Floc, &f); CHKERRQ(ierr); 24262e9c006SJeremy L Thompson } 243ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 24462e9c006SJeremy L Thompson CeedVectorSetArray(forceCeed, appCtx->memTypeRequested, CEED_USE_POINTER, f); 245ccaff030SJeremy L Thompson } 246ccaff030SJeremy L Thompson 247ccaff030SJeremy L Thompson // -- Restriction and prolongation QFunction 248ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 249ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 250ccaff030SJeremy L Thompson &qfRestrict); 251ccaff030SJeremy L Thompson CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 252ccaff030SJeremy L Thompson &qfProlong); 253ccaff030SJeremy L Thompson } 254ccaff030SJeremy L Thompson 255ccaff030SJeremy L Thompson // -- Setup libCEED objects 256ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 257ccaff030SJeremy L Thompson // ---- Setup residual evaluator and geometric information 258ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 259ccaff030SJeremy L Thompson { 260ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine); 261ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 262a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 263a3c02c40SJeremy L Thompson levelCeed, appCtx, phys, ceedData, fineLevel, 264a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 265a3c02c40SJeremy L Thompson forceCeed, qfRestrict, qfProlong); 266ccaff030SJeremy L Thompson CHKERRQ(ierr); 267ccaff030SJeremy L Thompson } 268ccaff030SJeremy L Thompson // ---- Setup Jacobian evaluator and prolongation/restriction 269e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 270ccaff030SJeremy L Thompson if (level != fineLevel) { 271ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 272ccaff030SJeremy L Thompson } 273ccaff030SJeremy L Thompson 2741c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 275ccaff030SJeremy L Thompson bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine); 276ccaff030SJeremy L Thompson Ceed levelCeed = highOrder ? ceedFine : ceed; 277ccaff030SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys, 278ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 279ccaff030SJeremy L Thompson Ulocsz[level], forceCeed, qfRestrict, 280ccaff030SJeremy L Thompson qfProlong); CHKERRQ(ierr); 281ccaff030SJeremy L Thompson } 282ccaff030SJeremy L Thompson 283ccaff030SJeremy L Thompson // Performance logging 284ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 285ccaff030SJeremy L Thompson 286ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 287ccaff030SJeremy L Thompson // Setup global forcing vector 288ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 289ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 290ccaff030SJeremy L Thompson 291ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 29262e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 293ccaff030SJeremy L Thompson ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 29462e9c006SJeremy L Thompson } else { 29562e9c006SJeremy L Thompson ierr = VecCUDARestoreArray(Floc, &f); CHKERRQ(ierr); 29662e9c006SJeremy L Thompson } 297ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 298ccaff030SJeremy L Thompson CHKERRQ(ierr); 299ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 300ccaff030SJeremy L Thompson } 301ccaff030SJeremy L Thompson 302ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 303ccaff030SJeremy L Thompson // Print problem summary 304ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 305ccaff030SJeremy L Thompson if (!appCtx->testMode) { 306ccaff030SJeremy L Thompson const char *usedresource; 307ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 308ccaff030SJeremy L Thompson 309ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 31017fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 311ccaff030SJeremy L Thompson " libCEED:\n" 31262e9c006SJeremy L Thompson " libCEED Backend : %s\n" 31362e9c006SJeremy L Thompson " libCEED Backend MemType : %s\n" 31462e9c006SJeremy L Thompson " libCEED User Requested MemType : %s\n", 31562e9c006SJeremy L Thompson usedresource, CeedMemTypes[memTypeBackend], 31662e9c006SJeremy L Thompson (appCtx->setMemTypeRequest) ? 31762e9c006SJeremy L Thompson CeedMemTypes[appCtx->memTypeRequested] : "none"); 31862e9c006SJeremy L Thompson CHKERRQ(ierr); 319ccaff030SJeremy L Thompson 320ccaff030SJeremy L Thompson if (ceedFine) { 321ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 322ccaff030SJeremy L Thompson " libCEED Backend - high order : %s\n", 323ccaff030SJeremy L Thompson appCtx->ceedResourceFine); CHKERRQ(ierr); 324ccaff030SJeremy L Thompson } 325ccaff030SJeremy L Thompson 32662e9c006SJeremy L Thompson VecType vecType; 32762e9c006SJeremy L Thompson ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 32862e9c006SJeremy L Thompson ierr = PetscPrintf(comm, 32962e9c006SJeremy L Thompson " PETSc:\n" 33062e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 33162e9c006SJeremy L Thompson vecType); CHKERRQ(ierr); 33262e9c006SJeremy L Thompson 333ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 334ccaff030SJeremy L Thompson " Problem:\n" 335ccaff030SJeremy L Thompson " Problem Name : %s\n" 336ccaff030SJeremy L Thompson " Forcing Function : %s\n" 337ccaff030SJeremy L Thompson " Mesh:\n" 338ccaff030SJeremy L Thompson " File : %s\n" 339ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 340ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 341ccaff030SJeremy L Thompson " Global nodes : %D\n" 342ccaff030SJeremy L Thompson " Owned nodes : %D\n" 343ccaff030SJeremy L Thompson " DoF per node : %D\n" 344ccaff030SJeremy L Thompson " Multigrid:\n" 345ccaff030SJeremy L Thompson " Type : %s\n" 346ccaff030SJeremy L Thompson " Number of Levels : %d\n", 347ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 348ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 349ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 350ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 351ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 352f892e6d1Sjeremylt (appCtx->degree == 1 && 353f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 354f892e6d1Sjeremylt "Algebraic multigrid" : 355ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 356f892e6d1Sjeremylt (appCtx->degree == 1 || 357f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 358f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 359ccaff030SJeremy L Thompson 360ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 361e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 362ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 363ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 364ccaff030SJeremy L Thompson " Level %D (%s):\n" 365ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 366ccaff030SJeremy L Thompson " Global Nodes : %D\n" 367ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 368ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 369ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 370ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 371ccaff030SJeremy L Thompson CHKERRQ(ierr); 372ccaff030SJeremy L Thompson } 373ccaff030SJeremy L Thompson } 374ccaff030SJeremy L Thompson } 375ccaff030SJeremy L Thompson 376ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 377ccaff030SJeremy L Thompson // Setup SNES 378ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 379ccaff030SJeremy L Thompson // Performance logging 380ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 381ccaff030SJeremy L Thompson CHKERRQ(ierr); 382ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 383ccaff030SJeremy L Thompson 384ccaff030SJeremy L Thompson // Create SNES 385ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 386ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 387ccaff030SJeremy L Thompson 388ccaff030SJeremy L Thompson // -- Jacobian evaluators 389ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 390ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 391e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 392ccaff030SJeremy L Thompson // -- Jacobian context for level 393ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 394ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 395f7b4142eSJeremy L Thompson Uloc[level], ceedData[level], ceed, phys, 396f7b4142eSJeremy L Thompson physSmoother, jacobCtx[level]); CHKERRQ(ierr); 397ccaff030SJeremy L Thompson 398ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 399ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 400ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 401ccaff030SJeremy L Thompson CHKERRQ(ierr); 402ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 403ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 404ccaff030SJeremy L Thompson CHKERRQ(ierr); 405ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 406ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 407a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 408a3658badSJeremy L Thompson ierr = MatShellSetVecType(jacobMat[level], VECCUDA); CHKERRQ(ierr); 409a3658badSJeremy L Thompson } 410ccaff030SJeremy L Thompson } 411e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 412ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 413ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 414ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 415ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 416ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 417ccaff030SJeremy L Thompson 418ccaff030SJeremy L Thompson // -- Residual evaluation function 419ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr); 420ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 421ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 422ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 423f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 424ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 425ccaff030SJeremy L Thompson 426ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 427ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 428ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 429e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 430ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 431ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 43262e9c006SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, appCtx, levelDMs[level-1], 43362e9c006SJeremy L Thompson levelDMs[level], Ug[level], Uloc[level-1], 43462e9c006SJeremy L Thompson Uloc[level], ceedData[level-1], 43562e9c006SJeremy L Thompson ceedData[level], ceed, 436ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 437ccaff030SJeremy L Thompson 438ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 439ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 440ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 441ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 442ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 443ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 444ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 445ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 446ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 447ccaff030SJeremy L Thompson CHKERRQ(ierr); 448a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 449a3658badSJeremy L Thompson ierr = MatShellSetVecType(prolongRestrMat[level], VECCUDA); CHKERRQ(ierr); 450a3658badSJeremy L Thompson } 451ccaff030SJeremy L Thompson } 452ccaff030SJeremy L Thompson 453ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 454ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 455ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 456e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 457e3e3df41Sjeremylt // -- Jacobian Matrix 458e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 459e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 460e3e3df41Sjeremylt 461e3e3df41Sjeremylt if (appCtx->degree > 1) { 462ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 463ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 464ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 465ccaff030SJeremy L Thompson 466e3e3df41Sjeremylt // -- Jacobian function 467ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 468ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 469ccaff030SJeremy L Thompson 470ccaff030SJeremy L Thompson // -- Residual evaluation function 471ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 472ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 473ccaff030SJeremy L Thompson CHKERRQ(ierr); 474ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 475ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 476ccaff030SJeremy L Thompson 477ccaff030SJeremy L Thompson // -- Update formJacobCtx 478ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 479ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 480ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 481e3e3df41Sjeremylt } 482e3e3df41Sjeremylt } 483e3e3df41Sjeremylt 484e3e3df41Sjeremylt // Set Jacobian function 485e3e3df41Sjeremylt if (appCtx->degree > 1) { 486e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 487e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 488e3e3df41Sjeremylt } else { 489e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 490e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 491e3e3df41Sjeremylt CHKERRQ(ierr); 492e3e3df41Sjeremylt } 493ccaff030SJeremy L Thompson 494ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 495ccaff030SJeremy L Thompson // Setup KSP 496ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 497ccaff030SJeremy L Thompson { 498ccaff030SJeremy L Thompson PC pc; 499ccaff030SJeremy L Thompson KSP ksp; 500ccaff030SJeremy L Thompson 501ccaff030SJeremy L Thompson // -- KSP 502ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 503ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 504ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 505ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 506ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 507ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 508ccaff030SJeremy L Thompson 509ccaff030SJeremy L Thompson // -- Preconditioning 510ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 511ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 512ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 513ccaff030SJeremy L Thompson 514ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 515ccaff030SJeremy L Thompson // ---- No Multigrid 516ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 517ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 518f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 519f892e6d1Sjeremylt // ---- AMG for degree 1 520f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 521ccaff030SJeremy L Thompson } else { 522ccaff030SJeremy L Thompson // ---- PCMG 523ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 524ccaff030SJeremy L Thompson 525ccaff030SJeremy L Thompson // ------ PCMG levels 526ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 527e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 528ccaff030SJeremy L Thompson // -------- Smoother 529ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 530ccaff030SJeremy L Thompson PC pcSmoother; 531ccaff030SJeremy L Thompson 532ccaff030SJeremy L Thompson // ---------- Smoother KSP 533ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 534ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 535ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 536ccaff030SJeremy L Thompson 537ccaff030SJeremy L Thompson // ---------- Chebyshev options 538ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 539ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 540ccaff030SJeremy L Thompson CHKERRQ(ierr); 541ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 542ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 543ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 544ccaff030SJeremy L Thompson CHKERRQ(ierr); 545ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 546ccaff030SJeremy L Thompson CHKERRQ(ierr); 547ccaff030SJeremy L Thompson 548ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 549ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 550ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 551ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 552ccaff030SJeremy L Thompson 553ccaff030SJeremy L Thompson // -------- Work vector 554ccaff030SJeremy L Thompson if (level != fineLevel) { 555ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 556ccaff030SJeremy L Thompson } 557ccaff030SJeremy L Thompson 558ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 559ccaff030SJeremy L Thompson if (level > 0) { 560ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 561ccaff030SJeremy L Thompson CHKERRQ(ierr); 562ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 563ccaff030SJeremy L Thompson CHKERRQ(ierr); 564ccaff030SJeremy L Thompson } 565ccaff030SJeremy L Thompson } 566ccaff030SJeremy L Thompson 567ccaff030SJeremy L Thompson // ------ PCMG coarse solve 568ccaff030SJeremy L Thompson KSP kspCoarse; 569ccaff030SJeremy L Thompson PC pcCoarse; 570ccaff030SJeremy L Thompson 571ccaff030SJeremy L Thompson // -------- Coarse KSP 572ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 573ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 574ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 575ccaff030SJeremy L Thompson CHKERRQ(ierr); 576ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 577ccaff030SJeremy L Thompson 578ccaff030SJeremy L Thompson // -------- Coarse preconditioner 579ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 580ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 581ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 582ccaff030SJeremy L Thompson 583ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 584ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 585ccaff030SJeremy L Thompson 586ccaff030SJeremy L Thompson // ------ PCMG options 587ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 588ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 589ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 590ccaff030SJeremy L Thompson } 591ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 592ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 593ccaff030SJeremy L Thompson } 594038d0cd7Sjeremylt { 595038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 596481a4840SJed Brown SNESLineSearch linesearch; 597481a4840SJed Brown 598481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 599481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 600481a4840SJed Brown } 601481a4840SJed Brown 602ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson 604ccaff030SJeremy L Thompson // Performance logging 605ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 606ccaff030SJeremy L Thompson 607ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 608ccaff030SJeremy L Thompson // Set initial guess 609ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 610f81c27eaSJed Brown ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 611ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 612ccaff030SJeremy L Thompson 613ccaff030SJeremy L Thompson // View solution 614ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 615ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 616ccaff030SJeremy L Thompson } 617ccaff030SJeremy L Thompson 618ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 619ccaff030SJeremy L Thompson // Solve SNES 620ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 6215f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 6225f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 6235f24f028Sjeremylt CHKERRQ(ierr); 6245f24f028Sjeremylt 625ccaff030SJeremy L Thompson // Performance logging 626ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 627ccaff030SJeremy L Thompson CHKERRQ(ierr); 628ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 629ccaff030SJeremy L Thompson 630ccaff030SJeremy L Thompson // Timing 631ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 632ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 633ccaff030SJeremy L Thompson 634ccaff030SJeremy L Thompson // Solve for each load increment 6355f24f028Sjeremylt PetscInt increment; 6365f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 6375f24f028Sjeremylt // -- Log increment count 6385f24f028Sjeremylt if (snesMonitor) { 6395f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 6405f24f028Sjeremylt CHKERRQ(ierr); 6415f24f028Sjeremylt } 6425f24f028Sjeremylt 643ccaff030SJeremy L Thompson // -- Scale the problem 644ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 645ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 646ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 647ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 648ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 649ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 650ccaff030SJeremy L Thompson } 651ccaff030SJeremy L Thompson 652ccaff030SJeremy L Thompson // -- Solve 653ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 654ccaff030SJeremy L Thompson 655ccaff030SJeremy L Thompson // -- View solution 656560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 657ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 658ccaff030SJeremy L Thompson } 659ccaff030SJeremy L Thompson 660ccaff030SJeremy L Thompson // -- Update SNES iteration count 661ccaff030SJeremy L Thompson PetscInt its; 662ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 663ccaff030SJeremy L Thompson snesIts += its; 6643951731eSjeremylt 6653951731eSjeremylt // -- Check for divergence 6663951731eSjeremylt SNESConvergedReason reason; 6673951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6683951731eSjeremylt if (reason < 0) 6693951731eSjeremylt break; 670ccaff030SJeremy L Thompson } 671ccaff030SJeremy L Thompson 672ccaff030SJeremy L Thompson // Timing 673ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 674ccaff030SJeremy L Thompson 675ccaff030SJeremy L Thompson // Performance logging 676ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 677ccaff030SJeremy L Thompson 678ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 679ccaff030SJeremy L Thompson // Output summary 680ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 681ccaff030SJeremy L Thompson if (!appCtx->testMode) { 682ccaff030SJeremy L Thompson // -- SNES 683ccaff030SJeremy L Thompson SNESType snesType; 684ccaff030SJeremy L Thompson SNESConvergedReason reason; 685ccaff030SJeremy L Thompson PetscReal rnorm; 686ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 687ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 688ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 689ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 690ccaff030SJeremy L Thompson " SNES:\n" 691ccaff030SJeremy L Thompson " SNES Type : %s\n" 692ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 693ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 6945f24f028Sjeremylt " Completed Load Increments : %d\n" 695ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 696ccaff030SJeremy L Thompson " Final rnorm : %e\n", 697ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 6985f24f028Sjeremylt appCtx->numIncrements, increment - 1, 6995f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 700ccaff030SJeremy L Thompson 701ccaff030SJeremy L Thompson // -- KSP 702ccaff030SJeremy L Thompson KSP ksp; 703ccaff030SJeremy L Thompson KSPType kspType; 704ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 705ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 706ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 707ccaff030SJeremy L Thompson " Linear Solver:\n" 708ccaff030SJeremy L Thompson " KSP Type : %s\n", 709ccaff030SJeremy L Thompson kspType); CHKERRQ(ierr); 710ccaff030SJeremy L Thompson 711ccaff030SJeremy L Thompson // -- PC 712ccaff030SJeremy L Thompson PC pc; 713e3e3df41Sjeremylt PCType pcType; 714ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 715e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 716e3e3df41Sjeremylt ierr = PetscPrintf(comm, 717e3e3df41Sjeremylt " PC Type : %s\n", 718e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 719e3e3df41Sjeremylt 7202b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 721e3e3df41Sjeremylt PCMGType pcmgType; 722ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 723ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 724ccaff030SJeremy L Thompson " P-Multigrid:\n" 725ccaff030SJeremy L Thompson " PCMG Type : %s\n" 726ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 727ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 728ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 729ccaff030SJeremy L Thompson 730ccaff030SJeremy L Thompson // -- Coarse Solve 731ccaff030SJeremy L Thompson KSP kspCoarse; 732ccaff030SJeremy L Thompson PC pcCoarse; 733ccaff030SJeremy L Thompson PCType pcType; 734ccaff030SJeremy L Thompson 735ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 736ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 737ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 738ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 739ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 740ccaff030SJeremy L Thompson " Coarse Solve:\n" 741ccaff030SJeremy L Thompson " KSP Type : %s\n" 742ccaff030SJeremy L Thompson " PC Type : %s\n", 743ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 744ccaff030SJeremy L Thompson } 745ccaff030SJeremy L Thompson } 746ccaff030SJeremy L Thompson 747ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 748ccaff030SJeremy L Thompson // Compute solve time 749ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 750ccaff030SJeremy L Thompson if (!appCtx->testMode) { 751ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 752ccaff030SJeremy L Thompson CHKERRQ(ierr); 753ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 754ccaff030SJeremy L Thompson CHKERRQ(ierr); 755ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 756ccaff030SJeremy L Thompson " Performance:\n" 757ccaff030SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n", 758ccaff030SJeremy L Thompson maxTime, minTime); CHKERRQ(ierr); 759ccaff030SJeremy L Thompson } 760ccaff030SJeremy L Thompson 761ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 762ccaff030SJeremy L Thompson // Compute error 763ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 764ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 765ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 766ccaff030SJeremy L Thompson const CeedScalar *truearray; 767ccaff030SJeremy L Thompson Vec errorVec, trueVec; 768ccaff030SJeremy L Thompson 769ccaff030SJeremy L Thompson // -- Work vectors 770ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 771ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 772ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 773ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 774ccaff030SJeremy L Thompson 775ccaff030SJeremy L Thompson // -- Assemble global true solution vector 77662e9c006SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, 77762e9c006SJeremy L Thompson appCtx->memTypeRequested, &truearray); 77862e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 77962e9c006SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 78062e9c006SJeremy L Thompson CHKERRQ(ierr); 78162e9c006SJeremy L Thompson } else { 78262e9c006SJeremy L Thompson ierr = VecCUDAPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 78362e9c006SJeremy L Thompson CHKERRQ(ierr); 78462e9c006SJeremy L Thompson } 785ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 786ccaff030SJeremy L Thompson CHKERRQ(ierr); 78762e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 788ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 78962e9c006SJeremy L Thompson } else { 79062e9c006SJeremy L Thompson ierr = VecCUDAResetArray(resCtx->Yloc); CHKERRQ(ierr); 79162e9c006SJeremy L Thompson } 792ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 793ccaff030SJeremy L Thompson 794ccaff030SJeremy L Thompson // -- Compute L2 error 795ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 796ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 797ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 798ccaff030SJeremy L Thompson l2Error /= l2Unorm; 799ccaff030SJeremy L Thompson 800ccaff030SJeremy L Thompson // -- Output 801ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 802ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 803ccaff030SJeremy L Thompson " L2 Error : %e\n", 804ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 805ccaff030SJeremy L Thompson } 806ccaff030SJeremy L Thompson 807ccaff030SJeremy L Thompson // -- Cleanup 808ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 809ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 8102d93065eSjeremylt } 8112d93065eSjeremylt 8122d93065eSjeremylt // --------------------------------------------------------------------------- 8132d93065eSjeremylt // Compute energy 8142d93065eSjeremylt // --------------------------------------------------------------------------- 8152d93065eSjeremylt if (!appCtx->testMode) { 8162d93065eSjeremylt // -- Compute L2 error 8172d93065eSjeremylt CeedScalar energy; 818a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 819a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 8202d93065eSjeremylt 8212d93065eSjeremylt // -- Output 8222d93065eSjeremylt ierr = PetscPrintf(comm, 8232d93065eSjeremylt " Strain Energy : %e\n", 8242d93065eSjeremylt energy); CHKERRQ(ierr); 825ccaff030SJeremy L Thompson } 826ccaff030SJeremy L Thompson 827ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 8285c25879aSJeremy L Thompson // Output diagnostic quantities 8295c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 8305c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 8315c25879aSJeremy L Thompson // -- Setup context 8325c25879aSJeremy L Thompson UserMult diagnosticCtx; 8335c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 8345c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 8355c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 8365c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 8375c25879aSJeremy L Thompson 8385c25879aSJeremy L Thompson // -- Compute and output 8395c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 8405c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 8415c25879aSJeremy L Thompson CHKERRQ(ierr); 8425c25879aSJeremy L Thompson 8435c25879aSJeremy L Thompson // -- Cleanup 8445c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 8455c25879aSJeremy L Thompson } 8465c25879aSJeremy L Thompson 8475c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 848ccaff030SJeremy L Thompson // Free objects 849ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 850ccaff030SJeremy L Thompson // Data in arrays per level 851e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 852ccaff030SJeremy L Thompson // Vectors 853ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 854ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 855ccaff030SJeremy L Thompson 856ccaff030SJeremy L Thompson // Jacobian matrix and data 857ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 858ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 859ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 860ccaff030SJeremy L Thompson 861ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 862ccaff030SJeremy L Thompson if (level > 0) { 863ccaff030SJeremy L Thompson ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr); 864ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 865ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 866ccaff030SJeremy L Thompson } 867ccaff030SJeremy L Thompson 868ccaff030SJeremy L Thompson // DM 869ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 870ccaff030SJeremy L Thompson 871ccaff030SJeremy L Thompson // libCEED objects 872ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 873ccaff030SJeremy L Thompson } 874ccaff030SJeremy L Thompson 875ccaff030SJeremy L Thompson // Arrays 876ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 877ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 878ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 879ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 880ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 881ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 882ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 883ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 884ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 885ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 886ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 887ccaff030SJeremy L Thompson 888ccaff030SJeremy L Thompson // libCEED objects 889ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfRestrict); 890ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qfProlong); 891ccaff030SJeremy L Thompson CeedDestroy(&ceed); 892ccaff030SJeremy L Thompson CeedDestroy(&ceedFine); 893ccaff030SJeremy L Thompson 894ccaff030SJeremy L Thompson // PETSc objects 895ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 896ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 897ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 898ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 899ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 900ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 901ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 902ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 903ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 904a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 905a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 906ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 907ccaff030SJeremy L Thompson 908ccaff030SJeremy L Thompson // Structs 909ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 910ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 911ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 912ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 913ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 914f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 915ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 916ccaff030SJeremy L Thompson 917ccaff030SJeremy L Thompson return PetscFinalize(); 918ccaff030SJeremy L Thompson } 919