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 // 30eccc2849SRezgar Shakeri // ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms 31eccc2849SRezgar Shakeri // ./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 SS-NH -forcing none -ceed /cpu/self 32eccc2849SRezgar Shakeri // ./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 FSInitial-NH1 -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 49d1d35e2fSjeremylt AppCtx app_ctx; // Contains problem options 50ccaff030SJeremy L Thompson Physics phys; // Contains physical constants 51d1d35e2fSjeremylt Physics phys_smoother = NULL; // Separate context if nu_smoother set 52ccaff030SJeremy L Thompson Units units; // Contains units scaling 53ccaff030SJeremy L Thompson // PETSc objects 54d1d35e2fSjeremylt PetscLogStage stage_dm_setup, stage_libceed_setup, 55d1d35e2fSjeremylt stage_snes_setup, stage_snes_solve; 56d1d35e2fSjeremylt DM dm_orig; // Distributed DM to clone 57d1d35e2fSjeremylt DM dm_energy, dm_diagnostic; // DMs for postprocessing 58d1d35e2fSjeremylt DM *level_dms; 59d1d35e2fSjeremylt Vec U, *U_g, *U_loc; // U: solution, R: residual, F: forcing 60d1d35e2fSjeremylt Vec R, R_loc, F, F_loc; // g: global, loc: local 61d1d35e2fSjeremylt Vec neumann_bcs = NULL, bcs_loc = NULL; 62d1d35e2fSjeremylt SNES snes, snes_coarse = NULL; 63d1d35e2fSjeremylt Mat *jacob_mat, jacob_mat_coarse, *prolong_restr_mat; 64ccaff030SJeremy L Thompson // PETSc data 65d1d35e2fSjeremylt UserMult res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx; 66d1d35e2fSjeremylt FormJacobCtx form_jacob_ctx; 67d1d35e2fSjeremylt UserMultProlongRestr *prolong_restr_ctx; 68d1d35e2fSjeremylt PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V; 69ccaff030SJeremy L Thompson // libCEED objects 70d99fa3c5SJeremy L Thompson Ceed ceed; 71d1d35e2fSjeremylt CeedData *ceed_data; 72d1d35e2fSjeremylt CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL; 73ccaff030SJeremy L Thompson // Parameters 74d1d35e2fSjeremylt PetscInt num_comp_u = 3; // 3 DoFs in 3D 75d1d35e2fSjeremylt PetscInt num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic 76d1d35e2fSjeremylt PetscInt num_levels = 1, fine_level = 0; 77*bf0f51feSjeremylt PetscInt *U_g_size, *U_l_size, *U_loc_size; 78*bf0f51feSjeremylt PetscInt snes_its = 0, ksp_its = 0; 79ccaff030SJeremy L Thompson // Timing 80d1d35e2fSjeremylt double start_time, elapsed_time, min_time, max_time; 81ccaff030SJeremy L Thompson 82ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 83*bf0f51feSjeremylt if (ierr) return ierr; 84ccaff030SJeremy L Thompson 85ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 86ccaff030SJeremy L Thompson // Process command line options 87ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 88ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 89ccaff030SJeremy L Thompson 90ccaff030SJeremy L Thompson // -- Set mesh file, polynomial degree, problem type 91d1d35e2fSjeremylt ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr); 92d1d35e2fSjeremylt ierr = ProcessCommandLineOptions(comm, app_ctx); CHKERRQ(ierr); 93d1d35e2fSjeremylt num_levels = app_ctx->num_levels; 94d1d35e2fSjeremylt fine_level = num_levels - 1; 95ccaff030SJeremy L Thompson 96ccaff030SJeremy L Thompson // -- Set Poison's ratio, Young's Modulus 97ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr); 98ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 99ccaff030SJeremy L Thompson ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr); 100d1d35e2fSjeremylt if (fabs(app_ctx->nu_smoother) > 1E-14) { 101d1d35e2fSjeremylt ierr = PetscMalloc1(1, &phys_smoother); CHKERRQ(ierr); 102d1d35e2fSjeremylt ierr = PetscMemcpy(phys_smoother, phys, sizeof(*phys)); CHKERRQ(ierr); 103d1d35e2fSjeremylt phys_smoother->nu = app_ctx->nu_smoother; 104f7b4142eSJeremy L Thompson } 105ccaff030SJeremy L Thompson 106ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 107d1d35e2fSjeremylt // Initialize libCEED 10862e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 109d1d35e2fSjeremylt // Initialize backend 110d1d35e2fSjeremylt CeedInit(app_ctx->ceed_resource, &ceed); 11162e9c006SJeremy L Thompson 11262e9c006SJeremy L Thompson // Check preferred MemType 113d1d35e2fSjeremylt CeedMemType mem_type_backend; 114d1d35e2fSjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 11562e9c006SJeremy L Thompson 116777ff853SJeremy L Thompson // Wrap context in libCEED objects 117d1d35e2fSjeremylt CeedQFunctionContextCreate(ceed, &ctx_phys); 118d1d35e2fSjeremylt CeedQFunctionContextSetData(ctx_phys, CEED_MEM_HOST, CEED_USE_POINTER, 119777ff853SJeremy L Thompson sizeof(*phys), phys); 120d1d35e2fSjeremylt if (phys_smoother) { 121d1d35e2fSjeremylt CeedQFunctionContextCreate(ceed, &ctx_phys_smoother); 122d1d35e2fSjeremylt CeedQFunctionContextSetData(ctx_phys_smoother, CEED_MEM_HOST, CEED_USE_POINTER, 123d1d35e2fSjeremylt sizeof(*phys_smoother), phys_smoother); 124777ff853SJeremy L Thompson } 125777ff853SJeremy L Thompson 12662e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 127ccaff030SJeremy L Thompson // Setup DM 128ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 129ccaff030SJeremy L Thompson // Performance logging 130d1d35e2fSjeremylt ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup); 131ccaff030SJeremy L Thompson CHKERRQ(ierr); 132d1d35e2fSjeremylt ierr = PetscLogStagePush(stage_dm_setup); CHKERRQ(ierr); 133ccaff030SJeremy L Thompson 134ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 135d1d35e2fSjeremylt ierr = CreateDistributedDM(comm, app_ctx, &dm_orig); CHKERRQ(ierr); 136b68a8d79SJed Brown VecType vectype; 137d1d35e2fSjeremylt switch (mem_type_backend) { 138b68a8d79SJed Brown case CEED_MEM_HOST: vectype = VECSTANDARD; break; 139b68a8d79SJed Brown case CEED_MEM_DEVICE: { 140b68a8d79SJed Brown const char *resolved; 141b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 142b68a8d79SJed Brown if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 143b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 144b68a8d79SJed Brown else vectype = VECSTANDARD; 145b68a8d79SJed Brown } 146b68a8d79SJed Brown } 147d1d35e2fSjeremylt ierr = DMSetVecType(dm_orig, vectype); CHKERRQ(ierr); 148d1d35e2fSjeremylt ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr); 149ccaff030SJeremy L Thompson 150ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 151d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &level_dms); CHKERRQ(ierr); 152d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 153d1d35e2fSjeremylt ierr = DMClone(dm_orig, &level_dms[level]); CHKERRQ(ierr); 154d1d35e2fSjeremylt ierr = DMGetVecType(dm_orig, &vectype); CHKERRQ(ierr); 155d1d35e2fSjeremylt ierr = DMSetVecType(level_dms[level], vectype); CHKERRQ(ierr); 156d1d35e2fSjeremylt ierr = SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level], 157d1d35e2fSjeremylt PETSC_TRUE, num_comp_u); CHKERRQ(ierr); 158a3c02c40SJeremy L Thompson // -- Label field components for viewing 159a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 160a3c02c40SJeremy L Thompson PetscSection section; 161d1d35e2fSjeremylt ierr = DMGetLocalSection(level_dms[level], §ion); CHKERRQ(ierr); 162a3c02c40SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 163a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 164a3c02c40SJeremy L Thompson CHKERRQ(ierr); 165a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 166a3c02c40SJeremy L Thompson CHKERRQ(ierr); 167a3c02c40SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 168a3c02c40SJeremy L Thompson CHKERRQ(ierr); 169a3c02c40SJeremy L Thompson } 170a3c02c40SJeremy L Thompson 171a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 172d1d35e2fSjeremylt ierr = DMClone(dm_orig, &dm_energy); CHKERRQ(ierr); 173d1d35e2fSjeremylt ierr = SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level], 174d1d35e2fSjeremylt PETSC_FALSE, num_comp_e); CHKERRQ(ierr); 175d1d35e2fSjeremylt ierr = DMClone(dm_orig, &dm_diagnostic); CHKERRQ(ierr); 176d1d35e2fSjeremylt ierr = SetupDMByDegree(dm_diagnostic, app_ctx, 177d1d35e2fSjeremylt app_ctx->level_degrees[fine_level], 178d1d35e2fSjeremylt PETSC_FALSE, num_comp_u + num_comp_d); CHKERRQ(ierr); 179d1d35e2fSjeremylt ierr = DMSetVecType(dm_energy, vectype); CHKERRQ(ierr); 180d1d35e2fSjeremylt ierr = DMSetVecType(dm_diagnostic, 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; 185d1d35e2fSjeremylt ierr = DMGetLocalSection(dm_diagnostic, §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 209d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &U_g); CHKERRQ(ierr); 210d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &U_loc); CHKERRQ(ierr); 211d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &U_g_size); CHKERRQ(ierr); 212d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &U_l_size); CHKERRQ(ierr); 213d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &U_loc_size); CHKERRQ(ierr); 214ccaff030SJeremy L Thompson 215ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 216d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 217ccaff030SJeremy L Thompson // -- Create global unknown vector U 218d1d35e2fSjeremylt ierr = DMCreateGlobalVector(level_dms[level], &U_g[level]); CHKERRQ(ierr); 219d1d35e2fSjeremylt ierr = VecGetSize(U_g[level], &U_g_size[level]); CHKERRQ(ierr); 220ccaff030SJeremy L Thompson // Note: Local size for matShell 221d1d35e2fSjeremylt ierr = VecGetLocalSize(U_g[level], &U_l_size[level]); CHKERRQ(ierr); 222ccaff030SJeremy L Thompson 223d1d35e2fSjeremylt // -- Create local unknown vector U_loc 224d1d35e2fSjeremylt ierr = DMCreateLocalVector(level_dms[level], &U_loc[level]); CHKERRQ(ierr); 225ccaff030SJeremy L Thompson // Note: local size for libCEED 226d1d35e2fSjeremylt ierr = VecGetSize(U_loc[level], &U_loc_size[level]); CHKERRQ(ierr); 227ccaff030SJeremy L Thompson } 228ccaff030SJeremy L Thompson 229ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 230d1d35e2fSjeremylt ierr = VecDuplicate(U_g[fine_level], &U); CHKERRQ(ierr); 231d1d35e2fSjeremylt ierr = VecDuplicate(U_g[fine_level], &R); CHKERRQ(ierr); 232d1d35e2fSjeremylt ierr = VecDuplicate(U_g[fine_level], &F); CHKERRQ(ierr); 233d1d35e2fSjeremylt ierr = VecDuplicate(U_loc[fine_level], &R_loc); CHKERRQ(ierr); 234d1d35e2fSjeremylt ierr = VecDuplicate(U_loc[fine_level], &F_loc); 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 243d1d35e2fSjeremylt ierr = PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup); 244ccaff030SJeremy L Thompson CHKERRQ(ierr); 245d1d35e2fSjeremylt ierr = PetscLogStagePush(stage_libceed_setup); CHKERRQ(ierr); 246ccaff030SJeremy L Thompson 247ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 248d1d35e2fSjeremylt CeedVector force_ceed; 249ccaff030SJeremy L Thompson CeedScalar *f; 250d1d35e2fSjeremylt PetscMemType force_mem_type; 251d1d35e2fSjeremylt if (app_ctx->forcing_choice != FORCE_NONE) { 252d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(F_loc, &f, &force_mem_type); CHKERRQ(ierr); 253d1d35e2fSjeremylt CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed); 254d1d35e2fSjeremylt CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f); 255ccaff030SJeremy L Thompson } 256ccaff030SJeremy L Thompson 257fe394131Sjeremylt // -- Create libCEED local Neumann BCs vector 258d1d35e2fSjeremylt CeedVector neumann_ceed; 259fe394131Sjeremylt CeedScalar *n; 260d1d35e2fSjeremylt PetscMemType nummann_mem_type; 261d1d35e2fSjeremylt if (app_ctx->bc_traction_count > 0) { 262d1d35e2fSjeremylt ierr = VecDuplicate(U, &neumann_bcs); CHKERRQ(ierr); 263d1d35e2fSjeremylt ierr = VecDuplicate(U_loc[fine_level], &bcs_loc); CHKERRQ(ierr); 264d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type); CHKERRQ(ierr); 265d1d35e2fSjeremylt CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed); 266d1d35e2fSjeremylt CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type), 267fe394131Sjeremylt CEED_USE_POINTER, n); 268fe394131Sjeremylt } 269fe394131Sjeremylt 270ccaff030SJeremy L Thompson // -- Setup libCEED objects 271d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr); 272d99fa3c5SJeremy L Thompson // ---- Setup residual, Jacobian evaluator and geometric information 273d1d35e2fSjeremylt ierr = PetscCalloc1(1, &ceed_data[fine_level]); CHKERRQ(ierr); 274ccaff030SJeremy L Thompson { 275d1d35e2fSjeremylt ierr = SetupLibceedFineLevel(level_dms[fine_level], dm_energy, dm_diagnostic, 276d1d35e2fSjeremylt ceed, app_ctx, ctx_phys, ceed_data, fine_level, 277d1d35e2fSjeremylt num_comp_u, U_g_size[fine_level], U_loc_size[fine_level], 278d1d35e2fSjeremylt force_ceed, neumann_ceed); 279ccaff030SJeremy L Thompson CHKERRQ(ierr); 280ccaff030SJeremy L Thompson } 281d99fa3c5SJeremy L Thompson // ---- Setup coarse Jacobian evaluator and prolongation/restriction 282d1d35e2fSjeremylt for (PetscInt level = num_levels - 2; level >= 0; level--) { 283d1d35e2fSjeremylt ierr = PetscCalloc1(1, &ceed_data[level]); CHKERRQ(ierr); 284d99fa3c5SJeremy L Thompson 285d99fa3c5SJeremy L Thompson // Get global communication restriction 286d1d35e2fSjeremylt ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr); 287d1d35e2fSjeremylt ierr = VecSet(U_loc[level+1], 1.0); CHKERRQ(ierr); 288d1d35e2fSjeremylt ierr = DMLocalToGlobal(level_dms[level+1], U_loc[level+1], ADD_VALUES, 289d1d35e2fSjeremylt U_g[level+1]); CHKERRQ(ierr); 290d1d35e2fSjeremylt ierr = DMGlobalToLocal(level_dms[level+1], U_g[level+1], INSERT_VALUES, 291d1d35e2fSjeremylt U_loc[level+1]); CHKERRQ(ierr); 292d99fa3c5SJeremy L Thompson 293d99fa3c5SJeremy L Thompson // Place in libCEED array 294d99fa3c5SJeremy L Thompson const PetscScalar *m; 295d1d35e2fSjeremylt PetscMemType m_mem_type; 296d1d35e2fSjeremylt ierr = VecGetArrayReadAndMemType(U_loc[level+1], &m, &m_mem_type); 297d1d35e2fSjeremylt CHKERRQ(ierr); 298d1d35e2fSjeremylt CeedVectorSetArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type), 299d99fa3c5SJeremy L Thompson CEED_USE_POINTER, (CeedScalar *)m); 300ccaff030SJeremy L Thompson 3011c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 302d1d35e2fSjeremylt ierr = SetupLibceedLevel(level_dms[level], ceed, app_ctx, 303d1d35e2fSjeremylt ceed_data, level, num_comp_u, U_g_size[level], 304d1d35e2fSjeremylt U_loc_size[level], ceed_data[level+1]->x_ceed); 305777ff853SJeremy L Thompson CHKERRQ(ierr); 306d99fa3c5SJeremy L Thompson 307d99fa3c5SJeremy L Thompson // Restore PETSc vector 308d1d35e2fSjeremylt CeedVectorTakeArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type), 309d99fa3c5SJeremy L Thompson (CeedScalar **)&m); 310d1d35e2fSjeremylt ierr = VecRestoreArrayReadAndMemType(U_loc[level+1], &m); CHKERRQ(ierr); 311d1d35e2fSjeremylt ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr); 312d1d35e2fSjeremylt ierr = VecZeroEntries(U_loc[level+1]); CHKERRQ(ierr); 313ccaff030SJeremy L Thompson } 314ccaff030SJeremy L Thompson 315ccaff030SJeremy L Thompson // Performance logging 316ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 317ccaff030SJeremy L Thompson 318ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 319fe394131Sjeremylt // Setup global forcing and Neumann BC vectors 320ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 321ccaff030SJeremy L Thompson ierr = VecZeroEntries(F); CHKERRQ(ierr); 322ccaff030SJeremy L Thompson 323d1d35e2fSjeremylt if (app_ctx->forcing_choice != FORCE_NONE) { 324d1d35e2fSjeremylt CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL); 325d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(F_loc, &f); CHKERRQ(ierr); 326d1d35e2fSjeremylt ierr = DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F); 327ccaff030SJeremy L Thompson CHKERRQ(ierr); 328d1d35e2fSjeremylt CeedVectorDestroy(&force_ceed); 329ccaff030SJeremy L Thompson } 330ccaff030SJeremy L Thompson 331d1d35e2fSjeremylt if (app_ctx->bc_traction_count > 0) { 332d1d35e2fSjeremylt ierr = VecZeroEntries(neumann_bcs); CHKERRQ(ierr); 333d1d35e2fSjeremylt CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL); 334d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(bcs_loc, &n); CHKERRQ(ierr); 335d1d35e2fSjeremylt ierr = DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs); 336fe394131Sjeremylt CHKERRQ(ierr); 337d1d35e2fSjeremylt CeedVectorDestroy(&neumann_ceed); 338fe394131Sjeremylt } 339fe394131Sjeremylt 340ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 341ccaff030SJeremy L Thompson // Print problem summary 342ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 343d1d35e2fSjeremylt if (!app_ctx->test_mode) { 344ccaff030SJeremy L Thompson const char *usedresource; 345ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 346ccaff030SJeremy L Thompson 347ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 34817fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 349ccaff030SJeremy L Thompson " libCEED:\n" 35062e9c006SJeremy L Thompson " libCEED Backend : %s\n" 351b68a8d79SJed Brown " libCEED Backend MemType : %s\n", 352d1d35e2fSjeremylt usedresource, CeedMemTypes[mem_type_backend]); 35362e9c006SJeremy L Thompson CHKERRQ(ierr); 354ccaff030SJeremy L Thompson 35562e9c006SJeremy L Thompson VecType vecType; 35662e9c006SJeremy L Thompson ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 35762e9c006SJeremy L Thompson ierr = PetscPrintf(comm, 35862e9c006SJeremy L Thompson " PETSc:\n" 35962e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 36062e9c006SJeremy L Thompson vecType); CHKERRQ(ierr); 36162e9c006SJeremy L Thompson 362ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 363ccaff030SJeremy L Thompson " Problem:\n" 364ccaff030SJeremy L Thompson " Problem Name : %s\n" 365ccaff030SJeremy L Thompson " Forcing Function : %s\n" 366ccaff030SJeremy L Thompson " Mesh:\n" 367ccaff030SJeremy L Thompson " File : %s\n" 368ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 369ccaff030SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 370ccaff030SJeremy L Thompson " Global nodes : %D\n" 371ccaff030SJeremy L Thompson " Owned nodes : %D\n" 372ccaff030SJeremy L Thompson " DoF per node : %D\n" 373ccaff030SJeremy L Thompson " Multigrid:\n" 374ccaff030SJeremy L Thompson " Type : %s\n" 375ccaff030SJeremy L Thompson " Number of Levels : %d\n", 376d1d35e2fSjeremylt problemTypesForDisp[app_ctx->problem_choice], 377d1d35e2fSjeremylt forcing_types_for_disp[app_ctx->forcing_choice], 378d1d35e2fSjeremylt app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh", 379d1d35e2fSjeremylt app_ctx->degree + 1, app_ctx->degree + 1, 380d1d35e2fSjeremylt U_g_size[fine_level]/num_comp_u, U_l_size[fine_level]/num_comp_u, num_comp_u, 381d1d35e2fSjeremylt (app_ctx->degree == 1 && 382d1d35e2fSjeremylt app_ctx->multigrid_choice != MULTIGRID_NONE) ? 383f892e6d1Sjeremylt "Algebraic multigrid" : 384d1d35e2fSjeremylt multigrid_types_for_disp[app_ctx->multigrid_choice], 385d1d35e2fSjeremylt (app_ctx->degree == 1 || 386d1d35e2fSjeremylt app_ctx->multigrid_choice == MULTIGRID_NONE) ? 387d1d35e2fSjeremylt 0 : num_levels); CHKERRQ(ierr); 388ccaff030SJeremy L Thompson 389d1d35e2fSjeremylt if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 390e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 391d1d35e2fSjeremylt CeedInt level = i ? fine_level : 0; 392ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 393ccaff030SJeremy L Thompson " Level %D (%s):\n" 394ccaff030SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 395ccaff030SJeremy L Thompson " Global Nodes : %D\n" 396ccaff030SJeremy L Thompson " Owned Nodes : %D\n", 397ccaff030SJeremy L Thompson level, i ? "fine" : "coarse", 398d1d35e2fSjeremylt app_ctx->level_degrees[level] + 1, 399d1d35e2fSjeremylt U_g_size[level]/num_comp_u, U_l_size[level]/num_comp_u); 400ccaff030SJeremy L Thompson CHKERRQ(ierr); 401ccaff030SJeremy L Thompson } 402ccaff030SJeremy L Thompson } 403ccaff030SJeremy L Thompson } 404ccaff030SJeremy L Thompson 405ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 406ccaff030SJeremy L Thompson // Setup SNES 407ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 408ccaff030SJeremy L Thompson // Performance logging 409d1d35e2fSjeremylt ierr = PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup); 410ccaff030SJeremy L Thompson CHKERRQ(ierr); 411d1d35e2fSjeremylt ierr = PetscLogStagePush(stage_snes_setup); CHKERRQ(ierr); 412ccaff030SJeremy L Thompson 413ccaff030SJeremy L Thompson // Create SNES 414ccaff030SJeremy L Thompson ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 415d1d35e2fSjeremylt ierr = SNESSetDM(snes, level_dms[fine_level]); CHKERRQ(ierr); 416ccaff030SJeremy L Thompson 417ccaff030SJeremy L Thompson // -- Jacobian evaluators 418d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &jacob_ctx); CHKERRQ(ierr); 419d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &jacob_mat); CHKERRQ(ierr); 420d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 421ccaff030SJeremy L Thompson // -- Jacobian context for level 422d1d35e2fSjeremylt ierr = PetscMalloc1(1, &jacob_ctx[level]); CHKERRQ(ierr); 423d1d35e2fSjeremylt ierr = SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level], 424d1d35e2fSjeremylt U_loc[level], ceed_data[level], ceed, ctx_phys, 425d1d35e2fSjeremylt ctx_phys_smoother, jacob_ctx[level]); CHKERRQ(ierr); 426ccaff030SJeremy L Thompson 427ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 428d1d35e2fSjeremylt ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level], 429d1d35e2fSjeremylt U_g_size[level], jacob_ctx[level], &jacob_mat[level]); 430ccaff030SJeremy L Thompson CHKERRQ(ierr); 431d1d35e2fSjeremylt ierr = MatShellSetOperation(jacob_mat[level], MATOP_MULT, 432ccaff030SJeremy L Thompson (void (*)(void))ApplyJacobian_Ceed); 433ccaff030SJeremy L Thompson CHKERRQ(ierr); 434d1d35e2fSjeremylt ierr = MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL, 435ccaff030SJeremy L Thompson (void(*)(void))GetDiag_Ceed); 436d1d35e2fSjeremylt ierr = MatShellSetVecType(jacob_mat[level], vectype); CHKERRQ(ierr); 437ccaff030SJeremy L Thompson } 438e3e3df41Sjeremylt // Note: FormJacobian updates Jacobian matrices on each level 439ccaff030SJeremy L Thompson // and assembles the Jpre matrix, if needed 440d1d35e2fSjeremylt ierr = PetscMalloc1(1, &form_jacob_ctx); CHKERRQ(ierr); 441d1d35e2fSjeremylt form_jacob_ctx->jacob_ctx = jacob_ctx; 442d1d35e2fSjeremylt form_jacob_ctx->num_levels = num_levels; 443d1d35e2fSjeremylt form_jacob_ctx->jacob_mat = jacob_mat; 444ccaff030SJeremy L Thompson 445ccaff030SJeremy L Thompson // -- Residual evaluation function 446d1d35e2fSjeremylt ierr = PetscCalloc1(1, &res_ctx); CHKERRQ(ierr); 447d1d35e2fSjeremylt ierr = PetscMemcpy(res_ctx, jacob_ctx[fine_level], 448d1d35e2fSjeremylt sizeof(*jacob_ctx[fine_level])); CHKERRQ(ierr); 449d1d35e2fSjeremylt res_ctx->op = ceed_data[fine_level]->op_apply; 450d1d35e2fSjeremylt res_ctx->qf = ceed_data[fine_level]->qf_apply; 451d1d35e2fSjeremylt if (app_ctx->bc_traction_count > 0) 452d1d35e2fSjeremylt res_ctx->neumann_bcs = neumann_bcs; 453fe394131Sjeremylt else 454d1d35e2fSjeremylt res_ctx->neumann_bcs = NULL; 455d1d35e2fSjeremylt ierr = SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx); CHKERRQ(ierr); 456ccaff030SJeremy L Thompson 457ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 458d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &prolong_restr_ctx); CHKERRQ(ierr); 459d1d35e2fSjeremylt ierr = PetscMalloc1(num_levels, &prolong_restr_mat); CHKERRQ(ierr); 460d1d35e2fSjeremylt for (PetscInt level = 1; level < num_levels; level++) { 461ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 462d1d35e2fSjeremylt ierr = PetscMalloc1(1, &prolong_restr_ctx[level]); CHKERRQ(ierr); 463d1d35e2fSjeremylt ierr = SetupProlongRestrictCtx(comm, app_ctx, level_dms[level-1], 464d1d35e2fSjeremylt level_dms[level], U_g[level], U_loc[level-1], 465d1d35e2fSjeremylt U_loc[level], ceed_data[level-1], 466d1d35e2fSjeremylt ceed_data[level], ceed, 467d1d35e2fSjeremylt prolong_restr_ctx[level]); CHKERRQ(ierr); 468ccaff030SJeremy L Thompson 469ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 470d1d35e2fSjeremylt ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level-1], U_g_size[level], 471d1d35e2fSjeremylt U_g_size[level-1], prolong_restr_ctx[level], 472d1d35e2fSjeremylt &prolong_restr_mat[level]); CHKERRQ(ierr); 473ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 474d1d35e2fSjeremylt ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT, 475ccaff030SJeremy L Thompson (void (*)(void))Prolong_Ceed); 476d1d35e2fSjeremylt ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE, 477ccaff030SJeremy L Thompson (void (*)(void))Restrict_Ceed); 478ccaff030SJeremy L Thompson CHKERRQ(ierr); 479d1d35e2fSjeremylt ierr = MatShellSetVecType(prolong_restr_mat[level], vectype); CHKERRQ(ierr); 480ccaff030SJeremy L Thompson } 481ccaff030SJeremy L Thompson 482ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 483ccaff030SJeremy L Thompson // Setup dummy SNES for AMG coarse solve 484ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 485d1d35e2fSjeremylt if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 486e3e3df41Sjeremylt // -- Jacobian Matrix 487d1d35e2fSjeremylt ierr = DMSetMatType(level_dms[0], MATAIJ); CHKERRQ(ierr); 488d1d35e2fSjeremylt ierr = DMCreateMatrix(level_dms[0], &jacob_mat_coarse); CHKERRQ(ierr); 489e3e3df41Sjeremylt 490d1d35e2fSjeremylt if (app_ctx->degree > 1) { 491d1d35e2fSjeremylt ierr = SNESCreate(comm, &snes_coarse); CHKERRQ(ierr); 492d1d35e2fSjeremylt ierr = SNESSetDM(snes_coarse, level_dms[0]); CHKERRQ(ierr); 493d1d35e2fSjeremylt ierr = SNESSetSolution(snes_coarse, U_g[0]); CHKERRQ(ierr); 494ccaff030SJeremy L Thompson 495e3e3df41Sjeremylt // -- Jacobian function 496d1d35e2fSjeremylt ierr = SNESSetJacobian(snes_coarse, jacob_mat_coarse, jacob_mat_coarse, NULL, 497ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 498ccaff030SJeremy L Thompson 499ccaff030SJeremy L Thompson // -- Residual evaluation function 500d1d35e2fSjeremylt ierr = PetscMalloc1(1, &jacob_coarse_ctx); CHKERRQ(ierr); 501d1d35e2fSjeremylt ierr = PetscMemcpy(jacob_coarse_ctx, jacob_ctx[0], sizeof(*jacob_ctx[0])); 502ccaff030SJeremy L Thompson CHKERRQ(ierr); 503d1d35e2fSjeremylt ierr = SNESSetFunction(snes_coarse, U_g[0], ApplyJacobianCoarse_Ceed, 504d1d35e2fSjeremylt jacob_coarse_ctx); CHKERRQ(ierr); 505ccaff030SJeremy L Thompson 506d1d35e2fSjeremylt // -- Update form_jacob_ctx 507d1d35e2fSjeremylt form_jacob_ctx->u_coarse = U_g[0]; 508d1d35e2fSjeremylt form_jacob_ctx->snes_coarse = snes_coarse; 509d1d35e2fSjeremylt form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse; 510e3e3df41Sjeremylt } 511e3e3df41Sjeremylt } 512e3e3df41Sjeremylt 513e3e3df41Sjeremylt // Set Jacobian function 514d1d35e2fSjeremylt if (app_ctx->degree > 1) { 515d1d35e2fSjeremylt ierr = SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level], 516d1d35e2fSjeremylt FormJacobian, form_jacob_ctx); CHKERRQ(ierr); 517e3e3df41Sjeremylt } else { 518d1d35e2fSjeremylt ierr = SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse, 519e3e3df41Sjeremylt SNESComputeJacobianDefaultColor, NULL); 520e3e3df41Sjeremylt CHKERRQ(ierr); 521e3e3df41Sjeremylt } 522ccaff030SJeremy L Thompson 523ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 524ccaff030SJeremy L Thompson // Setup KSP 525ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 526ccaff030SJeremy L Thompson { 527ccaff030SJeremy L Thompson PC pc; 528ccaff030SJeremy L Thompson KSP ksp; 529ccaff030SJeremy L Thompson 530ccaff030SJeremy L Thompson // -- KSP 531ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 532ccaff030SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 533ccaff030SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 534ccaff030SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 535ccaff030SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 536ccaff030SJeremy L Thompson ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 537ccaff030SJeremy L Thompson 538ccaff030SJeremy L Thompson // -- Preconditioning 539ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 540d1d35e2fSjeremylt ierr = PCSetDM(pc, level_dms[fine_level]); CHKERRQ(ierr); 541ccaff030SJeremy L Thompson ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 542ccaff030SJeremy L Thompson 543d1d35e2fSjeremylt if (app_ctx->multigrid_choice == MULTIGRID_NONE) { 544ccaff030SJeremy L Thompson // ---- No Multigrid 545ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 546ccaff030SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 547d1d35e2fSjeremylt } else if (app_ctx->degree == 1) { 548f892e6d1Sjeremylt // ---- AMG for degree 1 549f892e6d1Sjeremylt ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 550ccaff030SJeremy L Thompson } else { 551ccaff030SJeremy L Thompson // ---- PCMG 552ccaff030SJeremy L Thompson ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 553ccaff030SJeremy L Thompson 554ccaff030SJeremy L Thompson // ------ PCMG levels 555d1d35e2fSjeremylt ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr); 556d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 557ccaff030SJeremy L Thompson // -------- Smoother 558d1d35e2fSjeremylt KSP ksp_smoother, ksp_est; 559d1d35e2fSjeremylt PC pc_smoother; 560ccaff030SJeremy L Thompson 561ccaff030SJeremy L Thompson // ---------- Smoother KSP 562d1d35e2fSjeremylt ierr = PCMGGetSmoother(pc, level, &ksp_smoother); CHKERRQ(ierr); 563d1d35e2fSjeremylt ierr = KSPSetDM(ksp_smoother, level_dms[level]); CHKERRQ(ierr); 564d1d35e2fSjeremylt ierr = KSPSetDMActive(ksp_smoother, PETSC_FALSE); CHKERRQ(ierr); 565ccaff030SJeremy L Thompson 566ccaff030SJeremy L Thompson // ---------- Chebyshev options 567d1d35e2fSjeremylt ierr = KSPSetType(ksp_smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 568d1d35e2fSjeremylt ierr = KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1); 569ccaff030SJeremy L Thompson CHKERRQ(ierr); 570d1d35e2fSjeremylt ierr = KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est); CHKERRQ(ierr); 571d1d35e2fSjeremylt ierr = KSPSetType(ksp_est, KSPCG); CHKERRQ(ierr); 572d1d35e2fSjeremylt ierr = KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE); 573ccaff030SJeremy L Thompson CHKERRQ(ierr); 574d1d35e2fSjeremylt ierr = KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]); 575ccaff030SJeremy L Thompson CHKERRQ(ierr); 576ccaff030SJeremy L Thompson 577ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 578d1d35e2fSjeremylt ierr = KSPGetPC(ksp_smoother, &pc_smoother); CHKERRQ(ierr); 579d1d35e2fSjeremylt ierr = PCSetType(pc_smoother, PCJACOBI); CHKERRQ(ierr); 580d1d35e2fSjeremylt ierr = PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 581ccaff030SJeremy L Thompson 582ccaff030SJeremy L Thompson // -------- Work vector 583d1d35e2fSjeremylt if (level != fine_level) { 584d1d35e2fSjeremylt ierr = PCMGSetX(pc, level, U_g[level]); CHKERRQ(ierr); 585ccaff030SJeremy L Thompson } 586ccaff030SJeremy L Thompson 587ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 588ccaff030SJeremy L Thompson if (level > 0) { 589d1d35e2fSjeremylt ierr = PCMGSetInterpolation(pc, level, prolong_restr_mat[level]); 590ccaff030SJeremy L Thompson CHKERRQ(ierr); 591d1d35e2fSjeremylt ierr = PCMGSetRestriction(pc, level, prolong_restr_mat[level]); 592ccaff030SJeremy L Thompson CHKERRQ(ierr); 593ccaff030SJeremy L Thompson } 594ccaff030SJeremy L Thompson } 595ccaff030SJeremy L Thompson 596ccaff030SJeremy L Thompson // ------ PCMG coarse solve 597d1d35e2fSjeremylt KSP ksp_coarse; 598d1d35e2fSjeremylt PC pc_coarse; 599ccaff030SJeremy L Thompson 600ccaff030SJeremy L Thompson // -------- Coarse KSP 601d1d35e2fSjeremylt ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr); 602d1d35e2fSjeremylt ierr = KSPSetType(ksp_coarse, KSPPREONLY); CHKERRQ(ierr); 603d1d35e2fSjeremylt ierr = KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse); 604ccaff030SJeremy L Thompson CHKERRQ(ierr); 605d1d35e2fSjeremylt ierr = KSPSetOptionsPrefix(ksp_coarse, "coarse_"); CHKERRQ(ierr); 606ccaff030SJeremy L Thompson 607ccaff030SJeremy L Thompson // -------- Coarse preconditioner 608d1d35e2fSjeremylt ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr); 609d1d35e2fSjeremylt ierr = PCSetType(pc_coarse, PCGAMG); CHKERRQ(ierr); 610d1d35e2fSjeremylt ierr = PCSetOptionsPrefix(pc_coarse, "coarse_"); CHKERRQ(ierr); 611ccaff030SJeremy L Thompson 612d1d35e2fSjeremylt ierr = KSPSetFromOptions(ksp_coarse); CHKERRQ(ierr); 613d1d35e2fSjeremylt ierr = PCSetFromOptions(pc_coarse); CHKERRQ(ierr); 614ccaff030SJeremy L Thompson 615ccaff030SJeremy L Thompson // ------ PCMG options 616ccaff030SJeremy L Thompson ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 617ccaff030SJeremy L Thompson ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 618d1d35e2fSjeremylt ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr); 619ccaff030SJeremy L Thompson } 620ccaff030SJeremy L Thompson ierr = KSPSetFromOptions(ksp); 621ccaff030SJeremy L Thompson ierr = PCSetFromOptions(pc); 622ccaff030SJeremy L Thompson } 623038d0cd7Sjeremylt { 624038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 625d1d35e2fSjeremylt SNESLineSearch line_search; 626481a4840SJed Brown 627d1d35e2fSjeremylt ierr = SNESGetLineSearch(snes, &line_search); CHKERRQ(ierr); 628d1d35e2fSjeremylt ierr = SNESLineSearchSetType(line_search, SNESLINESEARCHCP); CHKERRQ(ierr); 629481a4840SJed Brown } 630481a4840SJed Brown 631ccaff030SJeremy L Thompson ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 632ccaff030SJeremy L Thompson 633ccaff030SJeremy L Thompson // Performance logging 634ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 635ccaff030SJeremy L Thompson 636ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 637ccaff030SJeremy L Thompson // Set initial guess 638ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 639f81c27eaSJed Brown ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 640ccaff030SJeremy L Thompson ierr = VecSet(U, 0.0); CHKERRQ(ierr); 641ccaff030SJeremy L Thompson 642ccaff030SJeremy L Thompson // View solution 643d1d35e2fSjeremylt if (app_ctx->view_soln) { 644d1d35e2fSjeremylt ierr = ViewSolution(comm, app_ctx, U, 0, 0.0); CHKERRQ(ierr); 645ccaff030SJeremy L Thompson } 646ccaff030SJeremy L Thompson 647ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 648ccaff030SJeremy L Thompson // Solve SNES 649ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 650d1d35e2fSjeremylt PetscBool snes_monitor = PETSC_FALSE; 651d1d35e2fSjeremylt ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor); 6525f24f028Sjeremylt CHKERRQ(ierr); 6535f24f028Sjeremylt 654ccaff030SJeremy L Thompson // Performance logging 655d1d35e2fSjeremylt ierr = PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve); 656ccaff030SJeremy L Thompson CHKERRQ(ierr); 657d1d35e2fSjeremylt ierr = PetscLogStagePush(stage_snes_solve); CHKERRQ(ierr); 658ccaff030SJeremy L Thompson 659ccaff030SJeremy L Thompson // Timing 660ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 661d1d35e2fSjeremylt start_time = MPI_Wtime(); 662ccaff030SJeremy L Thompson 663ccaff030SJeremy L Thompson // Solve for each load increment 6645f24f028Sjeremylt PetscInt increment; 665d1d35e2fSjeremylt for (increment = 1; increment <= app_ctx->num_increments; increment++) { 6665f24f028Sjeremylt // -- Log increment count 667d1d35e2fSjeremylt if (snes_monitor) { 6685f24f028Sjeremylt ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 6695f24f028Sjeremylt CHKERRQ(ierr); 6705f24f028Sjeremylt } 6715f24f028Sjeremylt 672ccaff030SJeremy L Thompson // -- Scale the problem 673d1d35e2fSjeremylt PetscScalar load_increment = 1.0*increment / app_ctx->num_increments, 674d1d35e2fSjeremylt scalingFactor = load_increment / 675d1d35e2fSjeremylt (increment == 1 ? 1 : res_ctx->load_increment); 676d1d35e2fSjeremylt res_ctx->load_increment = load_increment; 677d1d35e2fSjeremylt if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) { 678ccaff030SJeremy L Thompson ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 679ccaff030SJeremy L Thompson } 680ccaff030SJeremy L Thompson 681ccaff030SJeremy L Thompson // -- Solve 682ccaff030SJeremy L Thompson ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 683ccaff030SJeremy L Thompson 684ccaff030SJeremy L Thompson // -- View solution 685d1d35e2fSjeremylt if (app_ctx->view_soln) { 686d1d35e2fSjeremylt ierr = ViewSolution(comm, app_ctx, U, increment, load_increment); CHKERRQ(ierr); 687ccaff030SJeremy L Thompson } 688ccaff030SJeremy L Thompson 689ccaff030SJeremy L Thompson // -- Update SNES iteration count 690ccaff030SJeremy L Thompson PetscInt its; 691ccaff030SJeremy L Thompson ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 692*bf0f51feSjeremylt snes_its += its; 6937418ca88SJeremy L Thompson ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr); 694*bf0f51feSjeremylt ksp_its += its; 6953951731eSjeremylt 6963951731eSjeremylt // -- Check for divergence 6973951731eSjeremylt SNESConvergedReason reason; 6983951731eSjeremylt ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 6993951731eSjeremylt if (reason < 0) 7003951731eSjeremylt break; 701ccaff030SJeremy L Thompson } 702ccaff030SJeremy L Thompson 703ccaff030SJeremy L Thompson // Timing 704d1d35e2fSjeremylt elapsed_time = MPI_Wtime() - start_time; 705ccaff030SJeremy L Thompson 706ccaff030SJeremy L Thompson // Performance logging 707ccaff030SJeremy L Thompson ierr = PetscLogStagePop(); 708ccaff030SJeremy L Thompson 709ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 710ccaff030SJeremy L Thompson // Output summary 711ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 712d1d35e2fSjeremylt if (!app_ctx->test_mode) { 713ccaff030SJeremy L Thompson // -- SNES 714d1d35e2fSjeremylt SNESType snes_type; 715ccaff030SJeremy L Thompson SNESConvergedReason reason; 716ccaff030SJeremy L Thompson PetscReal rnorm; 717d1d35e2fSjeremylt ierr = SNESGetType(snes, &snes_type); CHKERRQ(ierr); 718ccaff030SJeremy L Thompson ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 719ccaff030SJeremy L Thompson ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 720ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 721ccaff030SJeremy L Thompson " SNES:\n" 722ccaff030SJeremy L Thompson " SNES Type : %s\n" 723ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 724ccaff030SJeremy L Thompson " Number of Load Increments : %d\n" 7255f24f028Sjeremylt " Completed Load Increments : %d\n" 726ccaff030SJeremy L Thompson " Total SNES Iterations : %D\n" 727ccaff030SJeremy L Thompson " Final rnorm : %e\n", 728d1d35e2fSjeremylt snes_type, SNESConvergedReasons[reason], 729d1d35e2fSjeremylt app_ctx->num_increments, increment - 1, 730*bf0f51feSjeremylt snes_its, (double)rnorm); CHKERRQ(ierr); 731ccaff030SJeremy L Thompson 732ccaff030SJeremy L Thompson // -- KSP 733ccaff030SJeremy L Thompson KSP ksp; 734d1d35e2fSjeremylt KSPType ksp_type; 735ccaff030SJeremy L Thompson ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 736d1d35e2fSjeremylt ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 737ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 738ccaff030SJeremy L Thompson " Linear Solver:\n" 7397418ca88SJeremy L Thompson " KSP Type : %s\n" 7407418ca88SJeremy L Thompson " Total KSP Iterations : %D\n", 741*bf0f51feSjeremylt ksp_type, ksp_its); CHKERRQ(ierr); 742ccaff030SJeremy L Thompson 743ccaff030SJeremy L Thompson // -- PC 744ccaff030SJeremy L Thompson PC pc; 745d1d35e2fSjeremylt PCType pc_type; 746ccaff030SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 747d1d35e2fSjeremylt ierr = PCGetType(pc, &pc_type); CHKERRQ(ierr); 748e3e3df41Sjeremylt ierr = PetscPrintf(comm, 749e3e3df41Sjeremylt " PC Type : %s\n", 750d1d35e2fSjeremylt pc_type); CHKERRQ(ierr); 751e3e3df41Sjeremylt 752d1d35e2fSjeremylt if (!strcmp(pc_type, PCMG)) { 753d1d35e2fSjeremylt PCMGType pcmg_type; 754d1d35e2fSjeremylt ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr); 755ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 756ccaff030SJeremy L Thompson " P-Multigrid:\n" 757ccaff030SJeremy L Thompson " PCMG Type : %s\n" 758ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 759d1d35e2fSjeremylt PCMGTypes[pcmg_type], 760d1d35e2fSjeremylt PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr); 761ccaff030SJeremy L Thompson 762ccaff030SJeremy L Thompson // -- Coarse Solve 763d1d35e2fSjeremylt KSP ksp_coarse; 764d1d35e2fSjeremylt PC pc_coarse; 765d1d35e2fSjeremylt PCType pc_type; 766ccaff030SJeremy L Thompson 767d1d35e2fSjeremylt ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr); 768d1d35e2fSjeremylt ierr = KSPGetType(ksp_coarse, &ksp_type); CHKERRQ(ierr); 769d1d35e2fSjeremylt ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr); 770d1d35e2fSjeremylt ierr = PCGetType(pc_coarse, &pc_type); CHKERRQ(ierr); 771ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 772ccaff030SJeremy L Thompson " Coarse Solve:\n" 773ccaff030SJeremy L Thompson " KSP Type : %s\n" 774ccaff030SJeremy L Thompson " PC Type : %s\n", 775d1d35e2fSjeremylt ksp_type, pc_type); CHKERRQ(ierr); 776ccaff030SJeremy L Thompson } 777ccaff030SJeremy L Thompson } 778ccaff030SJeremy L Thompson 779ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 780ccaff030SJeremy L Thompson // Compute solve time 781ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 782d1d35e2fSjeremylt if (!app_ctx->test_mode) { 783d1d35e2fSjeremylt ierr = MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm); 784ccaff030SJeremy L Thompson CHKERRQ(ierr); 785d1d35e2fSjeremylt ierr = MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm); 786ccaff030SJeremy L Thompson CHKERRQ(ierr); 787ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 788ccaff030SJeremy L Thompson " Performance:\n" 7897418ca88SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n" 7907418ca88SJeremy L Thompson " DoFs/Sec in SNES : %g (%g) million\n", 791*bf0f51feSjeremylt max_time, min_time, 1e-6*U_g_size[fine_level]*ksp_its/max_time, 792*bf0f51feSjeremylt 1e-6*U_g_size[fine_level]*ksp_its/min_time); CHKERRQ(ierr); 793ccaff030SJeremy L Thompson } 794ccaff030SJeremy L Thompson 795ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 796ccaff030SJeremy L Thompson // Compute error 797ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 798d1d35e2fSjeremylt if (app_ctx->forcing_choice == FORCE_MMS) { 799d1d35e2fSjeremylt CeedScalar l2_error = 1., l2_U_norm = 1.; 800d1d35e2fSjeremylt const CeedScalar *true_array; 801d1d35e2fSjeremylt Vec error_vec, true_vec; 802ccaff030SJeremy L Thompson 803ccaff030SJeremy L Thompson // -- Work vectors 804d1d35e2fSjeremylt ierr = VecDuplicate(U, &error_vec); CHKERRQ(ierr); 805d1d35e2fSjeremylt ierr = VecSet(error_vec, 0.0); CHKERRQ(ierr); 806d1d35e2fSjeremylt ierr = VecDuplicate(U, &true_vec); CHKERRQ(ierr); 807d1d35e2fSjeremylt ierr = VecSet(true_vec, 0.0); CHKERRQ(ierr); 808ccaff030SJeremy L Thompson 809ccaff030SJeremy L Thompson // -- Assemble global true solution vector 810d1d35e2fSjeremylt CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln, 811d1d35e2fSjeremylt CEED_MEM_HOST, &true_array); 812d1d35e2fSjeremylt ierr = VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array); 81362e9c006SJeremy L Thompson CHKERRQ(ierr); 814d1d35e2fSjeremylt ierr = DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec); 815ccaff030SJeremy L Thompson CHKERRQ(ierr); 816d1d35e2fSjeremylt ierr = VecResetArray(res_ctx->Y_loc); CHKERRQ(ierr); 817d1d35e2fSjeremylt CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array); 818ccaff030SJeremy L Thompson 819ccaff030SJeremy L Thompson // -- Compute L2 error 820d1d35e2fSjeremylt ierr = VecWAXPY(error_vec, -1.0, U, true_vec); CHKERRQ(ierr); 821d1d35e2fSjeremylt ierr = VecNorm(error_vec, NORM_2, &l2_error); CHKERRQ(ierr); 822d1d35e2fSjeremylt ierr = VecNorm(U, NORM_2, &l2_U_norm); CHKERRQ(ierr); 823d1d35e2fSjeremylt l2_error /= l2_U_norm; 824ccaff030SJeremy L Thompson 825ccaff030SJeremy L Thompson // -- Output 826d1d35e2fSjeremylt if (!app_ctx->test_mode || l2_error > 0.05) { 827ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, 828ccaff030SJeremy L Thompson " L2 Error : %e\n", 829d1d35e2fSjeremylt l2_error); CHKERRQ(ierr); 830ccaff030SJeremy L Thompson } 831ccaff030SJeremy L Thompson 832ccaff030SJeremy L Thompson // -- Cleanup 833d1d35e2fSjeremylt ierr = VecDestroy(&error_vec); CHKERRQ(ierr); 834d1d35e2fSjeremylt ierr = VecDestroy(&true_vec); CHKERRQ(ierr); 8352d93065eSjeremylt } 8362d93065eSjeremylt 8372d93065eSjeremylt // --------------------------------------------------------------------------- 8382d93065eSjeremylt // Compute energy 8392d93065eSjeremylt // --------------------------------------------------------------------------- 840d1d35e2fSjeremylt if (!app_ctx->test_mode) { 8412d93065eSjeremylt // -- Compute L2 error 8422d93065eSjeremylt CeedScalar energy; 843d1d35e2fSjeremylt ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, 844a3c02c40SJeremy L Thompson U, &energy); CHKERRQ(ierr); 8452d93065eSjeremylt 8462d93065eSjeremylt // -- Output 8472d93065eSjeremylt ierr = PetscPrintf(comm, 84872d03b64SArash Mehraban " Strain Energy : %.12e\n", 8492d93065eSjeremylt energy); CHKERRQ(ierr); 850ccaff030SJeremy L Thompson } 851ccaff030SJeremy L Thompson 852ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 8535c25879aSJeremy L Thompson // Output diagnostic quantities 8545c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 855d1d35e2fSjeremylt if (app_ctx->view_soln || app_ctx->view_final_soln) { 8565c25879aSJeremy L Thompson // -- Setup context 857d1d35e2fSjeremylt UserMult diagnostic_ctx; 858d1d35e2fSjeremylt ierr = PetscMalloc1(1, &diagnostic_ctx); CHKERRQ(ierr); 859d1d35e2fSjeremylt ierr = PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)); CHKERRQ(ierr); 860d1d35e2fSjeremylt diagnostic_ctx->dm = dm_diagnostic; 861d1d35e2fSjeremylt diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic; 8625c25879aSJeremy L Thompson 8635c25879aSJeremy L Thompson // -- Compute and output 864d1d35e2fSjeremylt ierr = ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx, 865d1d35e2fSjeremylt app_ctx, U, 866d1d35e2fSjeremylt ceed_data[fine_level]->elem_restr_diagnostic); 8675c25879aSJeremy L Thompson CHKERRQ(ierr); 8685c25879aSJeremy L Thompson 8695c25879aSJeremy L Thompson // -- Cleanup 870d1d35e2fSjeremylt ierr = PetscFree(diagnostic_ctx); CHKERRQ(ierr); 8715c25879aSJeremy L Thompson } 8725c25879aSJeremy L Thompson 8735c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 874ccaff030SJeremy L Thompson // Free objects 875ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 876ccaff030SJeremy L Thompson // Data in arrays per level 877d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 878ccaff030SJeremy L Thompson // Vectors 879d1d35e2fSjeremylt ierr = VecDestroy(&U_g[level]); CHKERRQ(ierr); 880d1d35e2fSjeremylt ierr = VecDestroy(&U_loc[level]); CHKERRQ(ierr); 881ccaff030SJeremy L Thompson 882ccaff030SJeremy L Thompson // Jacobian matrix and data 883d1d35e2fSjeremylt ierr = VecDestroy(&jacob_ctx[level]->Y_loc); CHKERRQ(ierr); 884d1d35e2fSjeremylt ierr = MatDestroy(&jacob_mat[level]); CHKERRQ(ierr); 885d1d35e2fSjeremylt ierr = PetscFree(jacob_ctx[level]); CHKERRQ(ierr); 886ccaff030SJeremy L Thompson 887ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 888ccaff030SJeremy L Thompson if (level > 0) { 889d1d35e2fSjeremylt ierr = PetscFree(prolong_restr_ctx[level]); CHKERRQ(ierr); 890d1d35e2fSjeremylt ierr = MatDestroy(&prolong_restr_mat[level]); CHKERRQ(ierr); 891ccaff030SJeremy L Thompson } 892ccaff030SJeremy L Thompson 893ccaff030SJeremy L Thompson // DM 894d1d35e2fSjeremylt ierr = DMDestroy(&level_dms[level]); CHKERRQ(ierr); 895ccaff030SJeremy L Thompson 896ccaff030SJeremy L Thompson // libCEED objects 897d1d35e2fSjeremylt ierr = CeedDataDestroy(level, ceed_data[level]); CHKERRQ(ierr); 898ccaff030SJeremy L Thompson } 899ccaff030SJeremy L Thompson 900ccaff030SJeremy L Thompson // Arrays 901d1d35e2fSjeremylt ierr = PetscFree(U_g); CHKERRQ(ierr); 902d1d35e2fSjeremylt ierr = PetscFree(U_loc); CHKERRQ(ierr); 903d1d35e2fSjeremylt ierr = PetscFree(U_g_size); CHKERRQ(ierr); 904d1d35e2fSjeremylt ierr = PetscFree(U_l_size); CHKERRQ(ierr); 905d1d35e2fSjeremylt ierr = PetscFree(U_loc_size); CHKERRQ(ierr); 906d1d35e2fSjeremylt ierr = PetscFree(jacob_ctx); CHKERRQ(ierr); 907d1d35e2fSjeremylt ierr = PetscFree(jacob_mat); CHKERRQ(ierr); 908d1d35e2fSjeremylt ierr = PetscFree(prolong_restr_ctx); CHKERRQ(ierr); 909d1d35e2fSjeremylt ierr = PetscFree(prolong_restr_mat); CHKERRQ(ierr); 910d1d35e2fSjeremylt ierr = PetscFree(app_ctx->level_degrees); CHKERRQ(ierr); 911d1d35e2fSjeremylt ierr = PetscFree(ceed_data); CHKERRQ(ierr); 912ccaff030SJeremy L Thompson 913ccaff030SJeremy L Thompson // libCEED objects 914d1d35e2fSjeremylt CeedQFunctionContextDestroy(&ctx_phys); 915d1d35e2fSjeremylt CeedQFunctionContextDestroy(&ctx_phys_smoother); 916ccaff030SJeremy L Thompson CeedDestroy(&ceed); 917ccaff030SJeremy L Thompson 918ccaff030SJeremy L Thompson // PETSc objects 919ccaff030SJeremy L Thompson ierr = VecDestroy(&U); CHKERRQ(ierr); 920ccaff030SJeremy L Thompson ierr = VecDestroy(&R); CHKERRQ(ierr); 921d1d35e2fSjeremylt ierr = VecDestroy(&R_loc); CHKERRQ(ierr); 922ccaff030SJeremy L Thompson ierr = VecDestroy(&F); CHKERRQ(ierr); 923d1d35e2fSjeremylt ierr = VecDestroy(&F_loc); CHKERRQ(ierr); 924d1d35e2fSjeremylt ierr = VecDestroy(&neumann_bcs); CHKERRQ(ierr); 925d1d35e2fSjeremylt ierr = VecDestroy(&bcs_loc); CHKERRQ(ierr); 926d1d35e2fSjeremylt ierr = MatDestroy(&jacob_mat_coarse); CHKERRQ(ierr); 927ccaff030SJeremy L Thompson ierr = SNESDestroy(&snes); CHKERRQ(ierr); 928d1d35e2fSjeremylt ierr = SNESDestroy(&snes_coarse); CHKERRQ(ierr); 929d1d35e2fSjeremylt ierr = DMDestroy(&dm_orig); CHKERRQ(ierr); 930d1d35e2fSjeremylt ierr = DMDestroy(&dm_energy); CHKERRQ(ierr); 931d1d35e2fSjeremylt ierr = DMDestroy(&dm_diagnostic); CHKERRQ(ierr); 932d1d35e2fSjeremylt ierr = PetscFree(level_dms); CHKERRQ(ierr); 933ccaff030SJeremy L Thompson 934ccaff030SJeremy L Thompson // Structs 935d1d35e2fSjeremylt ierr = PetscFree(res_ctx); CHKERRQ(ierr); 936d1d35e2fSjeremylt ierr = PetscFree(form_jacob_ctx); CHKERRQ(ierr); 937d1d35e2fSjeremylt ierr = PetscFree(jacob_coarse_ctx); CHKERRQ(ierr); 938d1d35e2fSjeremylt ierr = PetscFree(app_ctx); CHKERRQ(ierr); 939ccaff030SJeremy L Thompson ierr = PetscFree(phys); CHKERRQ(ierr); 940d1d35e2fSjeremylt ierr = PetscFree(phys_smoother); CHKERRQ(ierr); 941ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 942ccaff030SJeremy L Thompson 943ccaff030SJeremy L Thompson return PetscFinalize(); 944ccaff030SJeremy L Thompson } 945