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 // 36da5d3694SJeremy 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 69d99fa3c5SJeremy L Thompson Ceed ceed; 70ccaff030SJeremy L Thompson CeedData *ceedData; 71*777ff853SJeremy L Thompson CeedQFunctionContext ctxPhys, ctxPhysSmoother = 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 777418ca88SJeremy L Thompson PetscInt snesIts = 0, kspIts = 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 11262e9c006SJeremy L Thompson // Check preferred MemType 11362e9c006SJeremy L Thompson CeedMemType memTypeBackend; 11462e9c006SJeremy L Thompson CeedGetPreferredMemType(ceed, &memTypeBackend); 11562e9c006SJeremy L Thompson if (!appCtx->setMemTypeRequest) 11662e9c006SJeremy L Thompson appCtx->memTypeRequested = memTypeBackend; 11762e9c006SJeremy L Thompson else if (!appCtx->petscHaveCuda && appCtx->memTypeRequested == CEED_MEM_DEVICE) 11862e9c006SJeremy L Thompson SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 11962e9c006SJeremy L Thompson "PETSc was not built with CUDA. " 12062e9c006SJeremy L Thompson "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 12162e9c006SJeremy L Thompson 122*777ff853SJeremy L Thompson // Wrap context in libCEED objects 123*777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhys); 124*777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhys, CEED_MEM_HOST, CEED_USE_POINTER, 125*777ff853SJeremy L Thompson sizeof(*phys), phys); 126*777ff853SJeremy L Thompson if (physSmoother) { 127*777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhysSmoother); 128*777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhysSmoother, CEED_MEM_HOST, CEED_USE_POINTER, 129*777ff853SJeremy L Thompson sizeof(*physSmoother), physSmoother); 130*777ff853SJeremy L Thompson } 131*777ff853SJeremy L Thompson 13262e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 133ccaff030SJeremy L Thompson // Setup DM 134ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 135ccaff030SJeremy L Thompson // Performance logging 136ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 137ccaff030SJeremy L Thompson CHKERRQ(ierr); 138ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 139ccaff030SJeremy L Thompson 140ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 141ccaff030SJeremy L Thompson ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 142ccaff030SJeremy L Thompson 143ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 144ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 145e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 146ccaff030SJeremy L Thompson ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 14762e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 14862e9c006SJeremy L Thompson ierr = DMSetVecType(levelDMs[level], VECCUDA); CHKERRQ(ierr); 14962e9c006SJeremy L Thompson } 150ccaff030SJeremy L Thompson ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 151a3c02c40SJeremy L Thompson PETSC_TRUE, ncompu); CHKERRQ(ierr); 152a3c02c40SJeremy L Thompson // -- Label field components for viewing 153a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 154a3c02c40SJeremy L Thompson PetscSection section; 155a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(levelDMs[level], §ion); CHKERRQ(ierr); 156a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 157a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 158a3c02c40SJeremy L Thompson CHKERRQ(ierr); 159a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 160a3c02c40SJeremy L Thompson CHKERRQ(ierr); 161a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 162a3c02c40SJeremy L Thompson CHKERRQ(ierr); 163a3c02c40SJeremy L Thompson } 164a3c02c40SJeremy L Thompson 165a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 166a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr); 167a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel], 1685c25879aSJeremy L Thompson PETSC_FALSE, ncompe); CHKERRQ(ierr); 169a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr); 170a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel], 1718ca58ff3SJeremy L Thompson PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr); 17262e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 17362e9c006SJeremy L Thompson ierr = DMSetVecType(dmEnergy, VECCUDA); CHKERRQ(ierr); 17462e9c006SJeremy L Thompson ierr = DMSetVecType(dmDiagnostic, VECCUDA); CHKERRQ(ierr); 17562e9c006SJeremy L Thompson } 176a3c02c40SJeremy L Thompson { 177a3c02c40SJeremy L Thompson // -- Label field components for viewing 178a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 179a3c02c40SJeremy L Thompson PetscSection section; 180a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(dmDiagnostic, §ion); CHKERRQ(ierr); 181a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 1828ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 1838ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1848ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 1858ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1868ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 1878ca58ff3SJeremy L Thompson CHKERRQ(ierr); 188379387d4SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure"); 1898ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1907ab18febSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"); 19113fdad4bSJeremy L Thompson CHKERRQ(ierr); 19213fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 19313fdad4bSJeremy L Thompson CHKERRQ(ierr); 19413fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 19513fdad4bSJeremy L Thompson CHKERRQ(ierr); 19613fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 197a3c02c40SJeremy L Thompson CHKERRQ(ierr); 198ccaff030SJeremy L Thompson } 199ccaff030SJeremy L Thompson 200ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 201ccaff030SJeremy L Thompson // Setup solution and work vectors 202ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 203ccaff030SJeremy L Thompson // Allocate arrays 204ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 205ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 206ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 207ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 208ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 209ccaff030SJeremy L Thompson 210ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 211e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 212ccaff030SJeremy L Thompson // -- Create global unknown vector U 213ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 214ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 215ccaff030SJeremy L Thompson // Note: Local size for matShell 216ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 217ccaff030SJeremy L Thompson 218ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 219ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 220ccaff030SJeremy L Thompson // Note: local size for libCEED 221ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 222ccaff030SJeremy L Thompson } 223ccaff030SJeremy L Thompson 224ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 225ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 226ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 227ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 228ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 229ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 230ccaff030SJeremy L Thompson 231ccaff030SJeremy L Thompson // Performance logging 232ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 233ccaff030SJeremy L Thompson 234ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 235ccaff030SJeremy L Thompson // Set up libCEED 236ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 237ccaff030SJeremy L Thompson // Performance logging 238ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 239ccaff030SJeremy L Thompson CHKERRQ(ierr); 240ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 241ccaff030SJeremy L Thompson 242ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 243ccaff030SJeremy L Thompson CeedVector forceCeed; 244ccaff030SJeremy L Thompson CeedScalar *f; 245ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 24662e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 247ccaff030SJeremy L Thompson ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 24862e9c006SJeremy L Thompson } else { 24962e9c006SJeremy L Thompson ierr = VecCUDAGetArray(Floc, &f); CHKERRQ(ierr); 25062e9c006SJeremy L Thompson } 251ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 25262e9c006SJeremy L Thompson CeedVectorSetArray(forceCeed, appCtx->memTypeRequested, CEED_USE_POINTER, f); 253ccaff030SJeremy L Thompson } 254ccaff030SJeremy L Thompson 255ccaff030SJeremy L Thompson // -- Setup libCEED objects 256ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 257d99fa3c5SJeremy L Thompson // ---- Setup residual, Jacobian evaluator and geometric information 258ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 259ccaff030SJeremy L Thompson { 260a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 261*777ff853SJeremy L Thompson ceed, appCtx, ctxPhys, ceedData, fineLevel, 262a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 263*777ff853SJeremy L Thompson forceCeed); 264ccaff030SJeremy L Thompson CHKERRQ(ierr); 265ccaff030SJeremy L Thompson } 266d99fa3c5SJeremy L Thompson // ---- Setup coarse Jacobian evaluator and prolongation/restriction 267d99fa3c5SJeremy L Thompson for (PetscInt level = numLevels - 2; level >= 0; level--) { 268ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 269d99fa3c5SJeremy L Thompson 270d99fa3c5SJeremy L Thompson // Get global communication restriction 271d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 272d99fa3c5SJeremy L Thompson ierr = VecSet(Uloc[level+1], 1.0); CHKERRQ(ierr); 273d99fa3c5SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[level+1], Uloc[level+1], ADD_VALUES, 274d99fa3c5SJeremy L Thompson Ug[level+1]); CHKERRQ(ierr); 275d99fa3c5SJeremy L Thompson ierr = DMGlobalToLocal(levelDMs[level+1], Ug[level+1], INSERT_VALUES, 276d99fa3c5SJeremy L Thompson Uloc[level+1]); CHKERRQ(ierr); 277d99fa3c5SJeremy L Thompson 278d99fa3c5SJeremy L Thompson // Place in libCEED array 279d99fa3c5SJeremy L Thompson const PetscScalar *m; 280d99fa3c5SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 281d99fa3c5SJeremy L Thompson ierr = VecGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 282d99fa3c5SJeremy L Thompson } else { 283d99fa3c5SJeremy L Thompson ierr = VecCUDAGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 284ccaff030SJeremy L Thompson } 285d99fa3c5SJeremy L Thompson CeedVectorSetArray(ceedData[level+1]->xceed, appCtx->memTypeRequested, 286d99fa3c5SJeremy L Thompson CEED_USE_POINTER, (CeedScalar *)m); 287ccaff030SJeremy L Thompson 2881c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 289*777ff853SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], ceed, appCtx, 290ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 291*777ff853SJeremy L Thompson Ulocsz[level], ceedData[level+1]->xceed); 292*777ff853SJeremy L Thompson CHKERRQ(ierr); 293d99fa3c5SJeremy L Thompson 294d99fa3c5SJeremy L Thompson // Restore PETSc vector 295d99fa3c5SJeremy L Thompson CeedVectorTakeArray(ceedData[level+1]->xceed, appCtx->memTypeRequested, 296d99fa3c5SJeremy L Thompson (CeedScalar **)&m); 297d99fa3c5SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 298d99fa3c5SJeremy L Thompson ierr = VecRestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 299d99fa3c5SJeremy L Thompson } else { 300d99fa3c5SJeremy L Thompson ierr = VecCUDARestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 301d99fa3c5SJeremy L Thompson } 302d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 303d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Uloc[level+1]); CHKERRQ(ierr); 304ccaff030SJeremy L Thompson } 305ccaff030SJeremy L Thompson 306ccaff030SJeremy L Thompson // Performance logging 307ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 308ccaff030SJeremy L Thompson 309ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 310ccaff030SJeremy L Thompson // Setup global forcing vector 311ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 312ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 313ccaff030SJeremy L Thompson 314ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 3156a6c615bSJeremy L Thompson CeedVectorTakeArray(forceCeed, appCtx->memTypeRequested, NULL); 31662e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 317ccaff030SJeremy L Thompson ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 31862e9c006SJeremy L Thompson } else { 31962e9c006SJeremy L Thompson ierr = VecCUDARestoreArray(Floc, &f); CHKERRQ(ierr); 32062e9c006SJeremy L Thompson } 321ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 322ccaff030SJeremy L Thompson CHKERRQ(ierr); 323ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 324ccaff030SJeremy L Thompson } 325ccaff030SJeremy L Thompson 326ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 327ccaff030SJeremy L Thompson // Print problem summary 328ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 329ccaff030SJeremy L Thompson if (!appCtx->testMode) { 330ccaff030SJeremy L Thompson const char *usedresource; 331ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 332ccaff030SJeremy L Thompson 333ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 33417fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 335ccaff030SJeremy L Thompson " libCEED:\n" 33662e9c006SJeremy L Thompson " libCEED Backend : %s\n" 33762e9c006SJeremy L Thompson " libCEED Backend MemType : %s\n" 33862e9c006SJeremy L Thompson " libCEED User Requested MemType : %s\n", 33962e9c006SJeremy L Thompson usedresource, CeedMemTypes[memTypeBackend], 34062e9c006SJeremy L Thompson (appCtx->setMemTypeRequest) ? 34162e9c006SJeremy L Thompson CeedMemTypes[appCtx->memTypeRequested] : "none"); 34262e9c006SJeremy L Thompson CHKERRQ(ierr); 343ccaff030SJeremy L Thompson 34462e9c006SJeremy L Thompson VecType vecType; 34562e9c006SJeremy L Thompson ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 34662e9c006SJeremy L Thompson ierr = PetscPrintf(comm, 34762e9c006SJeremy L Thompson " PETSc:\n" 34862e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 34962e9c006SJeremy L Thompson vecType); CHKERRQ(ierr); 35062e9c006SJeremy L Thompson 351ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 352ccaff030SJeremy L Thompson " Problem:\n" 353ccaff030SJeremy L Thompson " Problem Name : %s\n" 354ccaff030SJeremy L Thompson " Forcing Function : %s\n" 355ccaff030SJeremy L Thompson " Mesh:\n" 356ccaff030SJeremy L Thompson " File : %s\n" 357ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 358ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 359ccaff030SJeremy L Thompson " Global nodes : %D\n" 360ccaff030SJeremy L Thompson " Owned nodes : %D\n" 361ccaff030SJeremy L Thompson " DoF per node : %D\n" 362ccaff030SJeremy L Thompson " Multigrid:\n" 363ccaff030SJeremy L Thompson " Type : %s\n" 364ccaff030SJeremy L Thompson " Number of Levels : %d\n", 365ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 366ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 367ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 368ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 369ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 370f892e6d1Sjeremylt (appCtx->degree == 1 && 371f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 372f892e6d1Sjeremylt "Algebraic multigrid" : 373ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 374f892e6d1Sjeremylt (appCtx->degree == 1 || 375f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 376f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 377ccaff030SJeremy L Thompson 378ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 379e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 380ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 381ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 382ccaff030SJeremy L Thompson " Level %D (%s):\n" 383ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 384ccaff030SJeremy L Thompson " Global Nodes : %D\n" 385ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 386ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 387ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 388ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 389ccaff030SJeremy L Thompson CHKERRQ(ierr); 390ccaff030SJeremy L Thompson } 391ccaff030SJeremy L Thompson } 392ccaff030SJeremy L Thompson } 393ccaff030SJeremy L Thompson 394ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 395ccaff030SJeremy L Thompson // Setup SNES 396ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 397ccaff030SJeremy L Thompson // Performance logging 398ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 399ccaff030SJeremy L Thompson CHKERRQ(ierr); 400ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 401ccaff030SJeremy L Thompson 402ccaff030SJeremy L Thompson // Create SNES 403ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 404ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 405ccaff030SJeremy L Thompson 406ccaff030SJeremy L Thompson // -- Jacobian evaluators 407ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 408ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 409e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 410ccaff030SJeremy L Thompson // -- Jacobian context for level 411ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 412ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 413*777ff853SJeremy L Thompson Uloc[level], ceedData[level], ceed, ctxPhys, 414*777ff853SJeremy L Thompson ctxPhysSmoother, jacobCtx[level]); CHKERRQ(ierr); 415ccaff030SJeremy L Thompson 416ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 417ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 418ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 419ccaff030SJeremy L Thompson CHKERRQ(ierr); 420ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 421ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 422ccaff030SJeremy L Thompson CHKERRQ(ierr); 423ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 424ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 425a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 426a3658badSJeremy L Thompson ierr = MatShellSetVecType(jacobMat[level], VECCUDA); CHKERRQ(ierr); 427a3658badSJeremy L Thompson } 428ccaff030SJeremy L Thompson } 429e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 430ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 431ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 432ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 433ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 434ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 435ccaff030SJeremy L Thompson 436ccaff030SJeremy L Thompson // -- Residual evaluation function 437ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr); 438ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 439ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 440ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 441f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 442ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 443ccaff030SJeremy L Thompson 444ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 445ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 446ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 447e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 448ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 449ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 45062e9c006SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, appCtx, levelDMs[level-1], 45162e9c006SJeremy L Thompson levelDMs[level], Ug[level], Uloc[level-1], 45262e9c006SJeremy L Thompson Uloc[level], ceedData[level-1], 45362e9c006SJeremy L Thompson ceedData[level], ceed, 454ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 455ccaff030SJeremy L Thompson 456ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 457ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 458ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 459ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 460ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 461ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 462ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 463ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 464ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 465ccaff030SJeremy L Thompson CHKERRQ(ierr); 466a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 467a3658badSJeremy L Thompson ierr = MatShellSetVecType(prolongRestrMat[level], VECCUDA); CHKERRQ(ierr); 468a3658badSJeremy L Thompson } 469ccaff030SJeremy L Thompson } 470ccaff030SJeremy L Thompson 471ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 472ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 473ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 474e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 475e3e3df41Sjeremylt // -- Jacobian Matrix 476e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 477e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 478e3e3df41Sjeremylt 479e3e3df41Sjeremylt if (appCtx->degree > 1) { 480ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 481ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 482ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 483ccaff030SJeremy L Thompson 484e3e3df41Sjeremylt // -- Jacobian function 485ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 486ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 487ccaff030SJeremy L Thompson 488ccaff030SJeremy L Thompson // -- Residual evaluation function 489ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 490ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 491ccaff030SJeremy L Thompson CHKERRQ(ierr); 492ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 493ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 494ccaff030SJeremy L Thompson 495ccaff030SJeremy L Thompson // -- Update formJacobCtx 496ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 497ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 498ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 499e3e3df41Sjeremylt } 500e3e3df41Sjeremylt } 501e3e3df41Sjeremylt 502e3e3df41Sjeremylt // Set Jacobian function 503e3e3df41Sjeremylt if (appCtx->degree > 1) { 504e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 505e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 506e3e3df41Sjeremylt } else { 507e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 508e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 509e3e3df41Sjeremylt CHKERRQ(ierr); 510e3e3df41Sjeremylt } 511ccaff030SJeremy L Thompson 512ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 513ccaff030SJeremy L Thompson // Setup KSP 514ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 515ccaff030SJeremy L Thompson { 516ccaff030SJeremy L Thompson PC pc; 517ccaff030SJeremy L Thompson KSP ksp; 518ccaff030SJeremy L Thompson 519ccaff030SJeremy L Thompson // -- KSP 520ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 521ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 522ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 523ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 524ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 525ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 526ccaff030SJeremy L Thompson 527ccaff030SJeremy L Thompson // -- Preconditioning 528ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 529ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 530ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 531ccaff030SJeremy L Thompson 532ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 533ccaff030SJeremy L Thompson // ---- No Multigrid 534ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 535ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 536f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 537f892e6d1Sjeremylt // ---- AMG for degree 1 538f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 539ccaff030SJeremy L Thompson } else { 540ccaff030SJeremy L Thompson // ---- PCMG 541ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 542ccaff030SJeremy L Thompson 543ccaff030SJeremy L Thompson // ------ PCMG levels 544ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 545e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 546ccaff030SJeremy L Thompson // -------- Smoother 547ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 548ccaff030SJeremy L Thompson PC pcSmoother; 549ccaff030SJeremy L Thompson 550ccaff030SJeremy L Thompson // ---------- Smoother KSP 551ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 552ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 553ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 554ccaff030SJeremy L Thompson 555ccaff030SJeremy L Thompson // ---------- Chebyshev options 556ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 557ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 558ccaff030SJeremy L Thompson CHKERRQ(ierr); 559ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 560ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 561ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 562ccaff030SJeremy L Thompson CHKERRQ(ierr); 563ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 564ccaff030SJeremy L Thompson CHKERRQ(ierr); 565ccaff030SJeremy L Thompson 566ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 567ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 568ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 569ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 570ccaff030SJeremy L Thompson 571ccaff030SJeremy L Thompson // -------- Work vector 572ccaff030SJeremy L Thompson if (level != fineLevel) { 573ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 574ccaff030SJeremy L Thompson } 575ccaff030SJeremy L Thompson 576ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 577ccaff030SJeremy L Thompson if (level > 0) { 578ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 579ccaff030SJeremy L Thompson CHKERRQ(ierr); 580ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 581ccaff030SJeremy L Thompson CHKERRQ(ierr); 582ccaff030SJeremy L Thompson } 583ccaff030SJeremy L Thompson } 584ccaff030SJeremy L Thompson 585ccaff030SJeremy L Thompson // ------ PCMG coarse solve 586ccaff030SJeremy L Thompson KSP kspCoarse; 587ccaff030SJeremy L Thompson PC pcCoarse; 588ccaff030SJeremy L Thompson 589ccaff030SJeremy L Thompson // -------- Coarse KSP 590ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 591ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 592ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 593ccaff030SJeremy L Thompson CHKERRQ(ierr); 594ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 595ccaff030SJeremy L Thompson 596ccaff030SJeremy L Thompson // -------- Coarse preconditioner 597ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 598ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 599ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 600ccaff030SJeremy L Thompson 601ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 602ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson 604ccaff030SJeremy L Thompson // ------ PCMG options 605ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 606ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 607ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 608ccaff030SJeremy L Thompson } 609ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 610ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 611ccaff030SJeremy L Thompson } 612038d0cd7Sjeremylt { 613038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 614481a4840SJed Brown SNESLineSearch linesearch; 615481a4840SJed Brown 616481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 617481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 618481a4840SJed Brown } 619481a4840SJed Brown 620ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 621ccaff030SJeremy L Thompson 622ccaff030SJeremy L Thompson // Performance logging 623ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 624ccaff030SJeremy L Thompson 625ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 626ccaff030SJeremy L Thompson // Set initial guess 627ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 628f81c27eaSJed Brown ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 629ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 630ccaff030SJeremy L Thompson 631ccaff030SJeremy L Thompson // View solution 632ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 633ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 634ccaff030SJeremy L Thompson } 635ccaff030SJeremy L Thompson 636ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 637ccaff030SJeremy L Thompson // Solve SNES 638ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 6395f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 6405f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 6415f24f028Sjeremylt CHKERRQ(ierr); 6425f24f028Sjeremylt 643ccaff030SJeremy L Thompson // Performance logging 644ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 645ccaff030SJeremy L Thompson CHKERRQ(ierr); 646ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 647ccaff030SJeremy L Thompson 648ccaff030SJeremy L Thompson // Timing 649ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 650ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 651ccaff030SJeremy L Thompson 652ccaff030SJeremy L Thompson // Solve for each load increment 6535f24f028Sjeremylt PetscInt increment; 6545f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 6555f24f028Sjeremylt // -- Log increment count 6565f24f028Sjeremylt if (snesMonitor) { 6575f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 6585f24f028Sjeremylt CHKERRQ(ierr); 6595f24f028Sjeremylt } 6605f24f028Sjeremylt 661ccaff030SJeremy L Thompson // -- Scale the problem 662ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 663ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 664ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 665ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 666ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 667ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 668ccaff030SJeremy L Thompson } 669ccaff030SJeremy L Thompson 670ccaff030SJeremy L Thompson // -- Solve 671ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 672ccaff030SJeremy L Thompson 673ccaff030SJeremy L Thompson // -- View solution 674560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 675ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 676ccaff030SJeremy L Thompson } 677ccaff030SJeremy L Thompson 678ccaff030SJeremy L Thompson // -- Update SNES iteration count 679ccaff030SJeremy L Thompson PetscInt its; 680ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 681ccaff030SJeremy L Thompson snesIts += its; 6827418ca88SJeremy L Thompson ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr); 6837418ca88SJeremy L Thompson kspIts += its; 6843951731eSjeremylt 6853951731eSjeremylt // -- Check for divergence 6863951731eSjeremylt SNESConvergedReason reason; 6873951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6883951731eSjeremylt if (reason < 0) 6893951731eSjeremylt break; 690ccaff030SJeremy L Thompson } 691ccaff030SJeremy L Thompson 692ccaff030SJeremy L Thompson // Timing 693ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 694ccaff030SJeremy L Thompson 695ccaff030SJeremy L Thompson // Performance logging 696ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 697ccaff030SJeremy L Thompson 698ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 699ccaff030SJeremy L Thompson // Output summary 700ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 701ccaff030SJeremy L Thompson if (!appCtx->testMode) { 702ccaff030SJeremy L Thompson // -- SNES 703ccaff030SJeremy L Thompson SNESType snesType; 704ccaff030SJeremy L Thompson SNESConvergedReason reason; 705ccaff030SJeremy L Thompson PetscReal rnorm; 706ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 707ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 708ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 709ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 710ccaff030SJeremy L Thompson " SNES:\n" 711ccaff030SJeremy L Thompson " SNES Type : %s\n" 712ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 713ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 7145f24f028Sjeremylt " Completed Load Increments : %d\n" 715ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 716ccaff030SJeremy L Thompson " Final rnorm : %e\n", 717ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 7185f24f028Sjeremylt appCtx->numIncrements, increment - 1, 7195f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 720ccaff030SJeremy L Thompson 721ccaff030SJeremy L Thompson // -- KSP 722ccaff030SJeremy L Thompson KSP ksp; 723ccaff030SJeremy L Thompson KSPType kspType; 724ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 725ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 726ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 727ccaff030SJeremy L Thompson " Linear Solver:\n" 7287418ca88SJeremy L Thompson " KSP Type : %s\n" 7297418ca88SJeremy L Thompson " Total KSP Iterations : %D\n", 7307418ca88SJeremy L Thompson kspType, kspIts); CHKERRQ(ierr); 731ccaff030SJeremy L Thompson 732ccaff030SJeremy L Thompson // -- PC 733ccaff030SJeremy L Thompson PC pc; 734e3e3df41Sjeremylt PCType pcType; 735ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 736e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 737e3e3df41Sjeremylt ierr = PetscPrintf(comm, 738e3e3df41Sjeremylt " PC Type : %s\n", 739e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 740e3e3df41Sjeremylt 7412b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 742e3e3df41Sjeremylt PCMGType pcmgType; 743ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 744ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 745ccaff030SJeremy L Thompson " P-Multigrid:\n" 746ccaff030SJeremy L Thompson " PCMG Type : %s\n" 747ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 748ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 749ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 750ccaff030SJeremy L Thompson 751ccaff030SJeremy L Thompson // -- Coarse Solve 752ccaff030SJeremy L Thompson KSP kspCoarse; 753ccaff030SJeremy L Thompson PC pcCoarse; 754ccaff030SJeremy L Thompson PCType pcType; 755ccaff030SJeremy L Thompson 756ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 757ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 758ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 759ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 760ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 761ccaff030SJeremy L Thompson " Coarse Solve:\n" 762ccaff030SJeremy L Thompson " KSP Type : %s\n" 763ccaff030SJeremy L Thompson " PC Type : %s\n", 764ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 765ccaff030SJeremy L Thompson } 766ccaff030SJeremy L Thompson } 767ccaff030SJeremy L Thompson 768ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 769ccaff030SJeremy L Thompson // Compute solve time 770ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 771ccaff030SJeremy L Thompson if (!appCtx->testMode) { 772ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 773ccaff030SJeremy L Thompson CHKERRQ(ierr); 774ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 775ccaff030SJeremy L Thompson CHKERRQ(ierr); 776ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 777ccaff030SJeremy L Thompson " Performance:\n" 7787418ca88SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n" 7797418ca88SJeremy L Thompson " DoFs/Sec in SNES : %g (%g) million\n", 7807418ca88SJeremy L Thompson maxTime, minTime, 1e-6*Ugsz[fineLevel]*kspIts/maxTime, 7817418ca88SJeremy L Thompson 1e-6*Ugsz[fineLevel]*kspIts/minTime); CHKERRQ(ierr); 782ccaff030SJeremy L Thompson } 783ccaff030SJeremy L Thompson 784ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 785ccaff030SJeremy L Thompson // Compute error 786ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 787ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 788ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 789ccaff030SJeremy L Thompson const CeedScalar *truearray; 790ccaff030SJeremy L Thompson Vec errorVec, trueVec; 791ccaff030SJeremy L Thompson 792ccaff030SJeremy L Thompson // -- Work vectors 793ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 794ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 795ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 796ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 797ccaff030SJeremy L Thompson 798ccaff030SJeremy L Thompson // -- Assemble global true solution vector 79962e9c006SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, 80062e9c006SJeremy L Thompson appCtx->memTypeRequested, &truearray); 80162e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 80262e9c006SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 80362e9c006SJeremy L Thompson CHKERRQ(ierr); 80462e9c006SJeremy L Thompson } else { 80562e9c006SJeremy L Thompson ierr = VecCUDAPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 80662e9c006SJeremy L Thompson CHKERRQ(ierr); 80762e9c006SJeremy L Thompson } 808ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 809ccaff030SJeremy L Thompson CHKERRQ(ierr); 81062e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 811ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 81262e9c006SJeremy L Thompson } else { 81362e9c006SJeremy L Thompson ierr = VecCUDAResetArray(resCtx->Yloc); CHKERRQ(ierr); 81462e9c006SJeremy L Thompson } 815ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 816ccaff030SJeremy L Thompson 817ccaff030SJeremy L Thompson // -- Compute L2 error 818ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 819ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 821ccaff030SJeremy L Thompson l2Error /= l2Unorm; 822ccaff030SJeremy L Thompson 823ccaff030SJeremy L Thompson // -- Output 824ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 825ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 826ccaff030SJeremy L Thompson " L2 Error : %e\n", 827ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 828ccaff030SJeremy L Thompson } 829ccaff030SJeremy L Thompson 830ccaff030SJeremy L Thompson // -- Cleanup 831ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 832ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 8332d93065eSjeremylt } 8342d93065eSjeremylt 8352d93065eSjeremylt // --------------------------------------------------------------------------- 8362d93065eSjeremylt // Compute energy 8372d93065eSjeremylt // --------------------------------------------------------------------------- 8382d93065eSjeremylt if (!appCtx->testMode) { 8392d93065eSjeremylt // -- Compute L2 error 8402d93065eSjeremylt CeedScalar energy; 841a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 842a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 8432d93065eSjeremylt 8442d93065eSjeremylt // -- Output 8452d93065eSjeremylt ierr = PetscPrintf(comm, 8462d93065eSjeremylt " Strain Energy : %e\n", 8472d93065eSjeremylt energy); CHKERRQ(ierr); 848ccaff030SJeremy L Thompson } 849ccaff030SJeremy L Thompson 850ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 8515c25879aSJeremy L Thompson // Output diagnostic quantities 8525c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 8535c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 8545c25879aSJeremy L Thompson // -- Setup context 8555c25879aSJeremy L Thompson UserMult diagnosticCtx; 8565c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 8575c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 8585c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 8595c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 8605c25879aSJeremy L Thompson 8615c25879aSJeremy L Thompson // -- Compute and output 8625c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 8635c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 8645c25879aSJeremy L Thompson CHKERRQ(ierr); 8655c25879aSJeremy L Thompson 8665c25879aSJeremy L Thompson // -- Cleanup 8675c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 8685c25879aSJeremy L Thompson } 8695c25879aSJeremy L Thompson 8705c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 871ccaff030SJeremy L Thompson // Free objects 872ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 873ccaff030SJeremy L Thompson // Data in arrays per level 874e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 875ccaff030SJeremy L Thompson // Vectors 876ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 877ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 878ccaff030SJeremy L Thompson 879ccaff030SJeremy L Thompson // Jacobian matrix and data 880ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 881ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 882ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 883ccaff030SJeremy L Thompson 884ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 885ccaff030SJeremy L Thompson if (level > 0) { 886ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 887ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 888ccaff030SJeremy L Thompson } 889ccaff030SJeremy L Thompson 890ccaff030SJeremy L Thompson // DM 891ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 892ccaff030SJeremy L Thompson 893ccaff030SJeremy L Thompson // libCEED objects 894ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 895ccaff030SJeremy L Thompson } 896ccaff030SJeremy L Thompson 897ccaff030SJeremy L Thompson // Arrays 898ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 899ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 900ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 901ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 902ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 903ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 904ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 905ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 906ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 907ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 908ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 909ccaff030SJeremy L Thompson 910ccaff030SJeremy L Thompson // libCEED objects 911*777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhys); 912*777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhysSmoother); 913ccaff030SJeremy L Thompson CeedDestroy(&ceed); 914ccaff030SJeremy L Thompson 915ccaff030SJeremy L Thompson // PETSc objects 916ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 917ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 918ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 919ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 920ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 921ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 922ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 923ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 924ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 925a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 926a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 927ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 928ccaff030SJeremy L Thompson 929ccaff030SJeremy L Thompson // Structs 930ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 931ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 932ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 933ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 934ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 935f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 936ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 937ccaff030SJeremy L Thompson 938ccaff030SJeremy L Thompson return PetscFinalize(); 939ccaff030SJeremy L Thompson } 940