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 3228688798Sjeremylt // ./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/cuda 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 61fe394131Sjeremylt Vec NBCs = NULL, NBCsloc = NULL; 62e3e3df41Sjeremylt SNES snes, snesCoarse = NULL; 63ccaff030SJeremy L Thompson Mat *jacobMat, jacobMatCoarse, *prolongRestrMat; 64ccaff030SJeremy L Thompson // PETSc data 65e3e3df41Sjeremylt UserMult resCtx, jacobCoarseCtx = NULL, *jacobCtx; 66ccaff030SJeremy L Thompson FormJacobCtx formJacobCtx; 67ccaff030SJeremy L Thompson UserMultProlongRestr *prolongRestrCtx; 68ccaff030SJeremy L Thompson PCMGCycleType pcmgCycleType = PC_MG_CYCLE_V; 69ccaff030SJeremy L Thompson // libCEED objects 70d99fa3c5SJeremy L Thompson Ceed ceed; 71ccaff030SJeremy L Thompson CeedData *ceedData; 72777ff853SJeremy L Thompson CeedQFunctionContext ctxPhys, ctxPhysSmoother = NULL; 73ccaff030SJeremy L Thompson // Parameters 74ccaff030SJeremy L Thompson PetscInt ncompu = 3; // 3 DoFs in 3D 7513fdad4bSJeremy L Thompson PetscInt ncompe = 1, ncompd = 5; // 1 energy output, 5 diagnostic 76ccaff030SJeremy L Thompson PetscInt numLevels = 1, fineLevel = 0; 77ccaff030SJeremy L Thompson PetscInt *Ugsz, *Ulsz, *Ulocsz; // sz: size 787418ca88SJeremy L Thompson PetscInt snesIts = 0, kspIts = 0; 79ccaff030SJeremy L Thompson // Timing 80ccaff030SJeremy L Thompson double startTime, elapsedTime, minTime, maxTime; 81ccaff030SJeremy L Thompson 82ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 83ccaff030SJeremy L Thompson if (ierr) 84ccaff030SJeremy L Thompson return ierr; 85ccaff030SJeremy L Thompson 86ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 87ccaff030SJeremy L Thompson // Process command line options 88ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 89ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 90ccaff030SJeremy L Thompson 91ccaff030SJeremy L Thompson // -- Set mesh file, polynomial degree, problem type 92ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr); 93ccaff030SJeremy L Thompson ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr); 94ccaff030SJeremy L Thompson numLevels = appCtx->numLevels; 95ccaff030SJeremy L Thompson fineLevel = numLevels - 1; 96ccaff030SJeremy L Thompson 97ccaff030SJeremy L Thompson // -- Set Poison's ratio, Young's Modulus 98ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr); 99ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 100ccaff030SJeremy L Thompson ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr); 101f7b4142eSJeremy L Thompson if (fabs(appCtx->nuSmoother) > 1E-14) { 102f7b4142eSJeremy L Thompson ierr = PetscMalloc1(1, &physSmoother); CHKERRQ(ierr); 103f7b4142eSJeremy L Thompson ierr = PetscMemcpy(physSmoother, phys, sizeof(*phys)); CHKERRQ(ierr); 104f7b4142eSJeremy L Thompson physSmoother->nu = appCtx->nuSmoother; 105f7b4142eSJeremy L Thompson } 106ccaff030SJeremy L Thompson 107ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 10862e9c006SJeremy L Thompson // Initalize libCEED 10962e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 11062e9c006SJeremy L Thompson // Initalize backend 11162e9c006SJeremy L Thompson CeedInit(appCtx->ceedResource, &ceed); 11262e9c006SJeremy L Thompson 11362e9c006SJeremy L Thompson // Check preferred MemType 11462e9c006SJeremy L Thompson CeedMemType memTypeBackend; 11562e9c006SJeremy L Thompson CeedGetPreferredMemType(ceed, &memTypeBackend); 11662e9c006SJeremy L Thompson 117777ff853SJeremy L Thompson // Wrap context in libCEED objects 118777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhys); 119777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhys, CEED_MEM_HOST, CEED_USE_POINTER, 120777ff853SJeremy L Thompson sizeof(*phys), phys); 121777ff853SJeremy L Thompson if (physSmoother) { 122777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhysSmoother); 123777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhysSmoother, CEED_MEM_HOST, CEED_USE_POINTER, 124777ff853SJeremy L Thompson sizeof(*physSmoother), physSmoother); 125777ff853SJeremy L Thompson } 126777ff853SJeremy L Thompson 12762e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 128ccaff030SJeremy L Thompson // Setup DM 129ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 130ccaff030SJeremy L Thompson // Performance logging 131ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 132ccaff030SJeremy L Thompson CHKERRQ(ierr); 133ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 134ccaff030SJeremy L Thompson 135ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 136ccaff030SJeremy L Thompson ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 137*b68a8d79SJed Brown VecType vectype; 138*b68a8d79SJed Brown switch (memTypeBackend) { 139*b68a8d79SJed Brown case CEED_MEM_HOST: vectype = VECSTANDARD; break; 140*b68a8d79SJed Brown case CEED_MEM_DEVICE: { 141*b68a8d79SJed Brown const char *resolved; 142*b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 143*b68a8d79SJed Brown if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 144*b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 145*b68a8d79SJed Brown else vectype = VECSTANDARD; 146*b68a8d79SJed Brown } 147*b68a8d79SJed Brown } 148*b68a8d79SJed Brown ierr = DMSetVecType(dmOrig, vectype); CHKERRQ(ierr); 149*b68a8d79SJed Brown ierr = DMSetFromOptions(dmOrig); CHKERRQ(ierr); 150ccaff030SJeremy L Thompson 151ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 152ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 153e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 154ccaff030SJeremy L Thompson ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 155*b68a8d79SJed Brown ierr = DMGetVecType(dmOrig, &vectype); CHKERRQ(ierr); 156*b68a8d79SJed Brown ierr = DMSetVecType(levelDMs[level], vectype); CHKERRQ(ierr); 157ccaff030SJeremy L Thompson ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 158a3c02c40SJeremy L Thompson PETSC_TRUE, ncompu); CHKERRQ(ierr); 159a3c02c40SJeremy L Thompson // -- Label field components for viewing 160a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 161a3c02c40SJeremy L Thompson PetscSection section; 162a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(levelDMs[level], §ion); CHKERRQ(ierr); 163a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 164a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 165a3c02c40SJeremy L Thompson CHKERRQ(ierr); 166a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 167a3c02c40SJeremy L Thompson CHKERRQ(ierr); 168a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 169a3c02c40SJeremy L Thompson CHKERRQ(ierr); 170a3c02c40SJeremy L Thompson } 171a3c02c40SJeremy L Thompson 172a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 173a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr); 174a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel], 1755c25879aSJeremy L Thompson PETSC_FALSE, ncompe); CHKERRQ(ierr); 176a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr); 177a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel], 1788ca58ff3SJeremy L Thompson PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr); 179*b68a8d79SJed Brown ierr = DMSetVecType(dmEnergy, vectype); CHKERRQ(ierr); 180*b68a8d79SJed Brown ierr = DMSetVecType(dmDiagnostic, vectype); CHKERRQ(ierr); 181a3c02c40SJeremy L Thompson { 182a3c02c40SJeremy L Thompson // -- Label field components for viewing 183a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 184a3c02c40SJeremy L Thompson PetscSection section; 185a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(dmDiagnostic, §ion); CHKERRQ(ierr); 186a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 1878ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 1888ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1898ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 1908ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1918ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 1928ca58ff3SJeremy L Thompson CHKERRQ(ierr); 193379387d4SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure"); 1948ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1957ab18febSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"); 19613fdad4bSJeremy L Thompson CHKERRQ(ierr); 19713fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 19813fdad4bSJeremy L Thompson CHKERRQ(ierr); 19913fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 20013fdad4bSJeremy L Thompson CHKERRQ(ierr); 20113fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 202a3c02c40SJeremy L Thompson CHKERRQ(ierr); 203ccaff030SJeremy L Thompson } 204ccaff030SJeremy L Thompson 205ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 206ccaff030SJeremy L Thompson // Setup solution and work vectors 207ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 208ccaff030SJeremy L Thompson // Allocate arrays 209ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 210ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 211ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 212ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 213ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 214ccaff030SJeremy L Thompson 215ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 216e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 217ccaff030SJeremy L Thompson // -- Create global unknown vector U 218ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 219ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 220ccaff030SJeremy L Thompson // Note: Local size for matShell 221ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 222ccaff030SJeremy L Thompson 223ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 224ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 225ccaff030SJeremy L Thompson // Note: local size for libCEED 226ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 227ccaff030SJeremy L Thompson } 228ccaff030SJeremy L Thompson 229ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 230ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 231ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 232ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 233ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 234ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 235ccaff030SJeremy L Thompson 236ccaff030SJeremy L Thompson // Performance logging 237ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 238ccaff030SJeremy L Thompson 239ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 240ccaff030SJeremy L Thompson // Set up libCEED 241ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 242ccaff030SJeremy L Thompson // Performance logging 243ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 244ccaff030SJeremy L Thompson CHKERRQ(ierr); 245ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 246ccaff030SJeremy L Thompson 247ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 248ccaff030SJeremy L Thompson CeedVector forceCeed; 249ccaff030SJeremy L Thompson CeedScalar *f; 250*b68a8d79SJed Brown PetscMemType fmemtype; 251ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 252*b68a8d79SJed Brown ierr = VecGetArrayAndMemType(Floc, &f, &fmemtype); CHKERRQ(ierr); 253ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 254*b68a8d79SJed Brown CeedVectorSetArray(forceCeed, MemTypeP2C(fmemtype), CEED_USE_POINTER, f); 255ccaff030SJeremy L Thompson } 256ccaff030SJeremy L Thompson 257fe394131Sjeremylt // -- Create libCEED local Neumann BCs vector 258fe394131Sjeremylt CeedVector neumannCeed; 259fe394131Sjeremylt CeedScalar *n; 260*b68a8d79SJed Brown PetscMemType nmemtype; 261fe394131Sjeremylt if (appCtx->bcTractionCount > 0) { 262fe394131Sjeremylt ierr = VecDuplicate(U, &NBCs); CHKERRQ(ierr); 263fe394131Sjeremylt ierr = VecDuplicate(Uloc[fineLevel], &NBCsloc); CHKERRQ(ierr); 264*b68a8d79SJed Brown ierr = VecGetArrayAndMemType(NBCsloc, &n, &nmemtype); CHKERRQ(ierr); 265fe394131Sjeremylt CeedVectorCreate(ceed, Ulocsz[fineLevel], &neumannCeed); 266*b68a8d79SJed Brown CeedVectorSetArray(neumannCeed, MemTypeP2C(nmemtype), 267fe394131Sjeremylt CEED_USE_POINTER, n); 268fe394131Sjeremylt } 269fe394131Sjeremylt 270ccaff030SJeremy L Thompson // -- Setup libCEED objects 271ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 272d99fa3c5SJeremy L Thompson // ---- Setup residual, Jacobian evaluator and geometric information 273ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 274ccaff030SJeremy L Thompson { 275a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 276777ff853SJeremy L Thompson ceed, appCtx, ctxPhys, ceedData, fineLevel, 277a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 278fe394131Sjeremylt forceCeed, neumannCeed); 279ccaff030SJeremy L Thompson CHKERRQ(ierr); 280ccaff030SJeremy L Thompson } 281d99fa3c5SJeremy L Thompson // ---- Setup coarse Jacobian evaluator and prolongation/restriction 282d99fa3c5SJeremy L Thompson for (PetscInt level = numLevels - 2; level >= 0; level--) { 283ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 284d99fa3c5SJeremy L Thompson 285d99fa3c5SJeremy L Thompson // Get global communication restriction 286d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 287d99fa3c5SJeremy L Thompson ierr = VecSet(Uloc[level+1], 1.0); CHKERRQ(ierr); 288d99fa3c5SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[level+1], Uloc[level+1], ADD_VALUES, 289d99fa3c5SJeremy L Thompson Ug[level+1]); CHKERRQ(ierr); 290d99fa3c5SJeremy L Thompson ierr = DMGlobalToLocal(levelDMs[level+1], Ug[level+1], INSERT_VALUES, 291d99fa3c5SJeremy L Thompson Uloc[level+1]); CHKERRQ(ierr); 292d99fa3c5SJeremy L Thompson 293d99fa3c5SJeremy L Thompson // Place in libCEED array 294d99fa3c5SJeremy L Thompson const PetscScalar *m; 295*b68a8d79SJed Brown PetscMemType mmemtype; 296*b68a8d79SJed Brown ierr = VecGetArrayReadAndMemType(Uloc[level+1], &m, &mmemtype); CHKERRQ(ierr); 297*b68a8d79SJed Brown CeedVectorSetArray(ceedData[level+1]->xceed, MemTypeP2C(mmemtype), 298d99fa3c5SJeremy L Thompson CEED_USE_POINTER, (CeedScalar *)m); 299ccaff030SJeremy L Thompson 3001c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 301777ff853SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], ceed, appCtx, 302ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 303777ff853SJeremy L Thompson Ulocsz[level], ceedData[level+1]->xceed); 304777ff853SJeremy L Thompson CHKERRQ(ierr); 305d99fa3c5SJeremy L Thompson 306d99fa3c5SJeremy L Thompson // Restore PETSc vector 307*b68a8d79SJed Brown CeedVectorTakeArray(ceedData[level+1]->xceed, MemTypeP2C(mmemtype), 308d99fa3c5SJeremy L Thompson (CeedScalar **)&m); 309*b68a8d79SJed Brown ierr = VecRestoreArrayReadAndMemType(Uloc[level+1], &m); CHKERRQ(ierr); 310d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 311d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Uloc[level+1]); CHKERRQ(ierr); 312ccaff030SJeremy L Thompson } 313ccaff030SJeremy L Thompson 314ccaff030SJeremy L Thompson // Performance logging 315ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 316ccaff030SJeremy L Thompson 317ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 318fe394131Sjeremylt // Setup global forcing and Neumann BC vectors 319ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 320ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 321ccaff030SJeremy L Thompson 322ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 323*b68a8d79SJed Brown CeedVectorTakeArray(forceCeed, MemTypeP2C(fmemtype), NULL); 324*b68a8d79SJed Brown ierr = VecRestoreArrayAndMemType(Floc, &f); CHKERRQ(ierr); 325ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 326ccaff030SJeremy L Thompson CHKERRQ(ierr); 327ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 328ccaff030SJeremy L Thompson } 329ccaff030SJeremy L Thompson 330fe394131Sjeremylt if (appCtx->bcTractionCount > 0) { 331fe394131Sjeremylt ierr = VecZeroEntries(NBCs); CHKERRQ(ierr); 332*b68a8d79SJed Brown CeedVectorTakeArray(neumannCeed, MemTypeP2C(nmemtype), NULL); 333*b68a8d79SJed Brown ierr = VecRestoreArrayAndMemType(NBCsloc, &n); CHKERRQ(ierr); 334fe394131Sjeremylt ierr = DMLocalToGlobal(levelDMs[fineLevel], NBCsloc, ADD_VALUES, NBCs); 335fe394131Sjeremylt CHKERRQ(ierr); 336fe394131Sjeremylt CeedVectorDestroy(&neumannCeed); 337fe394131Sjeremylt } 338fe394131Sjeremylt 339ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 340ccaff030SJeremy L Thompson // Print problem summary 341ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 342ccaff030SJeremy L Thompson if (!appCtx->testMode) { 343ccaff030SJeremy L Thompson const char *usedresource; 344ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 345ccaff030SJeremy L Thompson 346ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 34717fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 348ccaff030SJeremy L Thompson " libCEED:\n" 34962e9c006SJeremy L Thompson " libCEED Backend : %s\n" 350*b68a8d79SJed Brown " libCEED Backend MemType : %s\n", 351*b68a8d79SJed Brown usedresource, CeedMemTypes[memTypeBackend]); 35262e9c006SJeremy L Thompson CHKERRQ(ierr); 353ccaff030SJeremy L Thompson 35462e9c006SJeremy L Thompson VecType vecType; 35562e9c006SJeremy L Thompson ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 35662e9c006SJeremy L Thompson ierr = PetscPrintf(comm, 35762e9c006SJeremy L Thompson " PETSc:\n" 35862e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 35962e9c006SJeremy L Thompson vecType); CHKERRQ(ierr); 36062e9c006SJeremy L Thompson 361ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 362ccaff030SJeremy L Thompson " Problem:\n" 363ccaff030SJeremy L Thompson " Problem Name : %s\n" 364ccaff030SJeremy L Thompson " Forcing Function : %s\n" 365ccaff030SJeremy L Thompson " Mesh:\n" 366ccaff030SJeremy L Thompson " File : %s\n" 367ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 368ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 369ccaff030SJeremy L Thompson " Global nodes : %D\n" 370ccaff030SJeremy L Thompson " Owned nodes : %D\n" 371ccaff030SJeremy L Thompson " DoF per node : %D\n" 372ccaff030SJeremy L Thompson " Multigrid:\n" 373ccaff030SJeremy L Thompson " Type : %s\n" 374ccaff030SJeremy L Thompson " Number of Levels : %d\n", 375ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 376ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 377ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 378ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 379ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 380f892e6d1Sjeremylt (appCtx->degree == 1 && 381f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 382f892e6d1Sjeremylt "Algebraic multigrid" : 383ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 384f892e6d1Sjeremylt (appCtx->degree == 1 || 385f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 386f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 387ccaff030SJeremy L Thompson 388ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 389e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 390ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 391ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 392ccaff030SJeremy L Thompson " Level %D (%s):\n" 393ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 394ccaff030SJeremy L Thompson " Global Nodes : %D\n" 395ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 396ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 397ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 398ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 399ccaff030SJeremy L Thompson CHKERRQ(ierr); 400ccaff030SJeremy L Thompson } 401ccaff030SJeremy L Thompson } 402ccaff030SJeremy L Thompson } 403ccaff030SJeremy L Thompson 404ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 405ccaff030SJeremy L Thompson // Setup SNES 406ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 407ccaff030SJeremy L Thompson // Performance logging 408ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 409ccaff030SJeremy L Thompson CHKERRQ(ierr); 410ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 411ccaff030SJeremy L Thompson 412ccaff030SJeremy L Thompson // Create SNES 413ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 414ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 415ccaff030SJeremy L Thompson 416ccaff030SJeremy L Thompson // -- Jacobian evaluators 417ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 418ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 419e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 420ccaff030SJeremy L Thompson // -- Jacobian context for level 421ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 422ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 423777ff853SJeremy L Thompson Uloc[level], ceedData[level], ceed, ctxPhys, 424777ff853SJeremy L Thompson ctxPhysSmoother, jacobCtx[level]); CHKERRQ(ierr); 425ccaff030SJeremy L Thompson 426ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 427ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 428ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 429ccaff030SJeremy L Thompson CHKERRQ(ierr); 430ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 431ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 432ccaff030SJeremy L Thompson CHKERRQ(ierr); 433ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 434ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 435*b68a8d79SJed Brown ierr = MatShellSetVecType(jacobMat[level], vectype); CHKERRQ(ierr); 436ccaff030SJeremy L Thompson } 437e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 438ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 439ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 440ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 441ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 442ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 443ccaff030SJeremy L Thompson 444ccaff030SJeremy L Thompson // -- Residual evaluation function 445fe394131Sjeremylt ierr = PetscCalloc1(1, &resCtx); CHKERRQ(ierr); 446ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 447ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 448ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 449f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 450fe394131Sjeremylt if (appCtx->bcTractionCount > 0) 451fe394131Sjeremylt resCtx->NBCs = NBCs; 452fe394131Sjeremylt else 453fe394131Sjeremylt resCtx->NBCs = NULL; 454ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 455ccaff030SJeremy L Thompson 456ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 457ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 458ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 459e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 460ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 461ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 46262e9c006SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, appCtx, levelDMs[level-1], 46362e9c006SJeremy L Thompson levelDMs[level], Ug[level], Uloc[level-1], 46462e9c006SJeremy L Thompson Uloc[level], ceedData[level-1], 46562e9c006SJeremy L Thompson ceedData[level], ceed, 466ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 467ccaff030SJeremy L Thompson 468ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 469ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 470ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 471ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 472ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 473ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 474ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 475ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 476ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 477ccaff030SJeremy L Thompson CHKERRQ(ierr); 478*b68a8d79SJed Brown ierr = MatShellSetVecType(prolongRestrMat[level], vectype); CHKERRQ(ierr); 479ccaff030SJeremy L Thompson } 480ccaff030SJeremy L Thompson 481ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 482ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 483ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 484e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 485e3e3df41Sjeremylt // -- Jacobian Matrix 486e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 487e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 488e3e3df41Sjeremylt 489e3e3df41Sjeremylt if (appCtx->degree > 1) { 490ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 491ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 492ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 493ccaff030SJeremy L Thompson 494e3e3df41Sjeremylt // -- Jacobian function 495ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 496ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 497ccaff030SJeremy L Thompson 498ccaff030SJeremy L Thompson // -- Residual evaluation function 499ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 500ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 501ccaff030SJeremy L Thompson CHKERRQ(ierr); 502ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 503ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 504ccaff030SJeremy L Thompson 505ccaff030SJeremy L Thompson // -- Update formJacobCtx 506ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 507ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 508ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 509e3e3df41Sjeremylt } 510e3e3df41Sjeremylt } 511e3e3df41Sjeremylt 512e3e3df41Sjeremylt // Set Jacobian function 513e3e3df41Sjeremylt if (appCtx->degree > 1) { 514e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 515e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 516e3e3df41Sjeremylt } else { 517e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 518e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 519e3e3df41Sjeremylt CHKERRQ(ierr); 520e3e3df41Sjeremylt } 521ccaff030SJeremy L Thompson 522ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 523ccaff030SJeremy L Thompson // Setup KSP 524ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 525ccaff030SJeremy L Thompson { 526ccaff030SJeremy L Thompson PC pc; 527ccaff030SJeremy L Thompson KSP ksp; 528ccaff030SJeremy L Thompson 529ccaff030SJeremy L Thompson // -- KSP 530ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 531ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 532ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 533ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 534ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 535ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 536ccaff030SJeremy L Thompson 537ccaff030SJeremy L Thompson // -- Preconditioning 538ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 539ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 540ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 541ccaff030SJeremy L Thompson 542ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 543ccaff030SJeremy L Thompson // ---- No Multigrid 544ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 545ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 546f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 547f892e6d1Sjeremylt // ---- AMG for degree 1 548f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 549ccaff030SJeremy L Thompson } else { 550ccaff030SJeremy L Thompson // ---- PCMG 551ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 552ccaff030SJeremy L Thompson 553ccaff030SJeremy L Thompson // ------ PCMG levels 554ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 555e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 556ccaff030SJeremy L Thompson // -------- Smoother 557ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 558ccaff030SJeremy L Thompson PC pcSmoother; 559ccaff030SJeremy L Thompson 560ccaff030SJeremy L Thompson // ---------- Smoother KSP 561ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 562ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 563ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 564ccaff030SJeremy L Thompson 565ccaff030SJeremy L Thompson // ---------- Chebyshev options 566ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 567ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 568ccaff030SJeremy L Thompson CHKERRQ(ierr); 569ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 570ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 571ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 572ccaff030SJeremy L Thompson CHKERRQ(ierr); 573ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 574ccaff030SJeremy L Thompson CHKERRQ(ierr); 575ccaff030SJeremy L Thompson 576ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 577ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 578ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 579ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 580ccaff030SJeremy L Thompson 581ccaff030SJeremy L Thompson // -------- Work vector 582ccaff030SJeremy L Thompson if (level != fineLevel) { 583ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 584ccaff030SJeremy L Thompson } 585ccaff030SJeremy L Thompson 586ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 587ccaff030SJeremy L Thompson if (level > 0) { 588ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 589ccaff030SJeremy L Thompson CHKERRQ(ierr); 590ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 591ccaff030SJeremy L Thompson CHKERRQ(ierr); 592ccaff030SJeremy L Thompson } 593ccaff030SJeremy L Thompson } 594ccaff030SJeremy L Thompson 595ccaff030SJeremy L Thompson // ------ PCMG coarse solve 596ccaff030SJeremy L Thompson KSP kspCoarse; 597ccaff030SJeremy L Thompson PC pcCoarse; 598ccaff030SJeremy L Thompson 599ccaff030SJeremy L Thompson // -------- Coarse KSP 600ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 601ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 602ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 603ccaff030SJeremy L Thompson CHKERRQ(ierr); 604ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 605ccaff030SJeremy L Thompson 606ccaff030SJeremy L Thompson // -------- Coarse preconditioner 607ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 608ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 609ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 610ccaff030SJeremy L Thompson 611ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 612ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 613ccaff030SJeremy L Thompson 614ccaff030SJeremy L Thompson // ------ PCMG options 615ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 616ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 617ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 618ccaff030SJeremy L Thompson } 619ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 620ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 621ccaff030SJeremy L Thompson } 622038d0cd7Sjeremylt { 623038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 624481a4840SJed Brown SNESLineSearch linesearch; 625481a4840SJed Brown 626481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 627481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 628481a4840SJed Brown } 629481a4840SJed Brown 630ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 631ccaff030SJeremy L Thompson 632ccaff030SJeremy L Thompson // Performance logging 633ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 634ccaff030SJeremy L Thompson 635ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 636ccaff030SJeremy L Thompson // Set initial guess 637ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 638f81c27eaSJed Brown ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 639ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 640ccaff030SJeremy L Thompson 641ccaff030SJeremy L Thompson // View solution 642ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 643ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 644ccaff030SJeremy L Thompson } 645ccaff030SJeremy L Thompson 646ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 647ccaff030SJeremy L Thompson // Solve SNES 648ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 6495f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 6505f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 6515f24f028Sjeremylt CHKERRQ(ierr); 6525f24f028Sjeremylt 653ccaff030SJeremy L Thompson // Performance logging 654ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 655ccaff030SJeremy L Thompson CHKERRQ(ierr); 656ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 657ccaff030SJeremy L Thompson 658ccaff030SJeremy L Thompson // Timing 659ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 660ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 661ccaff030SJeremy L Thompson 662ccaff030SJeremy L Thompson // Solve for each load increment 6635f24f028Sjeremylt PetscInt increment; 6645f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 6655f24f028Sjeremylt // -- Log increment count 6665f24f028Sjeremylt if (snesMonitor) { 6675f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 6685f24f028Sjeremylt CHKERRQ(ierr); 6695f24f028Sjeremylt } 6705f24f028Sjeremylt 671ccaff030SJeremy L Thompson // -- Scale the problem 672ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 673ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 674ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 675ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 676ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 677ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 678ccaff030SJeremy L Thompson } 679ccaff030SJeremy L Thompson 680ccaff030SJeremy L Thompson // -- Solve 681ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 682ccaff030SJeremy L Thompson 683ccaff030SJeremy L Thompson // -- View solution 684560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 685ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 686ccaff030SJeremy L Thompson } 687ccaff030SJeremy L Thompson 688ccaff030SJeremy L Thompson // -- Update SNES iteration count 689ccaff030SJeremy L Thompson PetscInt its; 690ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 691ccaff030SJeremy L Thompson snesIts += its; 6927418ca88SJeremy L Thompson ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr); 6937418ca88SJeremy L Thompson kspIts += its; 6943951731eSjeremylt 6953951731eSjeremylt // -- Check for divergence 6963951731eSjeremylt SNESConvergedReason reason; 6973951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6983951731eSjeremylt if (reason < 0) 6993951731eSjeremylt break; 700ccaff030SJeremy L Thompson } 701ccaff030SJeremy L Thompson 702ccaff030SJeremy L Thompson // Timing 703ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 704ccaff030SJeremy L Thompson 705ccaff030SJeremy L Thompson // Performance logging 706ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 707ccaff030SJeremy L Thompson 708ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 709ccaff030SJeremy L Thompson // Output summary 710ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 711ccaff030SJeremy L Thompson if (!appCtx->testMode) { 712ccaff030SJeremy L Thompson // -- SNES 713ccaff030SJeremy L Thompson SNESType snesType; 714ccaff030SJeremy L Thompson SNESConvergedReason reason; 715ccaff030SJeremy L Thompson PetscReal rnorm; 716ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 717ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 718ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 719ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 720ccaff030SJeremy L Thompson " SNES:\n" 721ccaff030SJeremy L Thompson " SNES Type : %s\n" 722ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 723ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 7245f24f028Sjeremylt " Completed Load Increments : %d\n" 725ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 726ccaff030SJeremy L Thompson " Final rnorm : %e\n", 727ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 7285f24f028Sjeremylt appCtx->numIncrements, increment - 1, 7295f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 730ccaff030SJeremy L Thompson 731ccaff030SJeremy L Thompson // -- KSP 732ccaff030SJeremy L Thompson KSP ksp; 733ccaff030SJeremy L Thompson KSPType kspType; 734ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 735ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 736ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 737ccaff030SJeremy L Thompson " Linear Solver:\n" 7387418ca88SJeremy L Thompson " KSP Type : %s\n" 7397418ca88SJeremy L Thompson " Total KSP Iterations : %D\n", 7407418ca88SJeremy L Thompson kspType, kspIts); CHKERRQ(ierr); 741ccaff030SJeremy L Thompson 742ccaff030SJeremy L Thompson // -- PC 743ccaff030SJeremy L Thompson PC pc; 744e3e3df41Sjeremylt PCType pcType; 745ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 746e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 747e3e3df41Sjeremylt ierr = PetscPrintf(comm, 748e3e3df41Sjeremylt " PC Type : %s\n", 749e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 750e3e3df41Sjeremylt 7512b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 752e3e3df41Sjeremylt PCMGType pcmgType; 753ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 754ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 755ccaff030SJeremy L Thompson " P-Multigrid:\n" 756ccaff030SJeremy L Thompson " PCMG Type : %s\n" 757ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 758ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 759ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 760ccaff030SJeremy L Thompson 761ccaff030SJeremy L Thompson // -- Coarse Solve 762ccaff030SJeremy L Thompson KSP kspCoarse; 763ccaff030SJeremy L Thompson PC pcCoarse; 764ccaff030SJeremy L Thompson PCType pcType; 765ccaff030SJeremy L Thompson 766ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 767ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 768ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 769ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 770ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 771ccaff030SJeremy L Thompson " Coarse Solve:\n" 772ccaff030SJeremy L Thompson " KSP Type : %s\n" 773ccaff030SJeremy L Thompson " PC Type : %s\n", 774ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 775ccaff030SJeremy L Thompson } 776ccaff030SJeremy L Thompson } 777ccaff030SJeremy L Thompson 778ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 779ccaff030SJeremy L Thompson // Compute solve time 780ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 781ccaff030SJeremy L Thompson if (!appCtx->testMode) { 782ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 783ccaff030SJeremy L Thompson CHKERRQ(ierr); 784ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 785ccaff030SJeremy L Thompson CHKERRQ(ierr); 786ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 787ccaff030SJeremy L Thompson " Performance:\n" 7887418ca88SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n" 7897418ca88SJeremy L Thompson " DoFs/Sec in SNES : %g (%g) million\n", 7907418ca88SJeremy L Thompson maxTime, minTime, 1e-6*Ugsz[fineLevel]*kspIts/maxTime, 7917418ca88SJeremy L Thompson 1e-6*Ugsz[fineLevel]*kspIts/minTime); CHKERRQ(ierr); 792ccaff030SJeremy L Thompson } 793ccaff030SJeremy L Thompson 794ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 795ccaff030SJeremy L Thompson // Compute error 796ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 797ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 798ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 799ccaff030SJeremy L Thompson const CeedScalar *truearray; 800ccaff030SJeremy L Thompson Vec errorVec, trueVec; 801ccaff030SJeremy L Thompson 802ccaff030SJeremy L Thompson // -- Work vectors 803ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 804ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 805ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 806ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 807ccaff030SJeremy L Thompson 808ccaff030SJeremy L Thompson // -- Assemble global true solution vector 80962e9c006SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, 810*b68a8d79SJed Brown CEED_MEM_HOST, &truearray); 81162e9c006SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 81262e9c006SJeremy L Thompson CHKERRQ(ierr); 813ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 814ccaff030SJeremy L Thompson CHKERRQ(ierr); 815ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 817ccaff030SJeremy L Thompson 818ccaff030SJeremy L Thompson // -- Compute L2 error 819ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 821ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 822ccaff030SJeremy L Thompson l2Error /= l2Unorm; 823ccaff030SJeremy L Thompson 824ccaff030SJeremy L Thompson // -- Output 825ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 826ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 827ccaff030SJeremy L Thompson " L2 Error : %e\n", 828ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 829ccaff030SJeremy L Thompson } 830ccaff030SJeremy L Thompson 831ccaff030SJeremy L Thompson // -- Cleanup 832ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 833ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 8342d93065eSjeremylt } 8352d93065eSjeremylt 8362d93065eSjeremylt // --------------------------------------------------------------------------- 8372d93065eSjeremylt // Compute energy 8382d93065eSjeremylt // --------------------------------------------------------------------------- 8392d93065eSjeremylt if (!appCtx->testMode) { 8402d93065eSjeremylt // -- Compute L2 error 8412d93065eSjeremylt CeedScalar energy; 842a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 843a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 8442d93065eSjeremylt 8452d93065eSjeremylt // -- Output 8462d93065eSjeremylt ierr = PetscPrintf(comm, 8472d93065eSjeremylt " Strain Energy : %e\n", 8482d93065eSjeremylt energy); CHKERRQ(ierr); 849ccaff030SJeremy L Thompson } 850ccaff030SJeremy L Thompson 851ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 8525c25879aSJeremy L Thompson // Output diagnostic quantities 8535c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 8545c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 8555c25879aSJeremy L Thompson // -- Setup context 8565c25879aSJeremy L Thompson UserMult diagnosticCtx; 8575c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 8585c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 8595c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 8605c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 8615c25879aSJeremy L Thompson 8625c25879aSJeremy L Thompson // -- Compute and output 8635c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 8645c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 8655c25879aSJeremy L Thompson CHKERRQ(ierr); 8665c25879aSJeremy L Thompson 8675c25879aSJeremy L Thompson // -- Cleanup 8685c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 8695c25879aSJeremy L Thompson } 8705c25879aSJeremy L Thompson 8715c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 872ccaff030SJeremy L Thompson // Free objects 873ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 874ccaff030SJeremy L Thompson // Data in arrays per level 875e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 876ccaff030SJeremy L Thompson // Vectors 877ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 878ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 879ccaff030SJeremy L Thompson 880ccaff030SJeremy L Thompson // Jacobian matrix and data 881ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 882ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 883ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 884ccaff030SJeremy L Thompson 885ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 886ccaff030SJeremy L Thompson if (level > 0) { 887ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 888ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 889ccaff030SJeremy L Thompson } 890ccaff030SJeremy L Thompson 891ccaff030SJeremy L Thompson // DM 892ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 893ccaff030SJeremy L Thompson 894ccaff030SJeremy L Thompson // libCEED objects 895ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 896ccaff030SJeremy L Thompson } 897ccaff030SJeremy L Thompson 898ccaff030SJeremy L Thompson // Arrays 899ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 900ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 901ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 902ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 903ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 904ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 905ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 906ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 907ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 908ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 909ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 910ccaff030SJeremy L Thompson 911ccaff030SJeremy L Thompson // libCEED objects 912777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhys); 913777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhysSmoother); 914ccaff030SJeremy L Thompson CeedDestroy(&ceed); 915ccaff030SJeremy L Thompson 916ccaff030SJeremy L Thompson // PETSc objects 917ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 918ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 919ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 920ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 921ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 922fe394131Sjeremylt ierr = VecDestroy(&NBCs); CHKERRQ(ierr); 923fe394131Sjeremylt ierr = VecDestroy(&NBCsloc); CHKERRQ(ierr); 924ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 925ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 926ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 927ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 928a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 929a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 930ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 931ccaff030SJeremy L Thompson 932ccaff030SJeremy L Thompson // Structs 933ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 934ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 935ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 936ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 937ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 938f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 939ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 940ccaff030SJeremy L Thompson 941ccaff030SJeremy L Thompson return PetscFinalize(); 942ccaff030SJeremy L Thompson } 943