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 61*fe394131Sjeremylt 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 if (!appCtx->setMemTypeRequest) 11762e9c006SJeremy L Thompson appCtx->memTypeRequested = memTypeBackend; 11862e9c006SJeremy L Thompson else if (!appCtx->petscHaveCuda && appCtx->memTypeRequested == CEED_MEM_DEVICE) 11962e9c006SJeremy L Thompson SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 12062e9c006SJeremy L Thompson "PETSc was not built with CUDA. " 12162e9c006SJeremy L Thompson "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 12262e9c006SJeremy L Thompson 123777ff853SJeremy L Thompson // Wrap context in libCEED objects 124777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhys); 125777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhys, CEED_MEM_HOST, CEED_USE_POINTER, 126777ff853SJeremy L Thompson sizeof(*phys), phys); 127777ff853SJeremy L Thompson if (physSmoother) { 128777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxPhysSmoother); 129777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxPhysSmoother, CEED_MEM_HOST, CEED_USE_POINTER, 130777ff853SJeremy L Thompson sizeof(*physSmoother), physSmoother); 131777ff853SJeremy L Thompson } 132777ff853SJeremy L Thompson 13362e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 134ccaff030SJeremy L Thompson // Setup DM 135ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 136ccaff030SJeremy L Thompson // Performance logging 137ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 138ccaff030SJeremy L Thompson CHKERRQ(ierr); 139ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 140ccaff030SJeremy L Thompson 141ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 142ccaff030SJeremy L Thompson ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 143ccaff030SJeremy L Thompson 144ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 145ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 146e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 147ccaff030SJeremy L Thompson ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 14862e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 14962e9c006SJeremy L Thompson ierr = DMSetVecType(levelDMs[level], VECCUDA); CHKERRQ(ierr); 15062e9c006SJeremy L Thompson } 151ccaff030SJeremy L Thompson ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 152a3c02c40SJeremy L Thompson PETSC_TRUE, ncompu); CHKERRQ(ierr); 153a3c02c40SJeremy L Thompson // -- Label field components for viewing 154a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 155a3c02c40SJeremy L Thompson PetscSection section; 156a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(levelDMs[level], §ion); CHKERRQ(ierr); 157a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 158a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 159a3c02c40SJeremy L Thompson CHKERRQ(ierr); 160a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 161a3c02c40SJeremy L Thompson CHKERRQ(ierr); 162a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 163a3c02c40SJeremy L Thompson CHKERRQ(ierr); 164a3c02c40SJeremy L Thompson } 165a3c02c40SJeremy L Thompson 166a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 167a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmEnergy); CHKERRQ(ierr); 168a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmEnergy, appCtx, appCtx->levelDegrees[fineLevel], 1695c25879aSJeremy L Thompson PETSC_FALSE, ncompe); CHKERRQ(ierr); 170a3c02c40SJeremy L Thompson ierr = DMClone(dmOrig, &dmDiagnostic); CHKERRQ(ierr); 171a3c02c40SJeremy L Thompson ierr = SetupDMByDegree(dmDiagnostic, appCtx, appCtx->levelDegrees[fineLevel], 1728ca58ff3SJeremy L Thompson PETSC_FALSE, ncompu + ncompd); CHKERRQ(ierr); 17362e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 17462e9c006SJeremy L Thompson ierr = DMSetVecType(dmEnergy, VECCUDA); CHKERRQ(ierr); 17562e9c006SJeremy L Thompson ierr = DMSetVecType(dmDiagnostic, VECCUDA); CHKERRQ(ierr); 17662e9c006SJeremy L Thompson } 177a3c02c40SJeremy L Thompson { 178a3c02c40SJeremy L Thompson // -- Label field components for viewing 179a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 180a3c02c40SJeremy L Thompson PetscSection section; 181a3c02c40SJeremy L Thompson ierr = DMGetLocalSection(dmDiagnostic, §ion); CHKERRQ(ierr); 182a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 1838ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 1848ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1858ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 1868ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1878ca58ff3SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 1888ca58ff3SJeremy L Thompson CHKERRQ(ierr); 189379387d4SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure"); 1908ca58ff3SJeremy L Thompson CHKERRQ(ierr); 1917ab18febSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"); 19213fdad4bSJeremy L Thompson CHKERRQ(ierr); 19313fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 19413fdad4bSJeremy L Thompson CHKERRQ(ierr); 19513fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 19613fdad4bSJeremy L Thompson CHKERRQ(ierr); 19713fdad4bSJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 198a3c02c40SJeremy L Thompson CHKERRQ(ierr); 199ccaff030SJeremy L Thompson } 200ccaff030SJeremy L Thompson 201ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 202ccaff030SJeremy L Thompson // Setup solution and work vectors 203ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 204ccaff030SJeremy L Thompson // Allocate arrays 205ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 206ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 207ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 208ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 209ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 210ccaff030SJeremy L Thompson 211ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 212e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 213ccaff030SJeremy L Thompson // -- Create global unknown vector U 214ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 215ccaff030SJeremy L Thompson ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 216ccaff030SJeremy L Thompson // Note: Local size for matShell 217ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 218ccaff030SJeremy L Thompson 219ccaff030SJeremy L Thompson // -- Create local unknown vector Uloc 220ccaff030SJeremy L Thompson ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 221ccaff030SJeremy L Thompson // Note: local size for libCEED 222ccaff030SJeremy L Thompson ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 223ccaff030SJeremy L Thompson } 224ccaff030SJeremy L Thompson 225ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 226ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 227ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 228ccaff030SJeremy L Thompson ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 229ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 230ccaff030SJeremy L Thompson ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 231ccaff030SJeremy L Thompson 232ccaff030SJeremy L Thompson // Performance logging 233ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 234ccaff030SJeremy L Thompson 235ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 236ccaff030SJeremy L Thompson // Set up libCEED 237ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 238ccaff030SJeremy L Thompson // Performance logging 239ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 240ccaff030SJeremy L Thompson CHKERRQ(ierr); 241ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 242ccaff030SJeremy L Thompson 243ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 244ccaff030SJeremy L Thompson CeedVector forceCeed; 245ccaff030SJeremy L Thompson CeedScalar *f; 246ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 24762e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 248ccaff030SJeremy L Thompson ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 24962e9c006SJeremy L Thompson } else { 25062e9c006SJeremy L Thompson ierr = VecCUDAGetArray(Floc, &f); CHKERRQ(ierr); 25162e9c006SJeremy L Thompson } 252ccaff030SJeremy L Thompson CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 25362e9c006SJeremy L Thompson CeedVectorSetArray(forceCeed, appCtx->memTypeRequested, CEED_USE_POINTER, f); 254ccaff030SJeremy L Thompson } 255ccaff030SJeremy L Thompson 256*fe394131Sjeremylt // -- Create libCEED local Neumann BCs vector 257*fe394131Sjeremylt CeedVector neumannCeed; 258*fe394131Sjeremylt CeedScalar *n; 259*fe394131Sjeremylt if (appCtx->bcTractionCount > 0) { 260*fe394131Sjeremylt ierr = VecDuplicate(U, &NBCs); CHKERRQ(ierr); 261*fe394131Sjeremylt ierr = VecDuplicate(Uloc[fineLevel], &NBCsloc); CHKERRQ(ierr); 262*fe394131Sjeremylt if (appCtx->memTypeRequested == CEED_MEM_HOST) { 263*fe394131Sjeremylt ierr = VecGetArray(NBCsloc, &n); CHKERRQ(ierr); 264*fe394131Sjeremylt } else { 265*fe394131Sjeremylt ierr = VecCUDAGetArray(NBCsloc, &n); CHKERRQ(ierr); 266*fe394131Sjeremylt } 267*fe394131Sjeremylt CeedVectorCreate(ceed, Ulocsz[fineLevel], &neumannCeed); 268*fe394131Sjeremylt CeedVectorSetArray(neumannCeed, appCtx->memTypeRequested, 269*fe394131Sjeremylt CEED_USE_POINTER, n); 270*fe394131Sjeremylt } 271*fe394131Sjeremylt 272ccaff030SJeremy L Thompson // -- Setup libCEED objects 273ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 274d99fa3c5SJeremy L Thompson // ---- Setup residual, Jacobian evaluator and geometric information 275ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 276ccaff030SJeremy L Thompson { 277a3c02c40SJeremy L Thompson ierr = SetupLibceedFineLevel(levelDMs[fineLevel], dmEnergy, dmDiagnostic, 278777ff853SJeremy L Thompson ceed, appCtx, ctxPhys, ceedData, fineLevel, 279a3c02c40SJeremy L Thompson ncompu, Ugsz[fineLevel], Ulocsz[fineLevel], 280*fe394131Sjeremylt forceCeed, neumannCeed); 281ccaff030SJeremy L Thompson CHKERRQ(ierr); 282ccaff030SJeremy L Thompson } 283d99fa3c5SJeremy L Thompson // ---- Setup coarse Jacobian evaluator and prolongation/restriction 284d99fa3c5SJeremy L Thompson for (PetscInt level = numLevels - 2; level >= 0; level--) { 285ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 286d99fa3c5SJeremy L Thompson 287d99fa3c5SJeremy L Thompson // Get global communication restriction 288d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 289d99fa3c5SJeremy L Thompson ierr = VecSet(Uloc[level+1], 1.0); CHKERRQ(ierr); 290d99fa3c5SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[level+1], Uloc[level+1], ADD_VALUES, 291d99fa3c5SJeremy L Thompson Ug[level+1]); CHKERRQ(ierr); 292d99fa3c5SJeremy L Thompson ierr = DMGlobalToLocal(levelDMs[level+1], Ug[level+1], INSERT_VALUES, 293d99fa3c5SJeremy L Thompson Uloc[level+1]); CHKERRQ(ierr); 294d99fa3c5SJeremy L Thompson 295d99fa3c5SJeremy L Thompson // Place in libCEED array 296d99fa3c5SJeremy L Thompson const PetscScalar *m; 297d99fa3c5SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 298d99fa3c5SJeremy L Thompson ierr = VecGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 299d99fa3c5SJeremy L Thompson } else { 300d99fa3c5SJeremy L Thompson ierr = VecCUDAGetArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 301ccaff030SJeremy L Thompson } 302d99fa3c5SJeremy L Thompson CeedVectorSetArray(ceedData[level+1]->xceed, appCtx->memTypeRequested, 303d99fa3c5SJeremy L Thompson CEED_USE_POINTER, (CeedScalar *)m); 304ccaff030SJeremy L Thompson 3051c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 306777ff853SJeremy L Thompson ierr = SetupLibceedLevel(levelDMs[level], ceed, appCtx, 307ccaff030SJeremy L Thompson ceedData, level, ncompu, Ugsz[level], 308777ff853SJeremy L Thompson Ulocsz[level], ceedData[level+1]->xceed); 309777ff853SJeremy L Thompson CHKERRQ(ierr); 310d99fa3c5SJeremy L Thompson 311d99fa3c5SJeremy L Thompson // Restore PETSc vector 312d99fa3c5SJeremy L Thompson CeedVectorTakeArray(ceedData[level+1]->xceed, appCtx->memTypeRequested, 313d99fa3c5SJeremy L Thompson (CeedScalar **)&m); 314d99fa3c5SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 315d99fa3c5SJeremy L Thompson ierr = VecRestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 316d99fa3c5SJeremy L Thompson } else { 317d99fa3c5SJeremy L Thompson ierr = VecCUDARestoreArrayRead(Uloc[level+1], &m); CHKERRQ(ierr); 318d99fa3c5SJeremy L Thompson } 319d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Ug[level+1]); CHKERRQ(ierr); 320d99fa3c5SJeremy L Thompson ierr = VecZeroEntries(Uloc[level+1]); CHKERRQ(ierr); 321ccaff030SJeremy L Thompson } 322ccaff030SJeremy L Thompson 323ccaff030SJeremy L Thompson // Performance logging 324ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 325ccaff030SJeremy L Thompson 326ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 327*fe394131Sjeremylt // Setup global forcing and Neumann BC vectors 328ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 329ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 330ccaff030SJeremy L Thompson 331ccaff030SJeremy L Thompson if (appCtx->forcingChoice != FORCE_NONE) { 3326a6c615bSJeremy L Thompson CeedVectorTakeArray(forceCeed, appCtx->memTypeRequested, NULL); 33362e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 334ccaff030SJeremy L Thompson ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 33562e9c006SJeremy L Thompson } else { 33662e9c006SJeremy L Thompson ierr = VecCUDARestoreArray(Floc, &f); CHKERRQ(ierr); 33762e9c006SJeremy L Thompson } 338ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 339ccaff030SJeremy L Thompson CHKERRQ(ierr); 340ccaff030SJeremy L Thompson CeedVectorDestroy(&forceCeed); 341ccaff030SJeremy L Thompson } 342ccaff030SJeremy L Thompson 343*fe394131Sjeremylt if (appCtx->bcTractionCount > 0) { 344*fe394131Sjeremylt ierr = VecZeroEntries(NBCs); CHKERRQ(ierr); 345*fe394131Sjeremylt CeedVectorTakeArray(neumannCeed, appCtx->memTypeRequested, NULL); 346*fe394131Sjeremylt if (appCtx->memTypeRequested == CEED_MEM_HOST) { 347*fe394131Sjeremylt ierr = VecRestoreArray(NBCsloc, &n); CHKERRQ(ierr); 348*fe394131Sjeremylt } else { 349*fe394131Sjeremylt ierr = VecCUDARestoreArray(NBCsloc, &n); CHKERRQ(ierr); 350*fe394131Sjeremylt } 351*fe394131Sjeremylt ierr = DMLocalToGlobal(levelDMs[fineLevel], NBCsloc, ADD_VALUES, NBCs); 352*fe394131Sjeremylt CHKERRQ(ierr); 353*fe394131Sjeremylt CeedVectorDestroy(&neumannCeed); 354*fe394131Sjeremylt } 355*fe394131Sjeremylt 356ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 357ccaff030SJeremy L Thompson // Print problem summary 358ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 359ccaff030SJeremy L Thompson if (!appCtx->testMode) { 360ccaff030SJeremy L Thompson const char *usedresource; 361ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 362ccaff030SJeremy L Thompson 363ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 36417fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 365ccaff030SJeremy L Thompson " libCEED:\n" 36662e9c006SJeremy L Thompson " libCEED Backend : %s\n" 36762e9c006SJeremy L Thompson " libCEED Backend MemType : %s\n" 36862e9c006SJeremy L Thompson " libCEED User Requested MemType : %s\n", 36962e9c006SJeremy L Thompson usedresource, CeedMemTypes[memTypeBackend], 37062e9c006SJeremy L Thompson (appCtx->setMemTypeRequest) ? 37162e9c006SJeremy L Thompson CeedMemTypes[appCtx->memTypeRequested] : "none"); 37262e9c006SJeremy L Thompson CHKERRQ(ierr); 373ccaff030SJeremy L Thompson 37462e9c006SJeremy L Thompson VecType vecType; 37562e9c006SJeremy L Thompson ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 37662e9c006SJeremy L Thompson ierr = PetscPrintf(comm, 37762e9c006SJeremy L Thompson " PETSc:\n" 37862e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 37962e9c006SJeremy L Thompson vecType); CHKERRQ(ierr); 38062e9c006SJeremy L Thompson 381ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 382ccaff030SJeremy L Thompson " Problem:\n" 383ccaff030SJeremy L Thompson " Problem Name : %s\n" 384ccaff030SJeremy L Thompson " Forcing Function : %s\n" 385ccaff030SJeremy L Thompson " Mesh:\n" 386ccaff030SJeremy L Thompson " File : %s\n" 387ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 388ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 389ccaff030SJeremy L Thompson " Global nodes : %D\n" 390ccaff030SJeremy L Thompson " Owned nodes : %D\n" 391ccaff030SJeremy L Thompson " DoF per node : %D\n" 392ccaff030SJeremy L Thompson " Multigrid:\n" 393ccaff030SJeremy L Thompson " Type : %s\n" 394ccaff030SJeremy L Thompson " Number of Levels : %d\n", 395ccaff030SJeremy L Thompson problemTypesForDisp[appCtx->problemChoice], 396ccaff030SJeremy L Thompson forcingTypesForDisp[appCtx->forcingChoice], 397ccaff030SJeremy L Thompson appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 398ccaff030SJeremy L Thompson appCtx->degree + 1, appCtx->degree + 1, 399ccaff030SJeremy L Thompson Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 400f892e6d1Sjeremylt (appCtx->degree == 1 && 401f892e6d1Sjeremylt appCtx->multigridChoice != MULTIGRID_NONE) ? 402f892e6d1Sjeremylt "Algebraic multigrid" : 403ccaff030SJeremy L Thompson multigridTypesForDisp[appCtx->multigridChoice], 404f892e6d1Sjeremylt (appCtx->degree == 1 || 405f892e6d1Sjeremylt appCtx->multigridChoice == MULTIGRID_NONE) ? 406f892e6d1Sjeremylt 0 : numLevels); CHKERRQ(ierr); 407ccaff030SJeremy L Thompson 408ccaff030SJeremy L Thompson if (appCtx->multigridChoice != MULTIGRID_NONE) { 409e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 410ccaff030SJeremy L Thompson CeedInt level = i ? fineLevel : 0; 411ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 412ccaff030SJeremy L Thompson " Level %D (%s):\n" 413ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 414ccaff030SJeremy L Thompson " Global Nodes : %D\n" 415ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 416ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 417ccaff030SJeremy L Thompson appCtx->levelDegrees[level] + 1, 418ccaff030SJeremy L Thompson Ugsz[level]/ncompu, Ulsz[level]/ncompu); 419ccaff030SJeremy L Thompson CHKERRQ(ierr); 420ccaff030SJeremy L Thompson } 421ccaff030SJeremy L Thompson } 422ccaff030SJeremy L Thompson } 423ccaff030SJeremy L Thompson 424ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 425ccaff030SJeremy L Thompson // Setup SNES 426ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 427ccaff030SJeremy L Thompson // Performance logging 428ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 429ccaff030SJeremy L Thompson CHKERRQ(ierr); 430ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 431ccaff030SJeremy L Thompson 432ccaff030SJeremy L Thompson // Create SNES 433ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 434ccaff030SJeremy L Thompson ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 435ccaff030SJeremy L Thompson 436ccaff030SJeremy L Thompson // -- Jacobian evaluators 437ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 438ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 439e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 440ccaff030SJeremy L Thompson // -- Jacobian context for level 441ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 442ccaff030SJeremy L Thompson ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 443777ff853SJeremy L Thompson Uloc[level], ceedData[level], ceed, ctxPhys, 444777ff853SJeremy L Thompson ctxPhysSmoother, jacobCtx[level]); CHKERRQ(ierr); 445ccaff030SJeremy L Thompson 446ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 447ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 448ccaff030SJeremy L Thompson Ugsz[level], jacobCtx[level], &jacobMat[level]); 449ccaff030SJeremy L Thompson CHKERRQ(ierr); 450ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 451ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 452ccaff030SJeremy L Thompson CHKERRQ(ierr); 453ccaff030SJeremy L Thompson ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 454ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 455a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 456a3658badSJeremy L Thompson ierr = MatShellSetVecType(jacobMat[level], VECCUDA); CHKERRQ(ierr); 457a3658badSJeremy L Thompson } 458ccaff030SJeremy L Thompson } 459e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 460ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 461ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 462ccaff030SJeremy L Thompson formJacobCtx->jacobCtx = jacobCtx; 463ccaff030SJeremy L Thompson formJacobCtx->numLevels = numLevels; 464ccaff030SJeremy L Thompson formJacobCtx->jacobMat = jacobMat; 465ccaff030SJeremy L Thompson 466ccaff030SJeremy L Thompson // -- Residual evaluation function 467*fe394131Sjeremylt ierr = PetscCalloc1(1, &resCtx); CHKERRQ(ierr); 468ccaff030SJeremy L Thompson ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 469ccaff030SJeremy L Thompson sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 470ccaff030SJeremy L Thompson resCtx->op = ceedData[fineLevel]->opApply; 471f7b4142eSJeremy L Thompson resCtx->qf = ceedData[fineLevel]->qfApply; 472*fe394131Sjeremylt if (appCtx->bcTractionCount > 0) 473*fe394131Sjeremylt resCtx->NBCs = NBCs; 474*fe394131Sjeremylt else 475*fe394131Sjeremylt resCtx->NBCs = NULL; 476ccaff030SJeremy L Thompson ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 477ccaff030SJeremy L Thompson 478ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 479ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 480ccaff030SJeremy L Thompson ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 481e3e3df41Sjeremylt for (PetscInt level = 1; level < numLevels; level++) { 482ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 483ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 48462e9c006SJeremy L Thompson ierr = SetupProlongRestrictCtx(comm, appCtx, levelDMs[level-1], 48562e9c006SJeremy L Thompson levelDMs[level], Ug[level], Uloc[level-1], 48662e9c006SJeremy L Thompson Uloc[level], ceedData[level-1], 48762e9c006SJeremy L Thompson ceedData[level], ceed, 488ccaff030SJeremy L Thompson prolongRestrCtx[level]); CHKERRQ(ierr); 489ccaff030SJeremy L Thompson 490ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 491ccaff030SJeremy L Thompson ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 492ccaff030SJeremy L Thompson Ugsz[level-1], prolongRestrCtx[level], 493ccaff030SJeremy L Thompson &prolongRestrMat[level]); CHKERRQ(ierr); 494ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 495ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 496ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 497ccaff030SJeremy L Thompson ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 498ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 499ccaff030SJeremy L Thompson CHKERRQ(ierr); 500a3658badSJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_DEVICE) { 501a3658badSJeremy L Thompson ierr = MatShellSetVecType(prolongRestrMat[level], VECCUDA); CHKERRQ(ierr); 502a3658badSJeremy L Thompson } 503ccaff030SJeremy L Thompson } 504ccaff030SJeremy L Thompson 505ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 506ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 507ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 508e3e3df41Sjeremylt if (appCtx->multigridChoice != MULTIGRID_NONE) { 509e3e3df41Sjeremylt // -- Jacobian Matrix 510e3e3df41Sjeremylt ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 511e3e3df41Sjeremylt ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 512e3e3df41Sjeremylt 513e3e3df41Sjeremylt if (appCtx->degree > 1) { 514ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 515ccaff030SJeremy L Thompson ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 516ccaff030SJeremy L Thompson ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 517ccaff030SJeremy L Thompson 518e3e3df41Sjeremylt // -- Jacobian function 519ccaff030SJeremy L Thompson ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 520ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 521ccaff030SJeremy L Thompson 522ccaff030SJeremy L Thompson // -- Residual evaluation function 523ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 524ccaff030SJeremy L Thompson ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 525ccaff030SJeremy L Thompson CHKERRQ(ierr); 526ccaff030SJeremy L Thompson ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 527ccaff030SJeremy L Thompson jacobCoarseCtx); CHKERRQ(ierr); 528ccaff030SJeremy L Thompson 529ccaff030SJeremy L Thompson // -- Update formJacobCtx 530ccaff030SJeremy L Thompson formJacobCtx->Ucoarse = Ug[0]; 531ccaff030SJeremy L Thompson formJacobCtx->snesCoarse = snesCoarse; 532ccaff030SJeremy L Thompson formJacobCtx->jacobMatCoarse = jacobMatCoarse; 533e3e3df41Sjeremylt } 534e3e3df41Sjeremylt } 535e3e3df41Sjeremylt 536e3e3df41Sjeremylt // Set Jacobian function 537e3e3df41Sjeremylt if (appCtx->degree > 1) { 538e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 539e3e3df41Sjeremylt FormJacobian, formJacobCtx); CHKERRQ(ierr); 540e3e3df41Sjeremylt } else { 541e3e3df41Sjeremylt ierr = SNESSetJacobian(snes, jacobMat[0], jacobMatCoarse, 542e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 543e3e3df41Sjeremylt CHKERRQ(ierr); 544e3e3df41Sjeremylt } 545ccaff030SJeremy L Thompson 546ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 547ccaff030SJeremy L Thompson // Setup KSP 548ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 549ccaff030SJeremy L Thompson { 550ccaff030SJeremy L Thompson PC pc; 551ccaff030SJeremy L Thompson KSP ksp; 552ccaff030SJeremy L Thompson 553ccaff030SJeremy L Thompson // -- KSP 554ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 555ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 556ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 557ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 558ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 559ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 560ccaff030SJeremy L Thompson 561ccaff030SJeremy L Thompson // -- Preconditioning 562ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 563ccaff030SJeremy L Thompson ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 564ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 565ccaff030SJeremy L Thompson 566ccaff030SJeremy L Thompson if (appCtx->multigridChoice == MULTIGRID_NONE) { 567ccaff030SJeremy L Thompson // ---- No Multigrid 568ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 569ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 570f892e6d1Sjeremylt } else if (appCtx->degree == 1) { 571f892e6d1Sjeremylt // ---- AMG for degree 1 572f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 573ccaff030SJeremy L Thompson } else { 574ccaff030SJeremy L Thompson // ---- PCMG 575ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 576ccaff030SJeremy L Thompson 577ccaff030SJeremy L Thompson // ------ PCMG levels 578ccaff030SJeremy L Thompson ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 579e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 580ccaff030SJeremy L Thompson // -------- Smoother 581ccaff030SJeremy L Thompson KSP kspSmoother, kspEst; 582ccaff030SJeremy L Thompson PC pcSmoother; 583ccaff030SJeremy L Thompson 584ccaff030SJeremy L Thompson // ---------- Smoother KSP 585ccaff030SJeremy L Thompson ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 586ccaff030SJeremy L Thompson ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 587ccaff030SJeremy L Thompson ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 588ccaff030SJeremy L Thompson 589ccaff030SJeremy L Thompson // ---------- Chebyshev options 590ccaff030SJeremy L Thompson ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 591ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 592ccaff030SJeremy L Thompson CHKERRQ(ierr); 593ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 594ccaff030SJeremy L Thompson ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 595ccaff030SJeremy L Thompson ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 596ccaff030SJeremy L Thompson CHKERRQ(ierr); 597ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 598ccaff030SJeremy L Thompson CHKERRQ(ierr); 599ccaff030SJeremy L Thompson 600ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 601ccaff030SJeremy L Thompson ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 602ccaff030SJeremy L Thompson ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 604ccaff030SJeremy L Thompson 605ccaff030SJeremy L Thompson // -------- Work vector 606ccaff030SJeremy L Thompson if (level != fineLevel) { 607ccaff030SJeremy L Thompson ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 608ccaff030SJeremy L Thompson } 609ccaff030SJeremy L Thompson 610ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 611ccaff030SJeremy L Thompson if (level > 0) { 612ccaff030SJeremy L Thompson ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 613ccaff030SJeremy L Thompson CHKERRQ(ierr); 614ccaff030SJeremy L Thompson ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 615ccaff030SJeremy L Thompson CHKERRQ(ierr); 616ccaff030SJeremy L Thompson } 617ccaff030SJeremy L Thompson } 618ccaff030SJeremy L Thompson 619ccaff030SJeremy L Thompson // ------ PCMG coarse solve 620ccaff030SJeremy L Thompson KSP kspCoarse; 621ccaff030SJeremy L Thompson PC pcCoarse; 622ccaff030SJeremy L Thompson 623ccaff030SJeremy L Thompson // -------- Coarse KSP 624ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 625ccaff030SJeremy L Thompson ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 626ccaff030SJeremy L Thompson ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 627ccaff030SJeremy L Thompson CHKERRQ(ierr); 628ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 629ccaff030SJeremy L Thompson 630ccaff030SJeremy L Thompson // -------- Coarse preconditioner 631ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 632ccaff030SJeremy L Thompson ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 633ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 634ccaff030SJeremy L Thompson 635ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 636ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 637ccaff030SJeremy L Thompson 638ccaff030SJeremy L Thompson // ------ PCMG options 639ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 640ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 641ccaff030SJeremy L Thompson ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 642ccaff030SJeremy L Thompson } 643ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 644ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 645ccaff030SJeremy L Thompson } 646038d0cd7Sjeremylt { 647038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 648481a4840SJed Brown SNESLineSearch linesearch; 649481a4840SJed Brown 650481a4840SJed Brown ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 651481a4840SJed Brown ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 652481a4840SJed Brown } 653481a4840SJed Brown 654ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 655ccaff030SJeremy L Thompson 656ccaff030SJeremy L Thompson // Performance logging 657ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 658ccaff030SJeremy L Thompson 659ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 660ccaff030SJeremy L Thompson // Set initial guess 661ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 662f81c27eaSJed Brown ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 663ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 664ccaff030SJeremy L Thompson 665ccaff030SJeremy L Thompson // View solution 666ccaff030SJeremy L Thompson if (appCtx->viewSoln) { 667ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 668ccaff030SJeremy L Thompson } 669ccaff030SJeremy L Thompson 670ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 671ccaff030SJeremy L Thompson // Solve SNES 672ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 6735f24f028Sjeremylt PetscBool snesMonitor = PETSC_FALSE; 6745f24f028Sjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 6755f24f028Sjeremylt CHKERRQ(ierr); 6765f24f028Sjeremylt 677ccaff030SJeremy L Thompson // Performance logging 678ccaff030SJeremy L Thompson ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 679ccaff030SJeremy L Thompson CHKERRQ(ierr); 680ccaff030SJeremy L Thompson ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 681ccaff030SJeremy L Thompson 682ccaff030SJeremy L Thompson // Timing 683ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 684ccaff030SJeremy L Thompson startTime = MPI_Wtime(); 685ccaff030SJeremy L Thompson 686ccaff030SJeremy L Thompson // Solve for each load increment 6875f24f028Sjeremylt PetscInt increment; 6885f24f028Sjeremylt for (increment = 1; increment <= appCtx->numIncrements; increment++) { 6895f24f028Sjeremylt // -- Log increment count 6905f24f028Sjeremylt if (snesMonitor) { 6915f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 6925f24f028Sjeremylt CHKERRQ(ierr); 6935f24f028Sjeremylt } 6945f24f028Sjeremylt 695ccaff030SJeremy L Thompson // -- Scale the problem 696ccaff030SJeremy L Thompson PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 697ccaff030SJeremy L Thompson scalingFactor = loadIncrement / 698ccaff030SJeremy L Thompson (increment == 1 ? 1 : resCtx->loadIncrement); 699ccaff030SJeremy L Thompson resCtx->loadIncrement = loadIncrement; 700ccaff030SJeremy L Thompson if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 701ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 702ccaff030SJeremy L Thompson } 703ccaff030SJeremy L Thompson 704ccaff030SJeremy L Thompson // -- Solve 705ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 706ccaff030SJeremy L Thompson 707ccaff030SJeremy L Thompson // -- View solution 708560e6f00SJeremy L Thompson if (appCtx->viewSoln) { 709ccaff030SJeremy L Thompson ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 710ccaff030SJeremy L Thompson } 711ccaff030SJeremy L Thompson 712ccaff030SJeremy L Thompson // -- Update SNES iteration count 713ccaff030SJeremy L Thompson PetscInt its; 714ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 715ccaff030SJeremy L Thompson snesIts += its; 7167418ca88SJeremy L Thompson ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr); 7177418ca88SJeremy L Thompson kspIts += its; 7183951731eSjeremylt 7193951731eSjeremylt // -- Check for divergence 7203951731eSjeremylt SNESConvergedReason reason; 7213951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 7223951731eSjeremylt if (reason < 0) 7233951731eSjeremylt break; 724ccaff030SJeremy L Thompson } 725ccaff030SJeremy L Thompson 726ccaff030SJeremy L Thompson // Timing 727ccaff030SJeremy L Thompson elapsedTime = MPI_Wtime() - startTime; 728ccaff030SJeremy L Thompson 729ccaff030SJeremy L Thompson // Performance logging 730ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 731ccaff030SJeremy L Thompson 732ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 733ccaff030SJeremy L Thompson // Output summary 734ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 735ccaff030SJeremy L Thompson if (!appCtx->testMode) { 736ccaff030SJeremy L Thompson // -- SNES 737ccaff030SJeremy L Thompson SNESType snesType; 738ccaff030SJeremy L Thompson SNESConvergedReason reason; 739ccaff030SJeremy L Thompson PetscReal rnorm; 740ccaff030SJeremy L Thompson ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 741ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 742ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 743ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 744ccaff030SJeremy L Thompson " SNES:\n" 745ccaff030SJeremy L Thompson " SNES Type : %s\n" 746ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 747ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 7485f24f028Sjeremylt " Completed Load Increments : %d\n" 749ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 750ccaff030SJeremy L Thompson " Final rnorm : %e\n", 751ccaff030SJeremy L Thompson snesType, SNESConvergedReasons[reason], 7525f24f028Sjeremylt appCtx->numIncrements, increment - 1, 7535f24f028Sjeremylt snesIts, (double)rnorm); CHKERRQ(ierr); 754ccaff030SJeremy L Thompson 755ccaff030SJeremy L Thompson // -- KSP 756ccaff030SJeremy L Thompson KSP ksp; 757ccaff030SJeremy L Thompson KSPType kspType; 758ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 759ccaff030SJeremy L Thompson ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 760ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 761ccaff030SJeremy L Thompson " Linear Solver:\n" 7627418ca88SJeremy L Thompson " KSP Type : %s\n" 7637418ca88SJeremy L Thompson " Total KSP Iterations : %D\n", 7647418ca88SJeremy L Thompson kspType, kspIts); CHKERRQ(ierr); 765ccaff030SJeremy L Thompson 766ccaff030SJeremy L Thompson // -- PC 767ccaff030SJeremy L Thompson PC pc; 768e3e3df41Sjeremylt PCType pcType; 769ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 770e3e3df41Sjeremylt ierr = PCGetType(pc, &pcType); CHKERRQ(ierr); 771e3e3df41Sjeremylt ierr = PetscPrintf(comm, 772e3e3df41Sjeremylt " PC Type : %s\n", 773e3e3df41Sjeremylt pcType); CHKERRQ(ierr); 774e3e3df41Sjeremylt 7752b901a5eSjeremylt if (!strcmp(pcType, PCMG)) { 776e3e3df41Sjeremylt PCMGType pcmgType; 777ccaff030SJeremy L Thompson ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 778ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 779ccaff030SJeremy L Thompson " P-Multigrid:\n" 780ccaff030SJeremy L Thompson " PCMG Type : %s\n" 781ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 782ccaff030SJeremy L Thompson PCMGTypes[pcmgType], 783ccaff030SJeremy L Thompson PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 784ccaff030SJeremy L Thompson 785ccaff030SJeremy L Thompson // -- Coarse Solve 786ccaff030SJeremy L Thompson KSP kspCoarse; 787ccaff030SJeremy L Thompson PC pcCoarse; 788ccaff030SJeremy L Thompson PCType pcType; 789ccaff030SJeremy L Thompson 790ccaff030SJeremy L Thompson ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 791ccaff030SJeremy L Thompson ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 792ccaff030SJeremy L Thompson ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 793ccaff030SJeremy L Thompson ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 794ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 795ccaff030SJeremy L Thompson " Coarse Solve:\n" 796ccaff030SJeremy L Thompson " KSP Type : %s\n" 797ccaff030SJeremy L Thompson " PC Type : %s\n", 798ccaff030SJeremy L Thompson kspType, pcType); CHKERRQ(ierr); 799ccaff030SJeremy L Thompson } 800ccaff030SJeremy L Thompson } 801ccaff030SJeremy L Thompson 802ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 803ccaff030SJeremy L Thompson // Compute solve time 804ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 805ccaff030SJeremy L Thompson if (!appCtx->testMode) { 806ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 807ccaff030SJeremy L Thompson CHKERRQ(ierr); 808ccaff030SJeremy L Thompson ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 809ccaff030SJeremy L Thompson CHKERRQ(ierr); 810ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 811ccaff030SJeremy L Thompson " Performance:\n" 8127418ca88SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n" 8137418ca88SJeremy L Thompson " DoFs/Sec in SNES : %g (%g) million\n", 8147418ca88SJeremy L Thompson maxTime, minTime, 1e-6*Ugsz[fineLevel]*kspIts/maxTime, 8157418ca88SJeremy L Thompson 1e-6*Ugsz[fineLevel]*kspIts/minTime); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson } 817ccaff030SJeremy L Thompson 818ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 819ccaff030SJeremy L Thompson // Compute error 820ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 821ccaff030SJeremy L Thompson if (appCtx->forcingChoice == FORCE_MMS) { 822ccaff030SJeremy L Thompson CeedScalar l2Error = 1., l2Unorm = 1.; 823ccaff030SJeremy L Thompson const CeedScalar *truearray; 824ccaff030SJeremy L Thompson Vec errorVec, trueVec; 825ccaff030SJeremy L Thompson 826ccaff030SJeremy L Thompson // -- Work vectors 827ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 828ccaff030SJeremy L Thompson ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 829ccaff030SJeremy L Thompson ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 830ccaff030SJeremy L Thompson ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 831ccaff030SJeremy L Thompson 832ccaff030SJeremy L Thompson // -- Assemble global true solution vector 83362e9c006SJeremy L Thompson CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, 83462e9c006SJeremy L Thompson appCtx->memTypeRequested, &truearray); 83562e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 83662e9c006SJeremy L Thompson ierr = VecPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 83762e9c006SJeremy L Thompson CHKERRQ(ierr); 83862e9c006SJeremy L Thompson } else { 83962e9c006SJeremy L Thompson ierr = VecCUDAPlaceArray(resCtx->Yloc, (PetscScalar *)truearray); 84062e9c006SJeremy L Thompson CHKERRQ(ierr); 84162e9c006SJeremy L Thompson } 842ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 843ccaff030SJeremy L Thompson CHKERRQ(ierr); 84462e9c006SJeremy L Thompson if (appCtx->memTypeRequested == CEED_MEM_HOST) { 845ccaff030SJeremy L Thompson ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 84662e9c006SJeremy L Thompson } else { 84762e9c006SJeremy L Thompson ierr = VecCUDAResetArray(resCtx->Yloc); CHKERRQ(ierr); 84862e9c006SJeremy L Thompson } 849ccaff030SJeremy L Thompson CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 850ccaff030SJeremy L Thompson 851ccaff030SJeremy L Thompson // -- Compute L2 error 852ccaff030SJeremy L Thompson ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 853ccaff030SJeremy L Thompson ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 854ccaff030SJeremy L Thompson ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 855ccaff030SJeremy L Thompson l2Error /= l2Unorm; 856ccaff030SJeremy L Thompson 857ccaff030SJeremy L Thompson // -- Output 858ccaff030SJeremy L Thompson if (!appCtx->testMode || l2Error > 0.05) { 859ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 860ccaff030SJeremy L Thompson " L2 Error : %e\n", 861ccaff030SJeremy L Thompson l2Error); CHKERRQ(ierr); 862ccaff030SJeremy L Thompson } 863ccaff030SJeremy L Thompson 864ccaff030SJeremy L Thompson // -- Cleanup 865ccaff030SJeremy L Thompson ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 866ccaff030SJeremy L Thompson ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 8672d93065eSjeremylt } 8682d93065eSjeremylt 8692d93065eSjeremylt // --------------------------------------------------------------------------- 8702d93065eSjeremylt // Compute energy 8712d93065eSjeremylt // --------------------------------------------------------------------------- 8722d93065eSjeremylt if (!appCtx->testMode) { 8732d93065eSjeremylt // -- Compute L2 error 8742d93065eSjeremylt CeedScalar energy; 875a3c02c40SJeremy L Thompson ierr = ComputeStrainEnergy(dmEnergy, resCtx, ceedData[fineLevel]->opEnergy, 876a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 8772d93065eSjeremylt 8782d93065eSjeremylt // -- Output 8792d93065eSjeremylt ierr = PetscPrintf(comm, 8802d93065eSjeremylt " Strain Energy : %e\n", 8812d93065eSjeremylt energy); CHKERRQ(ierr); 882ccaff030SJeremy L Thompson } 883ccaff030SJeremy L Thompson 884ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 8855c25879aSJeremy L Thompson // Output diagnostic quantities 8865c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 8875c25879aSJeremy L Thompson if (appCtx->viewSoln || appCtx->viewFinalSoln) { 8885c25879aSJeremy L Thompson // -- Setup context 8895c25879aSJeremy L Thompson UserMult diagnosticCtx; 8905c25879aSJeremy L Thompson ierr = PetscMalloc1(1, &diagnosticCtx); CHKERRQ(ierr); 8915c25879aSJeremy L Thompson ierr = PetscMemcpy(diagnosticCtx, resCtx, sizeof(*resCtx)); CHKERRQ(ierr); 8925c25879aSJeremy L Thompson diagnosticCtx->dm = dmDiagnostic; 8935c25879aSJeremy L Thompson diagnosticCtx->op = ceedData[fineLevel]->opDiagnostic; 8945c25879aSJeremy L Thompson 8955c25879aSJeremy L Thompson // -- Compute and output 8965c25879aSJeremy L Thompson ierr = ViewDiagnosticQuantities(comm, levelDMs[fineLevel], diagnosticCtx, U, 8975c25879aSJeremy L Thompson ceedData[fineLevel]->ErestrictDiagnostic); 8985c25879aSJeremy L Thompson CHKERRQ(ierr); 8995c25879aSJeremy L Thompson 9005c25879aSJeremy L Thompson // -- Cleanup 9015c25879aSJeremy L Thompson ierr = PetscFree(diagnosticCtx); CHKERRQ(ierr); 9025c25879aSJeremy L Thompson } 9035c25879aSJeremy L Thompson 9045c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 905ccaff030SJeremy L Thompson // Free objects 906ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 907ccaff030SJeremy L Thompson // Data in arrays per level 908e3e3df41Sjeremylt for (PetscInt level = 0; level < numLevels; level++) { 909ccaff030SJeremy L Thompson // Vectors 910ccaff030SJeremy L Thompson ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 911ccaff030SJeremy L Thompson ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 912ccaff030SJeremy L Thompson 913ccaff030SJeremy L Thompson // Jacobian matrix and data 914ccaff030SJeremy L Thompson ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 915ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 916ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 917ccaff030SJeremy L Thompson 918ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 919ccaff030SJeremy L Thompson if (level > 0) { 920ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 921ccaff030SJeremy L Thompson ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 922ccaff030SJeremy L Thompson } 923ccaff030SJeremy L Thompson 924ccaff030SJeremy L Thompson // DM 925ccaff030SJeremy L Thompson ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 926ccaff030SJeremy L Thompson 927ccaff030SJeremy L Thompson // libCEED objects 928ccaff030SJeremy L Thompson ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 929ccaff030SJeremy L Thompson } 930ccaff030SJeremy L Thompson 931ccaff030SJeremy L Thompson // Arrays 932ccaff030SJeremy L Thompson ierr = PetscFree(Ug); CHKERRQ(ierr); 933ccaff030SJeremy L Thompson ierr = PetscFree(Uloc); CHKERRQ(ierr); 934ccaff030SJeremy L Thompson ierr = PetscFree(Ugsz); CHKERRQ(ierr); 935ccaff030SJeremy L Thompson ierr = PetscFree(Ulsz); CHKERRQ(ierr); 936ccaff030SJeremy L Thompson ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 937ccaff030SJeremy L Thompson ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 938ccaff030SJeremy L Thompson ierr = PetscFree(jacobMat); CHKERRQ(ierr); 939ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 940ccaff030SJeremy L Thompson ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 941ccaff030SJeremy L Thompson ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 942ccaff030SJeremy L Thompson ierr = PetscFree(ceedData); CHKERRQ(ierr); 943ccaff030SJeremy L Thompson 944ccaff030SJeremy L Thompson // libCEED objects 945777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhys); 946777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxPhysSmoother); 947ccaff030SJeremy L Thompson CeedDestroy(&ceed); 948ccaff030SJeremy L Thompson 949ccaff030SJeremy L Thompson // PETSc objects 950ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 951ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 952ccaff030SJeremy L Thompson ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 953ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 954ccaff030SJeremy L Thompson ierr = VecDestroy(&Floc); CHKERRQ(ierr); 955*fe394131Sjeremylt ierr = VecDestroy(&NBCs); CHKERRQ(ierr); 956*fe394131Sjeremylt ierr = VecDestroy(&NBCsloc); CHKERRQ(ierr); 957ccaff030SJeremy L Thompson ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 958ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 959ccaff030SJeremy L Thompson ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 960ccaff030SJeremy L Thompson ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 961a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmEnergy); CHKERRQ(ierr); 962a3c02c40SJeremy L Thompson ierr = DMDestroy(&dmDiagnostic); CHKERRQ(ierr); 963ccaff030SJeremy L Thompson ierr = PetscFree(levelDMs); CHKERRQ(ierr); 964ccaff030SJeremy L Thompson 965ccaff030SJeremy L Thompson // Structs 966ccaff030SJeremy L Thompson ierr = PetscFree(resCtx); CHKERRQ(ierr); 967ccaff030SJeremy L Thompson ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 968ccaff030SJeremy L Thompson ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 969ccaff030SJeremy L Thompson ierr = PetscFree(appCtx); CHKERRQ(ierr); 970ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 971f7b4142eSJeremy L Thompson ierr = PetscFree(physSmoother); CHKERRQ(ierr); 972ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 973ccaff030SJeremy L Thompson 974ccaff030SJeremy L Thompson return PetscFinalize(); 975ccaff030SJeremy L Thompson } 976