13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ccaff030SJeremy L Thompson // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5ccaff030SJeremy L Thompson // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7ccaff030SJeremy L Thompson 8ccaff030SJeremy L Thompson // libCEED + PETSc Example: Elasticity 9ccaff030SJeremy L Thompson // 10*ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve elasticity problems. 11ccaff030SJeremy L Thompson // 12ccaff030SJeremy L Thompson // The code uses higher level communication protocols in DMPlex. 13ccaff030SJeremy L Thompson // 14ccaff030SJeremy L Thompson // Build with: 15ccaff030SJeremy L Thompson // 16ccaff030SJeremy L Thompson // make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 17ccaff030SJeremy L Thompson // 18ccaff030SJeremy L Thompson // Sample runs: 19ccaff030SJeremy L Thompson // 20eccc2849SRezgar Shakeri // ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms 212b730f8bSJeremy L Thompson // ./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 222b730f8bSJeremy L Thompson // /cpu/self 232b730f8bSJeremy L Thompson // ./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 242b730f8bSJeremy L Thompson // -ceed /gpu/cuda 25ccaff030SJeremy L Thompson // 26ccaff030SJeremy L Thompson // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes 27ccaff030SJeremy L Thompson // 288355605fSKaren (Ren) Stengel //TESTARGS(name="solids-Linear-MMS") -ceed {ceed_resource} -test -degree 3 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3 298355605fSKaren (Ren) Stengel //TESTARGS(name="solids-NH1-1") -ceed {ceed_resource} -test -problem FSInitial-NH1 -E 2.8 -nu 0.4 -degree 2 -dm_plex_box_faces 2,2,2 -num_steps 1 -bc_clamp 6 -bc_traction 5 -bc_traction_5 0,0,-.5 -expect_final_strain_energy 2.124627916174e-01 308355605fSKaren (Ren) Stengel //TESTARGS(name="solids-MR1-1") -ceed {ceed_resource} -test -problem FSInitial-MR1 -mu_1 .5 -mu_2 .5 -nu 0.4 -degree 2 -dm_plex_box_faces 2,2,2 -num_steps 1 -bc_clamp 6 -bc_traction 5 -bc_traction_5 0,0,-.5 -expect_final_strain_energy 2.339138880207e-01 31ccaff030SJeremy L Thompson 32ccaff030SJeremy L Thompson /// @file 33ccaff030SJeremy L Thompson /// CEED elasticity example using PETSc with DMPlex 34ccaff030SJeremy L Thompson 35ccaff030SJeremy L Thompson const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n"; 36ccaff030SJeremy L Thompson 37ccaff030SJeremy L Thompson #include "elasticity.h" 38ccaff030SJeremy L Thompson 39ccaff030SJeremy L Thompson int main(int argc, char **argv) { 40ccaff030SJeremy L Thompson MPI_Comm comm; 41ccaff030SJeremy L Thompson // Context structs 42d1d35e2fSjeremylt AppCtx app_ctx; // Contains problem options 435754ecacSJeremy L Thompson ProblemFunctions problem_functions; // Setup functions for each problem 44ccaff030SJeremy L Thompson Units units; // Contains units scaling 45ccaff030SJeremy L Thompson // PETSc objects 462b730f8bSJeremy L Thompson PetscLogStage stage_dm_setup, stage_libceed_setup, stage_snes_setup, stage_snes_solve; 47d1d35e2fSjeremylt DM dm_orig; // Distributed DM to clone 48d1d35e2fSjeremylt DM dm_energy, dm_diagnostic; // DMs for postprocessing 49d1d35e2fSjeremylt DM *level_dms; 50d1d35e2fSjeremylt Vec U, *U_g, *U_loc; // U: solution, R: residual, F: forcing 51d1d35e2fSjeremylt Vec R, R_loc, F, F_loc; // g: global, loc: local 52d1d35e2fSjeremylt Vec neumann_bcs = NULL, bcs_loc = NULL; 5317f0843fSJeremy L Thompson SNES snes; 54d1d35e2fSjeremylt Mat *jacob_mat, jacob_mat_coarse, *prolong_restr_mat; 55ccaff030SJeremy L Thompson // PETSc data 56d1d35e2fSjeremylt UserMult res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx; 57d1d35e2fSjeremylt FormJacobCtx form_jacob_ctx; 58d1d35e2fSjeremylt UserMultProlongRestr *prolong_restr_ctx; 59d1d35e2fSjeremylt PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V; 60ccaff030SJeremy L Thompson // libCEED objects 61d99fa3c5SJeremy L Thompson Ceed ceed; 62d1d35e2fSjeremylt CeedData *ceed_data; 63d1d35e2fSjeremylt CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL; 64ccaff030SJeremy L Thompson // Parameters 65d1d35e2fSjeremylt PetscInt num_comp_u = 3; // 3 DoFs in 3D 66d1d35e2fSjeremylt PetscInt num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic 67d1d35e2fSjeremylt PetscInt num_levels = 1, fine_level = 0; 68bf0f51feSjeremylt PetscInt *U_g_size, *U_l_size, *U_loc_size; 69bf0f51feSjeremylt PetscInt snes_its = 0, ksp_its = 0; 70d1d35e2fSjeremylt double start_time, elapsed_time, min_time, max_time; 71ccaff030SJeremy L Thompson 722b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 73ccaff030SJeremy L Thompson 74ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 75ccaff030SJeremy L Thompson // Process command line options 76ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 77ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 78ccaff030SJeremy L Thompson 79ccaff030SJeremy L Thompson // -- Set mesh file, polynomial degree, problem type 802b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &app_ctx)); 812b730f8bSJeremy L Thompson PetscCall(ProcessCommandLineOptions(comm, app_ctx)); 822b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &problem_functions)); 832b730f8bSJeremy L Thompson PetscCall(RegisterProblems(problem_functions)); 84d1d35e2fSjeremylt num_levels = app_ctx->num_levels; 85d1d35e2fSjeremylt fine_level = num_levels - 1; 86ccaff030SJeremy L Thompson 87ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 88d1d35e2fSjeremylt // Initialize libCEED 8962e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 90d1d35e2fSjeremylt // Initialize backend 91d1d35e2fSjeremylt CeedInit(app_ctx->ceed_resource, &ceed); 9262e9c006SJeremy L Thompson 9362e9c006SJeremy L Thompson // Check preferred MemType 94d1d35e2fSjeremylt CeedMemType mem_type_backend; 95d1d35e2fSjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 965754ecacSJeremy L Thompson // Setup physics context and wrap in libCEED object 975754ecacSJeremy L Thompson { 985754ecacSJeremy L Thompson PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *); 992b730f8bSJeremy L Thompson PetscCall(PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name, &SetupPhysics)); 1002b730f8bSJeremy L Thompson if (!SetupPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found", app_ctx->name); 1012b730f8bSJeremy L Thompson PetscCall((*SetupPhysics)(comm, ceed, &units, &ctx_phys)); 1022b730f8bSJeremy L Thompson PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext, CeedQFunctionContext *); 1032b730f8bSJeremy L Thompson PetscCall(PetscFunctionListFind(problem_functions->setupSmootherPhysics, app_ctx->name, &SetupSmootherPhysics)); 1042b730f8bSJeremy L Thompson if (!SetupSmootherPhysics) SETERRQ(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found", app_ctx->name); 1052b730f8bSJeremy L Thompson PetscCall((*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother)); 106777ff853SJeremy L Thompson } 107777ff853SJeremy L Thompson 10862e9c006SJeremy L Thompson // --------------------------------------------------------------------------- 109ccaff030SJeremy L Thompson // Setup DM 110ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 111ccaff030SJeremy L Thompson // Performance logging 1122b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup)); 1132b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(stage_dm_setup)); 114ccaff030SJeremy L Thompson 115ccaff030SJeremy L Thompson // -- Create distributed DM from mesh file 1162b730f8bSJeremy L Thompson PetscCall(CreateDistributedDM(comm, app_ctx, &dm_orig)); 117b68a8d79SJed Brown VecType vectype; 118d1d35e2fSjeremylt switch (mem_type_backend) { 1192b730f8bSJeremy L Thompson case CEED_MEM_HOST: 1202b730f8bSJeremy L Thompson vectype = VECSTANDARD; 1212b730f8bSJeremy L Thompson break; 122b68a8d79SJed Brown case CEED_MEM_DEVICE: { 123b68a8d79SJed Brown const char *resolved; 124b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 125b68a8d79SJed Brown if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 126b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 127b68a8d79SJed Brown else vectype = VECSTANDARD; 128b68a8d79SJed Brown } 129b68a8d79SJed Brown } 1302b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_orig, vectype)); 1312b730f8bSJeremy L Thompson PetscCall(DMPlexDistributeSetDefault(dm_orig, PETSC_FALSE)); 1322b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(dm_orig)); 133ccaff030SJeremy L Thompson 134ccaff030SJeremy L Thompson // -- Setup DM by polynomial degree 1352b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &level_dms)); 136d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 1372b730f8bSJeremy L Thompson PetscCall(DMClone(dm_orig, &level_dms[level])); 1382b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm_orig, &vectype)); 1392b730f8bSJeremy L Thompson PetscCall(DMSetVecType(level_dms[level], vectype)); 1402b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level], PETSC_TRUE, num_comp_u)); 141a3c02c40SJeremy L Thompson // -- Label field components for viewing 142a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 143a3c02c40SJeremy L Thompson PetscSection section; 1442b730f8bSJeremy L Thompson PetscCall(DMGetLocalSection(level_dms[level], §ion)); 1452b730f8bSJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "Displacement")); 1462b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX")); 1472b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY")); 1482b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ")); 149a3c02c40SJeremy L Thompson } 150a3c02c40SJeremy L Thompson 151a3c02c40SJeremy L Thompson // -- Setup postprocessing DMs 1522b730f8bSJeremy L Thompson PetscCall(DMClone(dm_orig, &dm_energy)); 1532b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_e)); 1542b730f8bSJeremy L Thompson PetscCall(DMClone(dm_orig, &dm_diagnostic)); 1552b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm_diagnostic, app_ctx, app_ctx->level_degrees[fine_level], PETSC_FALSE, num_comp_u + num_comp_d)); 1562b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_energy, vectype)); 1572b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_diagnostic, vectype)); 158a3c02c40SJeremy L Thompson { 159a3c02c40SJeremy L Thompson // -- Label field components for viewing 160a3c02c40SJeremy L Thompson // Empty name for conserved field (because there is only one field) 161a3c02c40SJeremy L Thompson PetscSection section; 1622b730f8bSJeremy L Thompson PetscCall(DMGetLocalSection(dm_diagnostic, §ion)); 1632b730f8bSJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "Diagnostics")); 1642b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "DisplacementX")); 1652b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "DisplacementY")); 1662b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "DisplacementZ")); 1672b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "Pressure")); 1682b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain")); 1692b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 5, "TraceE2")); 1702b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 6, "detJ")); 1712b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity")); 172ccaff030SJeremy L Thompson } 173ccaff030SJeremy L Thompson 174ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 175ccaff030SJeremy L Thompson // Setup solution and work vectors 176ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 177ccaff030SJeremy L Thompson // Allocate arrays 1782b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &U_g)); 1792b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &U_loc)); 1802b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &U_g_size)); 1812b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &U_l_size)); 1822b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &U_loc_size)); 183ccaff030SJeremy L Thompson 184ccaff030SJeremy L Thompson // -- Setup solution vectors for each level 185d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 186ccaff030SJeremy L Thompson // -- Create global unknown vector U 1872b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(level_dms[level], &U_g[level])); 1882b730f8bSJeremy L Thompson PetscCall(VecGetSize(U_g[level], &U_g_size[level])); 189ccaff030SJeremy L Thompson // Note: Local size for matShell 1902b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(U_g[level], &U_l_size[level])); 191ccaff030SJeremy L Thompson 192d1d35e2fSjeremylt // -- Create local unknown vector U_loc 1932b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(level_dms[level], &U_loc[level])); 194ccaff030SJeremy L Thompson // Note: local size for libCEED 1952b730f8bSJeremy L Thompson PetscCall(VecGetSize(U_loc[level], &U_loc_size[level])); 196ccaff030SJeremy L Thompson } 197ccaff030SJeremy L Thompson 198ccaff030SJeremy L Thompson // -- Create residual and forcing vectors 1992b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_g[fine_level], &U)); 2002b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_g[fine_level], &R)); 2012b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_g[fine_level], &F)); 2022b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_loc[fine_level], &R_loc)); 2032b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_loc[fine_level], &F_loc)); 204ccaff030SJeremy L Thompson 205ccaff030SJeremy L Thompson // Performance logging 2062b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 207ccaff030SJeremy L Thompson 208ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 209ccaff030SJeremy L Thompson // Set up libCEED 210ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 211ccaff030SJeremy L Thompson // Performance logging 2122b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup)); 2132b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(stage_libceed_setup)); 214ccaff030SJeremy L Thompson 215ccaff030SJeremy L Thompson // -- Create libCEED local forcing vector 216d1d35e2fSjeremylt CeedVector force_ceed; 217ccaff030SJeremy L Thompson CeedScalar *f; 218d1d35e2fSjeremylt PetscMemType force_mem_type; 219d1d35e2fSjeremylt if (app_ctx->forcing_choice != FORCE_NONE) { 2202b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(F_loc, &f, &force_mem_type)); 221d1d35e2fSjeremylt CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed); 222d1d35e2fSjeremylt CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f); 223ccaff030SJeremy L Thompson } 224ccaff030SJeremy L Thompson 225fe394131Sjeremylt // -- Create libCEED local Neumann BCs vector 226d1d35e2fSjeremylt CeedVector neumann_ceed; 227fe394131Sjeremylt CeedScalar *n; 228d1d35e2fSjeremylt PetscMemType nummann_mem_type; 229d1d35e2fSjeremylt if (app_ctx->bc_traction_count > 0) { 2302b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U, &neumann_bcs)); 2312b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U_loc[fine_level], &bcs_loc)); 2322b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type)); 233d1d35e2fSjeremylt CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed); 2342b730f8bSJeremy L Thompson CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type), CEED_USE_POINTER, n); 235fe394131Sjeremylt } 236fe394131Sjeremylt 237ccaff030SJeremy L Thompson // -- Setup libCEED objects 2382b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &ceed_data)); 239d99fa3c5SJeremy L Thompson // ---- Setup residual, Jacobian evaluator and geometric information 2402b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &ceed_data[fine_level])); 241ccaff030SJeremy L Thompson { 2422b730f8bSJeremy L Thompson PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx, CeedQFunctionContext, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector, 2432b730f8bSJeremy L Thompson CeedVector, CeedData *); 2442b730f8bSJeremy L Thompson PetscCall(PetscFunctionListFind(problem_functions->setupLibceedFineLevel, app_ctx->name, &SetupLibceedFineLevel)); 2452b730f8bSJeremy L Thompson if (!SetupLibceedFineLevel) SETERRQ(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found", app_ctx->name); 2462b730f8bSJeremy L Thompson PetscCall((*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic, ceed, app_ctx, ctx_phys, fine_level, num_comp_u, 2472b730f8bSJeremy L Thompson U_g_size[fine_level], U_loc_size[fine_level], force_ceed, neumann_ceed, ceed_data)); 248ccaff030SJeremy L Thompson } 249d99fa3c5SJeremy L Thompson // ---- Setup coarse Jacobian evaluator and prolongation/restriction 250d1d35e2fSjeremylt for (PetscInt level = num_levels - 2; level >= 0; level--) { 2512b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &ceed_data[level])); 252d99fa3c5SJeremy L Thompson 253d99fa3c5SJeremy L Thompson // Get global communication restriction 2542b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(U_g[level + 1])); 2552b730f8bSJeremy L Thompson PetscCall(VecSet(U_loc[level + 1], 1.0)); 2562b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(level_dms[level + 1], U_loc[level + 1], ADD_VALUES, U_g[level + 1])); 2572b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(level_dms[level + 1], U_g[level + 1], INSERT_VALUES, U_loc[level + 1])); 258d99fa3c5SJeremy L Thompson 259d99fa3c5SJeremy L Thompson // Place in libCEED array 260d99fa3c5SJeremy L Thompson const PetscScalar *m; 261d1d35e2fSjeremylt PetscMemType m_mem_type; 2622b730f8bSJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(U_loc[level + 1], &m, &m_mem_type)); 2632b730f8bSJeremy L Thompson CeedVectorSetArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, (CeedScalar *)m); 264ccaff030SJeremy L Thompson 2651c8142b3Sjeremylt // Note: use high order ceed, if specified and degree > 4 2662b730f8bSJeremy L Thompson PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt, PetscInt, PetscInt, PetscInt, CeedVector, CeedData *); 2672b730f8bSJeremy L Thompson PetscCall(PetscFunctionListFind(problem_functions->setupLibceedLevel, app_ctx->name, &SetupLibceedLevel)); 2682b730f8bSJeremy L Thompson if (!SetupLibceedLevel) SETERRQ(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found", app_ctx->name); 2692b730f8bSJeremy L Thompson PetscCall((*SetupLibceedLevel)(level_dms[level], ceed, app_ctx, level, num_comp_u, U_g_size[level], U_loc_size[level], 2702b730f8bSJeremy L Thompson ceed_data[level + 1]->x_ceed, ceed_data)); 271d99fa3c5SJeremy L Thompson 272d99fa3c5SJeremy L Thompson // Restore PETSc vector 2732b730f8bSJeremy L Thompson CeedVectorTakeArray(ceed_data[level + 1]->x_ceed, MemTypeP2C(m_mem_type), (CeedScalar **)&m); 2742b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(U_loc[level + 1], &m)); 2752b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(U_g[level + 1])); 2762b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(U_loc[level + 1])); 277ccaff030SJeremy L Thompson } 278ccaff030SJeremy L Thompson 279ccaff030SJeremy L Thompson // Performance logging 2802b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 281ccaff030SJeremy L Thompson 282ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 283fe394131Sjeremylt // Setup global forcing and Neumann BC vectors 284ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 2852b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(F)); 286ccaff030SJeremy L Thompson 287d1d35e2fSjeremylt if (app_ctx->forcing_choice != FORCE_NONE) { 288d1d35e2fSjeremylt CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL); 2892b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(F_loc, &f)); 2902b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F)); 291d1d35e2fSjeremylt CeedVectorDestroy(&force_ceed); 292ccaff030SJeremy L Thompson } 293ccaff030SJeremy L Thompson 294d1d35e2fSjeremylt if (app_ctx->bc_traction_count > 0) { 2952b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(neumann_bcs)); 296d1d35e2fSjeremylt CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL); 2972b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(bcs_loc, &n)); 2982b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs)); 299d1d35e2fSjeremylt CeedVectorDestroy(&neumann_ceed); 300fe394131Sjeremylt } 301fe394131Sjeremylt 302ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 303ccaff030SJeremy L Thompson // Print problem summary 304ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 305d1d35e2fSjeremylt if (!app_ctx->test_mode) { 306ccaff030SJeremy L Thompson const char *usedresource; 307ccaff030SJeremy L Thompson CeedGetResource(ceed, &usedresource); 30877841947SLeila Ghaffari char hostname[PETSC_MAX_PATH_LEN]; 3092b730f8bSJeremy L Thompson PetscCall(PetscGetHostName(hostname, sizeof hostname)); 31077841947SLeila Ghaffari PetscInt comm_size; 3112b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(comm, &comm_size)); 312ccaff030SJeremy L Thompson 3132b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 31417fd649aSJeremy L Thompson "\n-- Elasticity Example - libCEED + PETSc --\n" 31577841947SLeila Ghaffari " MPI:\n" 31677841947SLeila Ghaffari " Hostname : %s\n" 31777841947SLeila Ghaffari " Total ranks : %d\n" 318ccaff030SJeremy L Thompson " libCEED:\n" 31962e9c006SJeremy L Thompson " libCEED Backend : %s\n" 320b68a8d79SJed Brown " libCEED Backend MemType : %s\n", 3212b730f8bSJeremy L Thompson hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend])); 322ccaff030SJeremy L Thompson 32362e9c006SJeremy L Thompson VecType vecType; 3242b730f8bSJeremy L Thompson PetscCall(VecGetType(U, &vecType)); 3252b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 32662e9c006SJeremy L Thompson " PETSc:\n" 32762e9c006SJeremy L Thompson " PETSc Vec Type : %s\n", 3282b730f8bSJeremy L Thompson vecType)); 32962e9c006SJeremy L Thompson 3302b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 331ccaff030SJeremy L Thompson " Problem:\n" 332ccaff030SJeremy L Thompson " Problem Name : %s\n" 333ccaff030SJeremy L Thompson " Forcing Function : %s\n" 334ccaff030SJeremy L Thompson " Mesh:\n" 335ccaff030SJeremy L Thompson " File : %s\n" 336990fdeb6SJeremy L Thompson " Number of 1D Basis Nodes (p) : %" CeedInt_FMT "\n" 337990fdeb6SJeremy L Thompson " Number of 1D Quadrature Points (q) : %" CeedInt_FMT "\n" 33808140895SJed Brown " Global nodes : %" PetscInt_FMT "\n" 33908140895SJed Brown " Owned nodes : %" PetscInt_FMT "\n" 34008140895SJed Brown " DoF per node : %" PetscInt_FMT "\n" 341ccaff030SJeremy L Thompson " Multigrid:\n" 342ccaff030SJeremy L Thompson " Type : %s\n" 343990fdeb6SJeremy L Thompson " Number of Levels : %" CeedInt_FMT "\n", 3442b730f8bSJeremy L Thompson app_ctx->name_for_disp, forcing_types_for_disp[app_ctx->forcing_choice], 3452b730f8bSJeremy L Thompson app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh", app_ctx->degree + 1, app_ctx->degree + 1, 3462b730f8bSJeremy L Thompson U_g_size[fine_level] / num_comp_u, U_l_size[fine_level] / num_comp_u, num_comp_u, 3472b730f8bSJeremy L Thompson (app_ctx->degree == 1 && app_ctx->multigrid_choice != MULTIGRID_NONE) ? "Algebraic multigrid" 3482b730f8bSJeremy L Thompson : multigrid_types_for_disp[app_ctx->multigrid_choice], 3492b730f8bSJeremy L Thompson (app_ctx->degree == 1 || app_ctx->multigrid_choice == MULTIGRID_NONE) ? 0 : num_levels)); 350ccaff030SJeremy L Thompson 351d1d35e2fSjeremylt if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 352e3e3df41Sjeremylt for (PetscInt i = 0; i < 2; i++) { 353d1d35e2fSjeremylt CeedInt level = i ? fine_level : 0; 3542b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 35508140895SJed Brown " Level %" PetscInt_FMT " (%s):\n" 356990fdeb6SJeremy L Thompson " Number of 1D Basis Nodes (p) : %" CeedInt_FMT "\n" 35708140895SJed Brown " Global Nodes : %" PetscInt_FMT "\n" 35808140895SJed Brown " Owned Nodes : %" PetscInt_FMT "\n", 3592b730f8bSJeremy L Thompson level, i ? "fine" : "coarse", app_ctx->level_degrees[level] + 1, U_g_size[level] / num_comp_u, 3602b730f8bSJeremy L Thompson U_l_size[level] / num_comp_u)); 361ccaff030SJeremy L Thompson } 362ccaff030SJeremy L Thompson } 363ccaff030SJeremy L Thompson } 364ccaff030SJeremy L Thompson 365ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 366ccaff030SJeremy L Thompson // Setup SNES 367ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 368ccaff030SJeremy L Thompson // Performance logging 3692b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup)); 3702b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(stage_snes_setup)); 371ccaff030SJeremy L Thompson 372ccaff030SJeremy L Thompson // Create SNES 3732b730f8bSJeremy L Thompson PetscCall(SNESCreate(comm, &snes)); 3742b730f8bSJeremy L Thompson PetscCall(SNESSetDM(snes, level_dms[fine_level])); 375ccaff030SJeremy L Thompson 376ccaff030SJeremy L Thompson // -- Jacobian evaluators 3772b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &jacob_ctx)); 3782b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &jacob_mat)); 379d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 380ccaff030SJeremy L Thompson // -- Jacobian context for level 3812b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &jacob_ctx[level])); 3822b730f8bSJeremy L Thompson PetscCall(SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level], U_loc[level], ceed_data[level], ceed, ctx_phys, ctx_phys_smoother, 3832b730f8bSJeremy L Thompson jacob_ctx[level])); 384ccaff030SJeremy L Thompson 385ccaff030SJeremy L Thompson // -- Form Action of Jacobian on delta_u 3862b730f8bSJeremy L Thompson PetscCall(MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level], U_g_size[level], jacob_ctx[level], &jacob_mat[level])); 3872b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_MULT, (void (*)(void))ApplyJacobian_Ceed)); 3882b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL, (void (*)(void))GetDiag_Ceed)); 3892b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(jacob_mat[level], vectype)); 390ccaff030SJeremy L Thompson } 391*ea61e9acSJeremy L Thompson // Note: FormJacobian updates Jacobian matrices on each level and assembles the Jpre matrix, if needed 3922b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &form_jacob_ctx)); 393d1d35e2fSjeremylt form_jacob_ctx->jacob_ctx = jacob_ctx; 394d1d35e2fSjeremylt form_jacob_ctx->num_levels = num_levels; 395d1d35e2fSjeremylt form_jacob_ctx->jacob_mat = jacob_mat; 396ccaff030SJeremy L Thompson 397ccaff030SJeremy L Thompson // -- Residual evaluation function 3982b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &res_ctx)); 3992b730f8bSJeremy L Thompson PetscCall(PetscMemcpy(res_ctx, jacob_ctx[fine_level], sizeof(*jacob_ctx[fine_level]))); 4005754ecacSJeremy L Thompson res_ctx->op = ceed_data[fine_level]->op_residual; 4015754ecacSJeremy L Thompson res_ctx->qf = ceed_data[fine_level]->qf_residual; 4022b730f8bSJeremy L Thompson if (app_ctx->bc_traction_count > 0) res_ctx->neumann_bcs = neumann_bcs; 4032b730f8bSJeremy L Thompson else res_ctx->neumann_bcs = NULL; 4042b730f8bSJeremy L Thompson PetscCall(SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx)); 405ccaff030SJeremy L Thompson 406ccaff030SJeremy L Thompson // -- Prolongation/Restriction evaluation 4072b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &prolong_restr_ctx)); 4082b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &prolong_restr_mat)); 409d1d35e2fSjeremylt for (PetscInt level = 1; level < num_levels; level++) { 410ccaff030SJeremy L Thompson // ---- Prolongation/restriction context for level 4112b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &prolong_restr_ctx[level])); 4122b730f8bSJeremy L Thompson PetscCall(SetupProlongRestrictCtx(comm, app_ctx, level_dms[level - 1], level_dms[level], U_g[level], U_loc[level - 1], U_loc[level], 4132b730f8bSJeremy L Thompson ceed_data[level - 1], ceed_data[level], ceed, prolong_restr_ctx[level])); 414ccaff030SJeremy L Thompson 415ccaff030SJeremy L Thompson // ---- Form Action of Jacobian on delta_u 4162b730f8bSJeremy L Thompson PetscCall(MatCreateShell(comm, U_l_size[level], U_l_size[level - 1], U_g_size[level], U_g_size[level - 1], prolong_restr_ctx[level], 4172b730f8bSJeremy L Thompson &prolong_restr_mat[level])); 418ccaff030SJeremy L Thompson // Note: In PCMG, restriction is the transpose of prolongation 4192b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT, (void (*)(void))Prolong_Ceed)); 4202b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE, (void (*)(void))Restrict_Ceed)); 4212b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(prolong_restr_mat[level], vectype)); 422ccaff030SJeremy L Thompson } 423ccaff030SJeremy L Thompson 424ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 42517f0843fSJeremy L Thompson // Setup for AMG coarse solve 426ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 427d1d35e2fSjeremylt if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 428e3e3df41Sjeremylt // -- Jacobian Matrix 4292b730f8bSJeremy L Thompson PetscCall(DMCreateMatrix(level_dms[0], &jacob_mat_coarse)); 430e3e3df41Sjeremylt 431d1d35e2fSjeremylt if (app_ctx->degree > 1) { 43217f0843fSJeremy L Thompson // -- Assemble sparsity pattern 4333047f789SJeremy L Thompson PetscCount num_entries; 4343047f789SJeremy L Thompson CeedInt *rows, *cols; 43517f0843fSJeremy L Thompson CeedVector coo_values; 4362b730f8bSJeremy L Thompson CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries, &rows, &cols); 43717f0843fSJeremy L Thompson ISLocalToGlobalMapping ltog_row, ltog_col; 4382b730f8bSJeremy L Thompson PetscCall(MatGetLocalToGlobalMapping(jacob_mat_coarse, <og_row, <og_col)); 4392b730f8bSJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows)); 4402b730f8bSJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols)); 4412b730f8bSJeremy L Thompson PetscCall(MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows, cols)); 44217f0843fSJeremy L Thompson free(rows); 44317f0843fSJeremy L Thompson free(cols); 44417f0843fSJeremy L Thompson CeedVectorCreate(ceed, num_entries, &coo_values); 445ccaff030SJeremy L Thompson 446d1d35e2fSjeremylt // -- Update form_jacob_ctx 44717f0843fSJeremy L Thompson form_jacob_ctx->coo_values = coo_values; 44817f0843fSJeremy L Thompson form_jacob_ctx->op_coarse = ceed_data[0]->op_jacobian; 449d1d35e2fSjeremylt form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse; 450e3e3df41Sjeremylt } 451e3e3df41Sjeremylt } 452e3e3df41Sjeremylt 453e3e3df41Sjeremylt // Set Jacobian function 454d1d35e2fSjeremylt if (app_ctx->degree > 1) { 4552b730f8bSJeremy L Thompson PetscCall(SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level], FormJacobian, form_jacob_ctx)); 456e3e3df41Sjeremylt } else { 4572b730f8bSJeremy L Thompson PetscCall(SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse, SNESComputeJacobianDefaultColor, NULL)); 458e3e3df41Sjeremylt } 459ccaff030SJeremy L Thompson 460ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 461ccaff030SJeremy L Thompson // Setup KSP 462ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 463ccaff030SJeremy L Thompson { 464ccaff030SJeremy L Thompson PC pc; 465ccaff030SJeremy L Thompson KSP ksp; 466ccaff030SJeremy L Thompson 467ccaff030SJeremy L Thompson // -- KSP 4682b730f8bSJeremy L Thompson PetscCall(SNESGetKSP(snes, &ksp)); 4692b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG)); 4702b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 4712b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 4722b730f8bSJeremy L Thompson PetscCall(KSPSetOptionsPrefix(ksp, "outer_")); 473ccaff030SJeremy L Thompson 474ccaff030SJeremy L Thompson // -- Preconditioning 4752b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 4762b730f8bSJeremy L Thompson PetscCall(PCSetDM(pc, level_dms[fine_level])); 4772b730f8bSJeremy L Thompson PetscCall(PCSetOptionsPrefix(pc, "outer_")); 478ccaff030SJeremy L Thompson 479d1d35e2fSjeremylt if (app_ctx->multigrid_choice == MULTIGRID_NONE) { 480ccaff030SJeremy L Thompson // ---- No Multigrid 4812b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCJACOBI)); 4822b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 483d1d35e2fSjeremylt } else if (app_ctx->degree == 1) { 484f892e6d1Sjeremylt // ---- AMG for degree 1 4852b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCGAMG)); 486ccaff030SJeremy L Thompson } else { 487ccaff030SJeremy L Thompson // ---- PCMG 4882b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCMG)); 489ccaff030SJeremy L Thompson 490ccaff030SJeremy L Thompson // ------ PCMG levels 4912b730f8bSJeremy L Thompson PetscCall(PCMGSetLevels(pc, num_levels, NULL)); 492d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 493ccaff030SJeremy L Thompson // -------- Smoother 494d1d35e2fSjeremylt KSP ksp_smoother, ksp_est; 495d1d35e2fSjeremylt PC pc_smoother; 496ccaff030SJeremy L Thompson 497ccaff030SJeremy L Thompson // ---------- Smoother KSP 4982b730f8bSJeremy L Thompson PetscCall(PCMGGetSmoother(pc, level, &ksp_smoother)); 4992b730f8bSJeremy L Thompson PetscCall(KSPSetDM(ksp_smoother, level_dms[level])); 5002b730f8bSJeremy L Thompson PetscCall(KSPSetDMActive(ksp_smoother, PETSC_FALSE)); 501ccaff030SJeremy L Thompson 502ccaff030SJeremy L Thompson // ---------- Chebyshev options 5032b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp_smoother, KSPCHEBYSHEV)); 5042b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1)); 5052b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est)); 5062b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp_est, KSPCG)); 5072b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE)); 5082b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level])); 509ccaff030SJeremy L Thompson 510ccaff030SJeremy L Thompson // ---------- Smoother preconditioner 5112b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp_smoother, &pc_smoother)); 5122b730f8bSJeremy L Thompson PetscCall(PCSetType(pc_smoother, PCJACOBI)); 5132b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL)); 514ccaff030SJeremy L Thompson 515ccaff030SJeremy L Thompson // -------- Work vector 516d1d35e2fSjeremylt if (level != fine_level) { 5172b730f8bSJeremy L Thompson PetscCall(PCMGSetX(pc, level, U_g[level])); 518ccaff030SJeremy L Thompson } 519ccaff030SJeremy L Thompson 520ccaff030SJeremy L Thompson // -------- Level prolongation/restriction operator 521ccaff030SJeremy L Thompson if (level > 0) { 5222b730f8bSJeremy L Thompson PetscCall(PCMGSetInterpolation(pc, level, prolong_restr_mat[level])); 5232b730f8bSJeremy L Thompson PetscCall(PCMGSetRestriction(pc, level, prolong_restr_mat[level])); 524ccaff030SJeremy L Thompson } 525ccaff030SJeremy L Thompson } 526ccaff030SJeremy L Thompson 527ccaff030SJeremy L Thompson // ------ PCMG coarse solve 528d1d35e2fSjeremylt KSP ksp_coarse; 529d1d35e2fSjeremylt PC pc_coarse; 530ccaff030SJeremy L Thompson 531ccaff030SJeremy L Thompson // -------- Coarse KSP 5322b730f8bSJeremy L Thompson PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse)); 5332b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp_coarse, KSPPREONLY)); 5342b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse)); 5352b730f8bSJeremy L Thompson PetscCall(KSPSetOptionsPrefix(ksp_coarse, "coarse_")); 536ccaff030SJeremy L Thompson 537ccaff030SJeremy L Thompson // -------- Coarse preconditioner 5382b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp_coarse, &pc_coarse)); 5392b730f8bSJeremy L Thompson PetscCall(PCSetType(pc_coarse, PCGAMG)); 5402b730f8bSJeremy L Thompson PetscCall(PCSetOptionsPrefix(pc_coarse, "coarse_")); 541ccaff030SJeremy L Thompson 5422b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp_coarse)); 5432b730f8bSJeremy L Thompson PetscCall(PCSetFromOptions(pc_coarse)); 544ccaff030SJeremy L Thompson 545ccaff030SJeremy L Thompson // ------ PCMG options 5462b730f8bSJeremy L Thompson PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE)); 5472b730f8bSJeremy L Thompson PetscCall(PCMGSetNumberSmooth(pc, 3)); 5482b730f8bSJeremy L Thompson PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type)); 549ccaff030SJeremy L Thompson } 5502b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp)); 5512b730f8bSJeremy L Thompson PetscCall(PCSetFromOptions(pc)); 552ccaff030SJeremy L Thompson } 553038d0cd7Sjeremylt { 554038d0cd7Sjeremylt // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 555d1d35e2fSjeremylt SNESLineSearch line_search; 556481a4840SJed Brown 5572b730f8bSJeremy L Thompson PetscCall(SNESGetLineSearch(snes, &line_search)); 5582b730f8bSJeremy L Thompson PetscCall(SNESLineSearchSetType(line_search, SNESLINESEARCHCP)); 559481a4840SJed Brown } 560481a4840SJed Brown 5612b730f8bSJeremy L Thompson PetscCall(SNESSetFromOptions(snes)); 562ccaff030SJeremy L Thompson 563ccaff030SJeremy L Thompson // Performance logging 5642b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 565ccaff030SJeremy L Thompson 566ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 567ccaff030SJeremy L Thompson // Set initial guess 568ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 5692b730f8bSJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)U, "")); 5702b730f8bSJeremy L Thompson PetscCall(VecSet(U, 0.0)); 571ccaff030SJeremy L Thompson 572ccaff030SJeremy L Thompson // View solution 573d1d35e2fSjeremylt if (app_ctx->view_soln) { 5742b730f8bSJeremy L Thompson PetscCall(ViewSolution(comm, app_ctx, U, 0, 0.0)); 575ccaff030SJeremy L Thompson } 576ccaff030SJeremy L Thompson 577ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 578ccaff030SJeremy L Thompson // Solve SNES 579ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 580d1d35e2fSjeremylt PetscBool snes_monitor = PETSC_FALSE; 5812b730f8bSJeremy L Thompson PetscCall(PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor)); 5825f24f028Sjeremylt 583ccaff030SJeremy L Thompson // Performance logging 5842b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve)); 5852b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(stage_snes_solve)); 586ccaff030SJeremy L Thompson 587ccaff030SJeremy L Thompson // Timing 5882b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)snes)); 589d1d35e2fSjeremylt start_time = MPI_Wtime(); 590ccaff030SJeremy L Thompson 591ccaff030SJeremy L Thompson // Solve for each load increment 5925f24f028Sjeremylt PetscInt increment; 593d1d35e2fSjeremylt for (increment = 1; increment <= app_ctx->num_increments; increment++) { 5945f24f028Sjeremylt // -- Log increment count 595d1d35e2fSjeremylt if (snes_monitor) { 5962b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, "%" PetscInt_FMT " Load Increment\n", increment - 1)); 5975f24f028Sjeremylt } 5985f24f028Sjeremylt 599ccaff030SJeremy L Thompson // -- Scale the problem 600d1d35e2fSjeremylt PetscScalar load_increment = 1.0 * increment / app_ctx->num_increments, 6012b730f8bSJeremy L Thompson scalingFactor = load_increment / (increment == 1 ? 1 : res_ctx->load_increment); 602d1d35e2fSjeremylt res_ctx->load_increment = load_increment; 603d1d35e2fSjeremylt if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) { 6042b730f8bSJeremy L Thompson PetscCall(VecScale(F, scalingFactor)); 605ccaff030SJeremy L Thompson } 606ccaff030SJeremy L Thompson 607ccaff030SJeremy L Thompson // -- Solve 6082b730f8bSJeremy L Thompson PetscCall(SNESSolve(snes, F, U)); 609ccaff030SJeremy L Thompson 610ccaff030SJeremy L Thompson // -- View solution 611d1d35e2fSjeremylt if (app_ctx->view_soln) { 6122b730f8bSJeremy L Thompson PetscCall(ViewSolution(comm, app_ctx, U, increment, load_increment)); 613ccaff030SJeremy L Thompson } 614ccaff030SJeremy L Thompson 615ccaff030SJeremy L Thompson // -- Update SNES iteration count 616ccaff030SJeremy L Thompson PetscInt its; 6172b730f8bSJeremy L Thompson PetscCall(SNESGetIterationNumber(snes, &its)); 618bf0f51feSjeremylt snes_its += its; 6192b730f8bSJeremy L Thompson PetscCall(SNESGetLinearSolveIterations(snes, &its)); 620bf0f51feSjeremylt ksp_its += its; 6213951731eSjeremylt 6223951731eSjeremylt // -- Check for divergence 6233951731eSjeremylt SNESConvergedReason reason; 6242b730f8bSJeremy L Thompson PetscCall(SNESGetConvergedReason(snes, &reason)); 6252b730f8bSJeremy L Thompson if (reason < 0) break; 6268355605fSKaren (Ren) Stengel if (app_ctx->energy_viewer) { 6278355605fSKaren (Ren) Stengel // -- Log strain energy for current load increment 6288355605fSKaren (Ren) Stengel CeedScalar energy; 6292b730f8bSJeremy L Thompson PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy)); 6308355605fSKaren (Ren) Stengel 6318355605fSKaren (Ren) Stengel if (!app_ctx->test_mode) { 6328355605fSKaren (Ren) Stengel // -- Output 6332b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Strain Energy : %.12e\n", energy)); 6348355605fSKaren (Ren) Stengel } 6352b730f8bSJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment, energy)); 6368355605fSKaren (Ren) Stengel } 637ccaff030SJeremy L Thompson } 638ccaff030SJeremy L Thompson 639ccaff030SJeremy L Thompson // Timing 640d1d35e2fSjeremylt elapsed_time = MPI_Wtime() - start_time; 641ccaff030SJeremy L Thompson 642ccaff030SJeremy L Thompson // Performance logging 6432b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 644ccaff030SJeremy L Thompson 645ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 646ccaff030SJeremy L Thompson // Output summary 647ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 648d1d35e2fSjeremylt if (!app_ctx->test_mode) { 649ccaff030SJeremy L Thompson // -- SNES 650d1d35e2fSjeremylt SNESType snes_type; 651ccaff030SJeremy L Thompson SNESConvergedReason reason; 652ccaff030SJeremy L Thompson PetscReal rnorm; 6532b730f8bSJeremy L Thompson PetscCall(SNESGetType(snes, &snes_type)); 6542b730f8bSJeremy L Thompson PetscCall(SNESGetConvergedReason(snes, &reason)); 6552b730f8bSJeremy L Thompson PetscCall(SNESGetFunctionNorm(snes, &rnorm)); 6562b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 657ccaff030SJeremy L Thompson " SNES:\n" 658ccaff030SJeremy L Thompson " SNES Type : %s\n" 659ccaff030SJeremy L Thompson " SNES Convergence : %s\n" 660990fdeb6SJeremy L Thompson " Number of Load Increments : %" PetscInt_FMT "\n" 661990fdeb6SJeremy L Thompson " Completed Load Increments : %" PetscInt_FMT "\n" 66208140895SJed Brown " Total SNES Iterations : %" PetscInt_FMT "\n" 663ccaff030SJeremy L Thompson " Final rnorm : %e\n", 6642b730f8bSJeremy L Thompson snes_type, SNESConvergedReasons[reason], app_ctx->num_increments, increment - 1, snes_its, (double)rnorm)); 665ccaff030SJeremy L Thompson 666ccaff030SJeremy L Thompson // -- KSP 667ccaff030SJeremy L Thompson KSP ksp; 668d1d35e2fSjeremylt KSPType ksp_type; 6692b730f8bSJeremy L Thompson PetscCall(SNESGetKSP(snes, &ksp)); 6702b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type)); 6712b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 672ccaff030SJeremy L Thompson " Linear Solver:\n" 6737418ca88SJeremy L Thompson " KSP Type : %s\n" 67408140895SJed Brown " Total KSP Iterations : %" PetscInt_FMT "\n", 6752b730f8bSJeremy L Thompson ksp_type, ksp_its)); 676ccaff030SJeremy L Thompson 677ccaff030SJeremy L Thompson // -- PC 678ccaff030SJeremy L Thompson PC pc; 679d1d35e2fSjeremylt PCType pc_type; 6802b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 6812b730f8bSJeremy L Thompson PetscCall(PCGetType(pc, &pc_type)); 6822b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " PC Type : %s\n", pc_type)); 683e3e3df41Sjeremylt 684d1d35e2fSjeremylt if (!strcmp(pc_type, PCMG)) { 685d1d35e2fSjeremylt PCMGType pcmg_type; 6862b730f8bSJeremy L Thompson PetscCall(PCMGGetType(pc, &pcmg_type)); 6872b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 688ccaff030SJeremy L Thompson " P-Multigrid:\n" 689ccaff030SJeremy L Thompson " PCMG Type : %s\n" 690ccaff030SJeremy L Thompson " PCMG Cycle Type : %s\n", 6912b730f8bSJeremy L Thompson PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type])); 692ccaff030SJeremy L Thompson 693ccaff030SJeremy L Thompson // -- Coarse Solve 694d1d35e2fSjeremylt KSP ksp_coarse; 695d1d35e2fSjeremylt PC pc_coarse; 696d1d35e2fSjeremylt PCType pc_type; 697ccaff030SJeremy L Thompson 6982b730f8bSJeremy L Thompson PetscCall(PCMGGetCoarseSolve(pc, &ksp_coarse)); 6992b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp_coarse, &ksp_type)); 7002b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp_coarse, &pc_coarse)); 7012b730f8bSJeremy L Thompson PetscCall(PCGetType(pc_coarse, &pc_type)); 7022b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 703ccaff030SJeremy L Thompson " Coarse Solve:\n" 704ccaff030SJeremy L Thompson " KSP Type : %s\n" 705ccaff030SJeremy L Thompson " PC Type : %s\n", 7062b730f8bSJeremy L Thompson ksp_type, pc_type)); 707ccaff030SJeremy L Thompson } 708ccaff030SJeremy L Thompson } 709ccaff030SJeremy L Thompson 710ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 711ccaff030SJeremy L Thompson // Compute solve time 712ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 713d1d35e2fSjeremylt if (!app_ctx->test_mode) { 7142b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm)); 7152b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm)); 7162b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 717ccaff030SJeremy L Thompson " Performance:\n" 7187418ca88SJeremy L Thompson " SNES Solve Time : %g (%g) sec\n" 7197418ca88SJeremy L Thompson " DoFs/Sec in SNES : %g (%g) million\n", 7202b730f8bSJeremy L Thompson max_time, min_time, 1e-6 * U_g_size[fine_level] * ksp_its / max_time, 1e-6 * U_g_size[fine_level] * ksp_its / min_time)); 721ccaff030SJeremy L Thompson } 722ccaff030SJeremy L Thompson 723ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 724ccaff030SJeremy L Thompson // Compute error 725ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 726d1d35e2fSjeremylt if (app_ctx->forcing_choice == FORCE_MMS) { 727d1d35e2fSjeremylt CeedScalar l2_error = 1., l2_U_norm = 1.; 728d1d35e2fSjeremylt const CeedScalar *true_array; 729d1d35e2fSjeremylt Vec error_vec, true_vec; 730ccaff030SJeremy L Thompson 731ccaff030SJeremy L Thompson // -- Work vectors 7322b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U, &error_vec)); 7332b730f8bSJeremy L Thompson PetscCall(VecSet(error_vec, 0.0)); 7342b730f8bSJeremy L Thompson PetscCall(VecDuplicate(U, &true_vec)); 7352b730f8bSJeremy L Thompson PetscCall(VecSet(true_vec, 0.0)); 736ccaff030SJeremy L Thompson 737ccaff030SJeremy L Thompson // -- Assemble global true solution vector 7382b730f8bSJeremy L Thompson CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln, CEED_MEM_HOST, &true_array); 7392b730f8bSJeremy L Thompson PetscCall(VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array)); 7402b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec)); 7412b730f8bSJeremy L Thompson PetscCall(VecResetArray(res_ctx->Y_loc)); 742d1d35e2fSjeremylt CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array); 743ccaff030SJeremy L Thompson 744ccaff030SJeremy L Thompson // -- Compute L2 error 7452b730f8bSJeremy L Thompson PetscCall(VecWAXPY(error_vec, -1.0, U, true_vec)); 7462b730f8bSJeremy L Thompson PetscCall(VecNorm(error_vec, NORM_2, &l2_error)); 7472b730f8bSJeremy L Thompson PetscCall(VecNorm(U, NORM_2, &l2_U_norm)); 748d1d35e2fSjeremylt l2_error /= l2_U_norm; 749ccaff030SJeremy L Thompson 750ccaff030SJeremy L Thompson // -- Output 751d1d35e2fSjeremylt if (!app_ctx->test_mode || l2_error > 0.05) { 7522b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " L2 Error : %e\n", l2_error)); 753ccaff030SJeremy L Thompson } 754ccaff030SJeremy L Thompson 755ccaff030SJeremy L Thompson // -- Cleanup 7562b730f8bSJeremy L Thompson PetscCall(VecDestroy(&error_vec)); 7572b730f8bSJeremy L Thompson PetscCall(VecDestroy(&true_vec)); 7582d93065eSjeremylt } 7592d93065eSjeremylt 7602d93065eSjeremylt // --------------------------------------------------------------------------- 7612d93065eSjeremylt // Compute energy 7622d93065eSjeremylt // --------------------------------------------------------------------------- 7638355605fSKaren (Ren) Stengel PetscReal energy; 7642b730f8bSJeremy L Thompson PetscCall(ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, U, &energy)); 7658355605fSKaren (Ren) Stengel if (!app_ctx->test_mode) { 7662d93065eSjeremylt // -- Output 7672b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Strain Energy : %.12e\n", energy)); 768ccaff030SJeremy L Thompson } 7692b730f8bSJeremy L Thompson PetscCall(RegressionTests_solids(app_ctx, energy)); 770ccaff030SJeremy L Thompson 771ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 7725c25879aSJeremy L Thompson // Output diagnostic quantities 7735c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 774d1d35e2fSjeremylt if (app_ctx->view_soln || app_ctx->view_final_soln) { 7755c25879aSJeremy L Thompson // -- Setup context 776d1d35e2fSjeremylt UserMult diagnostic_ctx; 7772b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &diagnostic_ctx)); 7782b730f8bSJeremy L Thompson PetscCall(PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx))); 779d1d35e2fSjeremylt diagnostic_ctx->dm = dm_diagnostic; 780d1d35e2fSjeremylt diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic; 7815c25879aSJeremy L Thompson 7825c25879aSJeremy L Thompson // -- Compute and output 7832b730f8bSJeremy L Thompson PetscCall(ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx, app_ctx, U, ceed_data[fine_level]->elem_restr_diagnostic)); 7845c25879aSJeremy L Thompson 7855c25879aSJeremy L Thompson // -- Cleanup 7862b730f8bSJeremy L Thompson PetscCall(PetscFree(diagnostic_ctx)); 7875c25879aSJeremy L Thompson } 7885c25879aSJeremy L Thompson 7895c25879aSJeremy L Thompson // --------------------------------------------------------------------------- 790ccaff030SJeremy L Thompson // Free objects 791ccaff030SJeremy L Thompson // --------------------------------------------------------------------------- 792ccaff030SJeremy L Thompson // Data in arrays per level 793d1d35e2fSjeremylt for (PetscInt level = 0; level < num_levels; level++) { 794ccaff030SJeremy L Thompson // Vectors 7952b730f8bSJeremy L Thompson PetscCall(VecDestroy(&U_g[level])); 7962b730f8bSJeremy L Thompson PetscCall(VecDestroy(&U_loc[level])); 797ccaff030SJeremy L Thompson 798ccaff030SJeremy L Thompson // Jacobian matrix and data 7992b730f8bSJeremy L Thompson PetscCall(VecDestroy(&jacob_ctx[level]->Y_loc)); 8002b730f8bSJeremy L Thompson PetscCall(MatDestroy(&jacob_mat[level])); 8012b730f8bSJeremy L Thompson PetscCall(PetscFree(jacob_ctx[level])); 802ccaff030SJeremy L Thompson 803ccaff030SJeremy L Thompson // Prolongation/Restriction matrix and data 804ccaff030SJeremy L Thompson if (level > 0) { 8052b730f8bSJeremy L Thompson PetscCall(PetscFree(prolong_restr_ctx[level])); 8062b730f8bSJeremy L Thompson PetscCall(MatDestroy(&prolong_restr_mat[level])); 807ccaff030SJeremy L Thompson } 808ccaff030SJeremy L Thompson 809ccaff030SJeremy L Thompson // DM 8102b730f8bSJeremy L Thompson PetscCall(DMDestroy(&level_dms[level])); 811ccaff030SJeremy L Thompson 812ccaff030SJeremy L Thompson // libCEED objects 8132b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(level, ceed_data[level])); 814ccaff030SJeremy L Thompson } 815ccaff030SJeremy L Thompson 8162b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&app_ctx->energy_viewer)); 8178355605fSKaren (Ren) Stengel 818ccaff030SJeremy L Thompson // Arrays 8192b730f8bSJeremy L Thompson PetscCall(PetscFree(U_g)); 8202b730f8bSJeremy L Thompson PetscCall(PetscFree(U_loc)); 8212b730f8bSJeremy L Thompson PetscCall(PetscFree(U_g_size)); 8222b730f8bSJeremy L Thompson PetscCall(PetscFree(U_l_size)); 8232b730f8bSJeremy L Thompson PetscCall(PetscFree(U_loc_size)); 8242b730f8bSJeremy L Thompson PetscCall(PetscFree(jacob_ctx)); 8252b730f8bSJeremy L Thompson PetscCall(PetscFree(jacob_mat)); 8262b730f8bSJeremy L Thompson PetscCall(PetscFree(prolong_restr_ctx)); 8272b730f8bSJeremy L Thompson PetscCall(PetscFree(prolong_restr_mat)); 8282b730f8bSJeremy L Thompson PetscCall(PetscFree(app_ctx->level_degrees)); 8292b730f8bSJeremy L Thompson PetscCall(PetscFree(ceed_data)); 830ccaff030SJeremy L Thompson 831ccaff030SJeremy L Thompson // libCEED objects 83217f0843fSJeremy L Thompson CeedVectorDestroy(&form_jacob_ctx->coo_values); 833d1d35e2fSjeremylt CeedQFunctionContextDestroy(&ctx_phys); 834d1d35e2fSjeremylt CeedQFunctionContextDestroy(&ctx_phys_smoother); 835ccaff030SJeremy L Thompson CeedDestroy(&ceed); 836ccaff030SJeremy L Thompson 837ccaff030SJeremy L Thompson // PETSc objects 8382b730f8bSJeremy L Thompson PetscCall(VecDestroy(&U)); 8392b730f8bSJeremy L Thompson PetscCall(VecDestroy(&R)); 8402b730f8bSJeremy L Thompson PetscCall(VecDestroy(&R_loc)); 8412b730f8bSJeremy L Thompson PetscCall(VecDestroy(&F)); 8422b730f8bSJeremy L Thompson PetscCall(VecDestroy(&F_loc)); 8432b730f8bSJeremy L Thompson PetscCall(VecDestroy(&neumann_bcs)); 8442b730f8bSJeremy L Thompson PetscCall(VecDestroy(&bcs_loc)); 8452b730f8bSJeremy L Thompson PetscCall(MatDestroy(&jacob_mat_coarse)); 8462b730f8bSJeremy L Thompson PetscCall(SNESDestroy(&snes)); 8472b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_orig)); 8482b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_energy)); 8492b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_diagnostic)); 8502b730f8bSJeremy L Thompson PetscCall(PetscFree(level_dms)); 851ccaff030SJeremy L Thompson 8525754ecacSJeremy L Thompson // -- Function list 8532b730f8bSJeremy L Thompson PetscCall(PetscFunctionListDestroy(&problem_functions->setupPhysics)); 8542b730f8bSJeremy L Thompson PetscCall(PetscFunctionListDestroy(&problem_functions->setupSmootherPhysics)); 8552b730f8bSJeremy L Thompson PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel)); 8562b730f8bSJeremy L Thompson PetscCall(PetscFunctionListDestroy(&problem_functions->setupLibceedLevel)); 8575754ecacSJeremy L Thompson 858ccaff030SJeremy L Thompson // Structs 8592b730f8bSJeremy L Thompson PetscCall(PetscFree(res_ctx)); 8602b730f8bSJeremy L Thompson PetscCall(PetscFree(form_jacob_ctx)); 8612b730f8bSJeremy L Thompson PetscCall(PetscFree(jacob_coarse_ctx)); 8622b730f8bSJeremy L Thompson PetscCall(PetscFree(app_ctx)); 8632b730f8bSJeremy L Thompson PetscCall(PetscFree(problem_functions)); 8642b730f8bSJeremy L Thompson PetscCall(PetscFree(units)); 865ccaff030SJeremy L Thompson 866ccaff030SJeremy L Thompson return PetscFinalize(); 867ccaff030SJeremy L Thompson } 868