1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 // libCEED + PETSc Example: Elasticity 18 // 19 // This example demonstrates a simple usage of libCEED with PETSc to solve 20 // elasticity problems. 21 // 22 // The code uses higher level communication protocols in DMPlex. 23 // 24 // Build with: 25 // 26 // make elasticity [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 27 // 28 // Sample runs: 29 // 30 // ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -problem Linear -forcing mms 31 // ./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 32 // ./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 33 // 34 // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes 35 // 36 //TESTARGS(name="solids-Linear-MMS") -ceed {ceed_resource} -test -degree 3 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3 37 //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 38 //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 39 40 /// @file 41 /// CEED elasticity example using PETSc with DMPlex 42 43 const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n"; 44 45 #include "elasticity.h" 46 47 int main(int argc, char **argv) { 48 PetscInt ierr; 49 MPI_Comm comm; 50 // Context structs 51 AppCtx app_ctx; // Contains problem options 52 ProblemFunctions problem_functions; // Setup functions for each problem 53 Units units; // Contains units scaling 54 // PETSc objects 55 PetscLogStage stage_dm_setup, stage_libceed_setup, 56 stage_snes_setup, stage_snes_solve; 57 DM dm_orig; // Distributed DM to clone 58 DM dm_energy, dm_diagnostic; // DMs for postprocessing 59 DM *level_dms; 60 Vec U, *U_g, *U_loc; // U: solution, R: residual, F: forcing 61 Vec R, R_loc, F, F_loc; // g: global, loc: local 62 Vec neumann_bcs = NULL, bcs_loc = NULL; 63 SNES snes; 64 Mat *jacob_mat, jacob_mat_coarse, *prolong_restr_mat; 65 // PETSc data 66 UserMult res_ctx, jacob_coarse_ctx = NULL, *jacob_ctx; 67 FormJacobCtx form_jacob_ctx; 68 UserMultProlongRestr *prolong_restr_ctx; 69 PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V; 70 // libCEED objects 71 Ceed ceed; 72 CeedData *ceed_data; 73 CeedQFunctionContext ctx_phys, ctx_phys_smoother = NULL; 74 // Parameters 75 PetscInt num_comp_u = 3; // 3 DoFs in 3D 76 PetscInt num_comp_e = 1, num_comp_d = 5; // 1 energy output, 5 diagnostic 77 PetscInt num_levels = 1, fine_level = 0; 78 PetscInt *U_g_size, *U_l_size, *U_loc_size; 79 PetscInt snes_its = 0, ksp_its = 0; 80 double start_time, elapsed_time, min_time, max_time; 81 82 ierr = PetscInitialize(&argc, &argv, NULL, help); 83 if (ierr) return ierr; 84 85 // --------------------------------------------------------------------------- 86 // Process command line options 87 // --------------------------------------------------------------------------- 88 comm = PETSC_COMM_WORLD; 89 90 // -- Set mesh file, polynomial degree, problem type 91 ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr); 92 ierr = ProcessCommandLineOptions(comm, app_ctx); CHKERRQ(ierr); 93 ierr = PetscCalloc1(1, &problem_functions); CHKERRQ(ierr); 94 ierr = RegisterProblems(problem_functions); CHKERRQ(ierr); 95 num_levels = app_ctx->num_levels; 96 fine_level = num_levels - 1; 97 98 // --------------------------------------------------------------------------- 99 // Initialize libCEED 100 // --------------------------------------------------------------------------- 101 // Initialize backend 102 CeedInit(app_ctx->ceed_resource, &ceed); 103 104 // Check preferred MemType 105 CeedMemType mem_type_backend; 106 CeedGetPreferredMemType(ceed, &mem_type_backend); 107 // Setup physics context and wrap in libCEED object 108 { 109 PetscErrorCode (*SetupPhysics)(MPI_Comm, Ceed, Units *, CeedQFunctionContext *); 110 ierr = PetscFunctionListFind(problem_functions->setupPhysics, app_ctx->name, 111 &SetupPhysics); CHKERRQ(ierr); 112 if (!SetupPhysics) 113 SETERRQ(PETSC_COMM_SELF, 1, "Physics setup for '%s' not found", 114 app_ctx->name); 115 ierr = (*SetupPhysics)(comm, ceed, &units, &ctx_phys); CHKERRQ(ierr); 116 PetscErrorCode (*SetupSmootherPhysics)(MPI_Comm, Ceed, CeedQFunctionContext, 117 CeedQFunctionContext *); 118 ierr = PetscFunctionListFind(problem_functions->setupSmootherPhysics, 119 app_ctx->name, &SetupSmootherPhysics); 120 CHKERRQ(ierr); 121 if (!SetupSmootherPhysics) 122 SETERRQ(PETSC_COMM_SELF, 1, "Smoother physics setup for '%s' not found", 123 app_ctx->name); 124 ierr = (*SetupSmootherPhysics)(comm, ceed, ctx_phys, &ctx_phys_smoother); 125 CHKERRQ(ierr); 126 } 127 128 // --------------------------------------------------------------------------- 129 // Setup DM 130 // --------------------------------------------------------------------------- 131 // Performance logging 132 ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stage_dm_setup); 133 CHKERRQ(ierr); 134 ierr = PetscLogStagePush(stage_dm_setup); CHKERRQ(ierr); 135 136 // -- Create distributed DM from mesh file 137 ierr = CreateDistributedDM(comm, app_ctx, &dm_orig); CHKERRQ(ierr); 138 VecType vectype; 139 switch (mem_type_backend) { 140 case CEED_MEM_HOST: vectype = VECSTANDARD; break; 141 case CEED_MEM_DEVICE: { 142 const char *resolved; 143 CeedGetResource(ceed, &resolved); 144 if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 145 else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 146 else vectype = VECSTANDARD; 147 } 148 } 149 ierr = DMSetVecType(dm_orig, vectype); CHKERRQ(ierr); 150 ierr = DMPlexDistributeSetDefault(dm_orig, PETSC_FALSE); CHKERRQ(ierr); 151 ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr); 152 153 // -- Setup DM by polynomial degree 154 ierr = PetscMalloc1(num_levels, &level_dms); CHKERRQ(ierr); 155 for (PetscInt level = 0; level < num_levels; level++) { 156 ierr = DMClone(dm_orig, &level_dms[level]); CHKERRQ(ierr); 157 ierr = DMGetVecType(dm_orig, &vectype); CHKERRQ(ierr); 158 ierr = DMSetVecType(level_dms[level], vectype); CHKERRQ(ierr); 159 ierr = SetupDMByDegree(level_dms[level], app_ctx, app_ctx->level_degrees[level], 160 PETSC_TRUE, num_comp_u); CHKERRQ(ierr); 161 // -- Label field components for viewing 162 // Empty name for conserved field (because there is only one field) 163 PetscSection section; 164 ierr = DMGetLocalSection(level_dms[level], §ion); CHKERRQ(ierr); 165 ierr = PetscSectionSetFieldName(section, 0, "Displacement"); CHKERRQ(ierr); 166 ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 167 CHKERRQ(ierr); 168 ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 169 CHKERRQ(ierr); 170 ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 171 CHKERRQ(ierr); 172 } 173 174 // -- Setup postprocessing DMs 175 ierr = DMClone(dm_orig, &dm_energy); CHKERRQ(ierr); 176 ierr = SetupDMByDegree(dm_energy, app_ctx, app_ctx->level_degrees[fine_level], 177 PETSC_FALSE, num_comp_e); CHKERRQ(ierr); 178 ierr = DMClone(dm_orig, &dm_diagnostic); CHKERRQ(ierr); 179 ierr = SetupDMByDegree(dm_diagnostic, app_ctx, 180 app_ctx->level_degrees[fine_level], 181 PETSC_FALSE, num_comp_u + num_comp_d); CHKERRQ(ierr); 182 ierr = DMSetVecType(dm_energy, vectype); CHKERRQ(ierr); 183 ierr = DMSetVecType(dm_diagnostic, vectype); CHKERRQ(ierr); 184 { 185 // -- Label field components for viewing 186 // Empty name for conserved field (because there is only one field) 187 PetscSection section; 188 ierr = DMGetLocalSection(dm_diagnostic, §ion); CHKERRQ(ierr); 189 ierr = PetscSectionSetFieldName(section, 0, "Diagnostics"); CHKERRQ(ierr); 190 ierr = PetscSectionSetComponentName(section, 0, 0, "DisplacementX"); 191 CHKERRQ(ierr); 192 ierr = PetscSectionSetComponentName(section, 0, 1, "DisplacementY"); 193 CHKERRQ(ierr); 194 ierr = PetscSectionSetComponentName(section, 0, 2, "DisplacementZ"); 195 CHKERRQ(ierr); 196 ierr = PetscSectionSetComponentName(section, 0, 3, "Pressure"); 197 CHKERRQ(ierr); 198 ierr = PetscSectionSetComponentName(section, 0, 4, "VolumentricStrain"); 199 CHKERRQ(ierr); 200 ierr = PetscSectionSetComponentName(section, 0, 5, "TraceE2"); 201 CHKERRQ(ierr); 202 ierr = PetscSectionSetComponentName(section, 0, 6, "detJ"); 203 CHKERRQ(ierr); 204 ierr = PetscSectionSetComponentName(section, 0, 7, "StrainEnergyDensity"); 205 CHKERRQ(ierr); 206 } 207 208 // --------------------------------------------------------------------------- 209 // Setup solution and work vectors 210 // --------------------------------------------------------------------------- 211 // Allocate arrays 212 ierr = PetscMalloc1(num_levels, &U_g); CHKERRQ(ierr); 213 ierr = PetscMalloc1(num_levels, &U_loc); CHKERRQ(ierr); 214 ierr = PetscMalloc1(num_levels, &U_g_size); CHKERRQ(ierr); 215 ierr = PetscMalloc1(num_levels, &U_l_size); CHKERRQ(ierr); 216 ierr = PetscMalloc1(num_levels, &U_loc_size); CHKERRQ(ierr); 217 218 // -- Setup solution vectors for each level 219 for (PetscInt level = 0; level < num_levels; level++) { 220 // -- Create global unknown vector U 221 ierr = DMCreateGlobalVector(level_dms[level], &U_g[level]); CHKERRQ(ierr); 222 ierr = VecGetSize(U_g[level], &U_g_size[level]); CHKERRQ(ierr); 223 // Note: Local size for matShell 224 ierr = VecGetLocalSize(U_g[level], &U_l_size[level]); CHKERRQ(ierr); 225 226 // -- Create local unknown vector U_loc 227 ierr = DMCreateLocalVector(level_dms[level], &U_loc[level]); CHKERRQ(ierr); 228 // Note: local size for libCEED 229 ierr = VecGetSize(U_loc[level], &U_loc_size[level]); CHKERRQ(ierr); 230 } 231 232 // -- Create residual and forcing vectors 233 ierr = VecDuplicate(U_g[fine_level], &U); CHKERRQ(ierr); 234 ierr = VecDuplicate(U_g[fine_level], &R); CHKERRQ(ierr); 235 ierr = VecDuplicate(U_g[fine_level], &F); CHKERRQ(ierr); 236 ierr = VecDuplicate(U_loc[fine_level], &R_loc); CHKERRQ(ierr); 237 ierr = VecDuplicate(U_loc[fine_level], &F_loc); CHKERRQ(ierr); 238 239 // Performance logging 240 ierr = PetscLogStagePop(); 241 242 // --------------------------------------------------------------------------- 243 // Set up libCEED 244 // --------------------------------------------------------------------------- 245 // Performance logging 246 ierr = PetscLogStageRegister("libCEED Setup Stage", &stage_libceed_setup); 247 CHKERRQ(ierr); 248 ierr = PetscLogStagePush(stage_libceed_setup); CHKERRQ(ierr); 249 250 // -- Create libCEED local forcing vector 251 CeedVector force_ceed; 252 CeedScalar *f; 253 PetscMemType force_mem_type; 254 if (app_ctx->forcing_choice != FORCE_NONE) { 255 ierr = VecGetArrayAndMemType(F_loc, &f, &force_mem_type); CHKERRQ(ierr); 256 CeedVectorCreate(ceed, U_loc_size[fine_level], &force_ceed); 257 CeedVectorSetArray(force_ceed, MemTypeP2C(force_mem_type), CEED_USE_POINTER, f); 258 } 259 260 // -- Create libCEED local Neumann BCs vector 261 CeedVector neumann_ceed; 262 CeedScalar *n; 263 PetscMemType nummann_mem_type; 264 if (app_ctx->bc_traction_count > 0) { 265 ierr = VecDuplicate(U, &neumann_bcs); CHKERRQ(ierr); 266 ierr = VecDuplicate(U_loc[fine_level], &bcs_loc); CHKERRQ(ierr); 267 ierr = VecGetArrayAndMemType(bcs_loc, &n, &nummann_mem_type); CHKERRQ(ierr); 268 CeedVectorCreate(ceed, U_loc_size[fine_level], &neumann_ceed); 269 CeedVectorSetArray(neumann_ceed, MemTypeP2C(nummann_mem_type), 270 CEED_USE_POINTER, n); 271 } 272 273 // -- Setup libCEED objects 274 ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr); 275 // ---- Setup residual, Jacobian evaluator and geometric information 276 ierr = PetscCalloc1(1, &ceed_data[fine_level]); CHKERRQ(ierr); 277 { 278 PetscErrorCode (*SetupLibceedFineLevel)(DM, DM, DM, Ceed, AppCtx, 279 CeedQFunctionContext, PetscInt, 280 PetscInt, PetscInt, PetscInt, 281 CeedVector, CeedVector, CeedData *); 282 ierr = PetscFunctionListFind(problem_functions->setupLibceedFineLevel, 283 app_ctx->name, &SetupLibceedFineLevel); 284 CHKERRQ(ierr); 285 if (!SetupLibceedFineLevel) 286 SETERRQ(PETSC_COMM_SELF, 1, "Fine grid setup for '%s' not found", 287 app_ctx->name); 288 ierr = (*SetupLibceedFineLevel)(level_dms[fine_level], dm_energy, dm_diagnostic, 289 ceed, app_ctx, ctx_phys, fine_level, 290 num_comp_u, U_g_size[fine_level], 291 U_loc_size[fine_level], 292 force_ceed, neumann_ceed, ceed_data); 293 CHKERRQ(ierr); 294 } 295 // ---- Setup coarse Jacobian evaluator and prolongation/restriction 296 for (PetscInt level = num_levels - 2; level >= 0; level--) { 297 ierr = PetscCalloc1(1, &ceed_data[level]); CHKERRQ(ierr); 298 299 // Get global communication restriction 300 ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr); 301 ierr = VecSet(U_loc[level+1], 1.0); CHKERRQ(ierr); 302 ierr = DMLocalToGlobal(level_dms[level+1], U_loc[level+1], ADD_VALUES, 303 U_g[level+1]); CHKERRQ(ierr); 304 ierr = DMGlobalToLocal(level_dms[level+1], U_g[level+1], INSERT_VALUES, 305 U_loc[level+1]); CHKERRQ(ierr); 306 307 // Place in libCEED array 308 const PetscScalar *m; 309 PetscMemType m_mem_type; 310 ierr = VecGetArrayReadAndMemType(U_loc[level+1], &m, &m_mem_type); 311 CHKERRQ(ierr); 312 CeedVectorSetArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type), 313 CEED_USE_POINTER, (CeedScalar *)m); 314 315 // Note: use high order ceed, if specified and degree > 4 316 PetscErrorCode (*SetupLibceedLevel)(DM, Ceed, AppCtx, PetscInt, 317 PetscInt, PetscInt, PetscInt, CeedVector, CeedData *); 318 ierr = PetscFunctionListFind(problem_functions->setupLibceedLevel, 319 app_ctx->name, &SetupLibceedLevel); 320 CHKERRQ(ierr); 321 if (!SetupLibceedLevel) 322 SETERRQ(PETSC_COMM_SELF, 1, "Coarse grid setup for '%s' not found", 323 app_ctx->name); 324 ierr = (*SetupLibceedLevel)(level_dms[level], ceed, app_ctx, 325 level, num_comp_u, U_g_size[level], 326 U_loc_size[level], ceed_data[level+1]->x_ceed, 327 ceed_data); 328 CHKERRQ(ierr); 329 330 // Restore PETSc vector 331 CeedVectorTakeArray(ceed_data[level+1]->x_ceed, MemTypeP2C(m_mem_type), 332 (CeedScalar **)&m); 333 ierr = VecRestoreArrayReadAndMemType(U_loc[level+1], &m); CHKERRQ(ierr); 334 ierr = VecZeroEntries(U_g[level+1]); CHKERRQ(ierr); 335 ierr = VecZeroEntries(U_loc[level+1]); CHKERRQ(ierr); 336 } 337 338 // Performance logging 339 ierr = PetscLogStagePop(); 340 341 // --------------------------------------------------------------------------- 342 // Setup global forcing and Neumann BC vectors 343 // --------------------------------------------------------------------------- 344 ierr = VecZeroEntries(F); CHKERRQ(ierr); 345 346 if (app_ctx->forcing_choice != FORCE_NONE) { 347 CeedVectorTakeArray(force_ceed, MemTypeP2C(force_mem_type), NULL); 348 ierr = VecRestoreArrayAndMemType(F_loc, &f); CHKERRQ(ierr); 349 ierr = DMLocalToGlobal(level_dms[fine_level], F_loc, ADD_VALUES, F); 350 CHKERRQ(ierr); 351 CeedVectorDestroy(&force_ceed); 352 } 353 354 if (app_ctx->bc_traction_count > 0) { 355 ierr = VecZeroEntries(neumann_bcs); CHKERRQ(ierr); 356 CeedVectorTakeArray(neumann_ceed, MemTypeP2C(nummann_mem_type), NULL); 357 ierr = VecRestoreArrayAndMemType(bcs_loc, &n); CHKERRQ(ierr); 358 ierr = DMLocalToGlobal(level_dms[fine_level], bcs_loc, ADD_VALUES, neumann_bcs); 359 CHKERRQ(ierr); 360 CeedVectorDestroy(&neumann_ceed); 361 } 362 363 // --------------------------------------------------------------------------- 364 // Print problem summary 365 // --------------------------------------------------------------------------- 366 if (!app_ctx->test_mode) { 367 const char *usedresource; 368 CeedGetResource(ceed, &usedresource); 369 char hostname[PETSC_MAX_PATH_LEN]; 370 ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr); 371 PetscInt comm_size; 372 ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 373 374 ierr = PetscPrintf(comm, 375 "\n-- Elasticity Example - libCEED + PETSc --\n" 376 " MPI:\n" 377 " Hostname : %s\n" 378 " Total ranks : %d\n" 379 " libCEED:\n" 380 " libCEED Backend : %s\n" 381 " libCEED Backend MemType : %s\n", 382 hostname, comm_size, usedresource, CeedMemTypes[mem_type_backend]); 383 CHKERRQ(ierr); 384 385 VecType vecType; 386 ierr = VecGetType(U, &vecType); CHKERRQ(ierr); 387 ierr = PetscPrintf(comm, 388 " PETSc:\n" 389 " PETSc Vec Type : %s\n", 390 vecType); CHKERRQ(ierr); 391 392 ierr = PetscPrintf(comm, 393 " Problem:\n" 394 " Problem Name : %s\n" 395 " Forcing Function : %s\n" 396 " Mesh:\n" 397 " File : %s\n" 398 " Number of 1D Basis Nodes (p) : %d\n" 399 " Number of 1D Quadrature Points (q) : %d\n" 400 " Global nodes : %D\n" 401 " Owned nodes : %D\n" 402 " DoF per node : %D\n" 403 " Multigrid:\n" 404 " Type : %s\n" 405 " Number of Levels : %d\n", 406 app_ctx->name_for_disp, 407 forcing_types_for_disp[app_ctx->forcing_choice], 408 app_ctx->mesh_file[0] ? app_ctx->mesh_file : "Box Mesh", 409 app_ctx->degree + 1, app_ctx->degree + 1, 410 U_g_size[fine_level]/num_comp_u, U_l_size[fine_level]/num_comp_u, 411 num_comp_u, 412 (app_ctx->degree == 1 && 413 app_ctx->multigrid_choice != MULTIGRID_NONE) ? 414 "Algebraic multigrid" : 415 multigrid_types_for_disp[app_ctx->multigrid_choice], 416 (app_ctx->degree == 1 || 417 app_ctx->multigrid_choice == MULTIGRID_NONE) ? 418 0 : num_levels); CHKERRQ(ierr); 419 420 if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 421 for (PetscInt i = 0; i < 2; i++) { 422 CeedInt level = i ? fine_level : 0; 423 ierr = PetscPrintf(comm, 424 " Level %D (%s):\n" 425 " Number of 1D Basis Nodes (p) : %d\n" 426 " Global Nodes : %D\n" 427 " Owned Nodes : %D\n", 428 level, i ? "fine" : "coarse", 429 app_ctx->level_degrees[level] + 1, 430 U_g_size[level]/num_comp_u, U_l_size[level]/num_comp_u); 431 CHKERRQ(ierr); 432 } 433 } 434 } 435 436 // --------------------------------------------------------------------------- 437 // Setup SNES 438 // --------------------------------------------------------------------------- 439 // Performance logging 440 ierr = PetscLogStageRegister("SNES Setup Stage", &stage_snes_setup); 441 CHKERRQ(ierr); 442 ierr = PetscLogStagePush(stage_snes_setup); CHKERRQ(ierr); 443 444 // Create SNES 445 ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 446 ierr = SNESSetDM(snes, level_dms[fine_level]); CHKERRQ(ierr); 447 448 // -- Jacobian evaluators 449 ierr = PetscMalloc1(num_levels, &jacob_ctx); CHKERRQ(ierr); 450 ierr = PetscMalloc1(num_levels, &jacob_mat); CHKERRQ(ierr); 451 for (PetscInt level = 0; level < num_levels; level++) { 452 // -- Jacobian context for level 453 ierr = PetscMalloc1(1, &jacob_ctx[level]); CHKERRQ(ierr); 454 ierr = SetupJacobianCtx(comm, app_ctx, level_dms[level], U_g[level], 455 U_loc[level], ceed_data[level], ceed, ctx_phys, 456 ctx_phys_smoother, jacob_ctx[level]); CHKERRQ(ierr); 457 458 // -- Form Action of Jacobian on delta_u 459 ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level], U_g_size[level], 460 U_g_size[level], jacob_ctx[level], &jacob_mat[level]); 461 CHKERRQ(ierr); 462 ierr = MatShellSetOperation(jacob_mat[level], MATOP_MULT, 463 (void (*)(void))ApplyJacobian_Ceed); 464 CHKERRQ(ierr); 465 ierr = MatShellSetOperation(jacob_mat[level], MATOP_GET_DIAGONAL, 466 (void(*)(void))GetDiag_Ceed); 467 ierr = MatShellSetVecType(jacob_mat[level], vectype); CHKERRQ(ierr); 468 } 469 // Note: FormJacobian updates Jacobian matrices on each level 470 // and assembles the Jpre matrix, if needed 471 ierr = PetscMalloc1(1, &form_jacob_ctx); CHKERRQ(ierr); 472 form_jacob_ctx->jacob_ctx = jacob_ctx; 473 form_jacob_ctx->num_levels = num_levels; 474 form_jacob_ctx->jacob_mat = jacob_mat; 475 476 // -- Residual evaluation function 477 ierr = PetscCalloc1(1, &res_ctx); CHKERRQ(ierr); 478 ierr = PetscMemcpy(res_ctx, jacob_ctx[fine_level], 479 sizeof(*jacob_ctx[fine_level])); CHKERRQ(ierr); 480 res_ctx->op = ceed_data[fine_level]->op_residual; 481 res_ctx->qf = ceed_data[fine_level]->qf_residual; 482 if (app_ctx->bc_traction_count > 0) 483 res_ctx->neumann_bcs = neumann_bcs; 484 else 485 res_ctx->neumann_bcs = NULL; 486 ierr = SNESSetFunction(snes, R, FormResidual_Ceed, res_ctx); CHKERRQ(ierr); 487 488 // -- Prolongation/Restriction evaluation 489 ierr = PetscMalloc1(num_levels, &prolong_restr_ctx); CHKERRQ(ierr); 490 ierr = PetscMalloc1(num_levels, &prolong_restr_mat); CHKERRQ(ierr); 491 for (PetscInt level = 1; level < num_levels; level++) { 492 // ---- Prolongation/restriction context for level 493 ierr = PetscMalloc1(1, &prolong_restr_ctx[level]); CHKERRQ(ierr); 494 ierr = SetupProlongRestrictCtx(comm, app_ctx, level_dms[level-1], 495 level_dms[level], U_g[level], U_loc[level-1], 496 U_loc[level], ceed_data[level-1], 497 ceed_data[level], ceed, 498 prolong_restr_ctx[level]); CHKERRQ(ierr); 499 500 // ---- Form Action of Jacobian on delta_u 501 ierr = MatCreateShell(comm, U_l_size[level], U_l_size[level-1], U_g_size[level], 502 U_g_size[level-1], prolong_restr_ctx[level], 503 &prolong_restr_mat[level]); CHKERRQ(ierr); 504 // Note: In PCMG, restriction is the transpose of prolongation 505 ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT, 506 (void (*)(void))Prolong_Ceed); 507 ierr = MatShellSetOperation(prolong_restr_mat[level], MATOP_MULT_TRANSPOSE, 508 (void (*)(void))Restrict_Ceed); 509 CHKERRQ(ierr); 510 ierr = MatShellSetVecType(prolong_restr_mat[level], vectype); CHKERRQ(ierr); 511 } 512 513 // --------------------------------------------------------------------------- 514 // Setup for AMG coarse solve 515 // --------------------------------------------------------------------------- 516 if (app_ctx->multigrid_choice != MULTIGRID_NONE) { 517 // -- Jacobian Matrix 518 ierr = DMCreateMatrix(level_dms[0], &jacob_mat_coarse); CHKERRQ(ierr); 519 520 if (app_ctx->degree > 1) { 521 // -- Assemble sparsity pattern 522 CeedInt num_entries, *rows, *cols; 523 CeedVector coo_values; 524 CeedOperatorLinearAssembleSymbolic(ceed_data[0]->op_jacobian, &num_entries, 525 &rows, &cols); 526 ISLocalToGlobalMapping ltog_row, ltog_col; 527 ierr = MatGetLocalToGlobalMapping(jacob_mat_coarse, <og_row, <og_col); 528 CHKERRQ(ierr); 529 ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows); 530 CHKERRQ(ierr); 531 ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols); 532 CHKERRQ(ierr); 533 ierr = MatSetPreallocationCOO(jacob_mat_coarse, num_entries, rows, cols); 534 CHKERRQ(ierr); 535 free(rows); 536 free(cols); 537 CeedVectorCreate(ceed, num_entries, &coo_values); 538 539 // -- Update form_jacob_ctx 540 form_jacob_ctx->coo_values = coo_values; 541 form_jacob_ctx->op_coarse = ceed_data[0]->op_jacobian; 542 form_jacob_ctx->jacob_mat_coarse = jacob_mat_coarse; 543 } 544 } 545 546 // Set Jacobian function 547 if (app_ctx->degree > 1) { 548 ierr = SNESSetJacobian(snes, jacob_mat[fine_level], jacob_mat[fine_level], 549 FormJacobian, form_jacob_ctx); CHKERRQ(ierr); 550 } else { 551 ierr = SNESSetJacobian(snes, jacob_mat[0], jacob_mat_coarse, 552 SNESComputeJacobianDefaultColor, NULL); 553 CHKERRQ(ierr); 554 } 555 556 // --------------------------------------------------------------------------- 557 // Setup KSP 558 // --------------------------------------------------------------------------- 559 { 560 PC pc; 561 KSP ksp; 562 563 // -- KSP 564 ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 565 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 566 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 567 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 568 PETSC_DEFAULT); CHKERRQ(ierr); 569 ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 570 571 // -- Preconditioning 572 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 573 ierr = PCSetDM(pc, level_dms[fine_level]); CHKERRQ(ierr); 574 ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 575 576 if (app_ctx->multigrid_choice == MULTIGRID_NONE) { 577 // ---- No Multigrid 578 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 579 ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 580 } else if (app_ctx->degree == 1) { 581 // ---- AMG for degree 1 582 ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 583 } else { 584 // ---- PCMG 585 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 586 587 // ------ PCMG levels 588 ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr); 589 for (PetscInt level = 0; level < num_levels; level++) { 590 // -------- Smoother 591 KSP ksp_smoother, ksp_est; 592 PC pc_smoother; 593 594 // ---------- Smoother KSP 595 ierr = PCMGGetSmoother(pc, level, &ksp_smoother); CHKERRQ(ierr); 596 ierr = KSPSetDM(ksp_smoother, level_dms[level]); CHKERRQ(ierr); 597 ierr = KSPSetDMActive(ksp_smoother, PETSC_FALSE); CHKERRQ(ierr); 598 599 // ---------- Chebyshev options 600 ierr = KSPSetType(ksp_smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 601 ierr = KSPChebyshevEstEigSet(ksp_smoother, 0, 0.1, 0, 1.1); 602 CHKERRQ(ierr); 603 ierr = KSPChebyshevEstEigGetKSP(ksp_smoother, &ksp_est); CHKERRQ(ierr); 604 ierr = KSPSetType(ksp_est, KSPCG); CHKERRQ(ierr); 605 ierr = KSPChebyshevEstEigSetUseNoisy(ksp_smoother, PETSC_TRUE); 606 CHKERRQ(ierr); 607 ierr = KSPSetOperators(ksp_smoother, jacob_mat[level], jacob_mat[level]); 608 CHKERRQ(ierr); 609 610 // ---------- Smoother preconditioner 611 ierr = KSPGetPC(ksp_smoother, &pc_smoother); CHKERRQ(ierr); 612 ierr = PCSetType(pc_smoother, PCJACOBI); CHKERRQ(ierr); 613 ierr = PCJacobiSetType(pc_smoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 614 615 // -------- Work vector 616 if (level != fine_level) { 617 ierr = PCMGSetX(pc, level, U_g[level]); CHKERRQ(ierr); 618 } 619 620 // -------- Level prolongation/restriction operator 621 if (level > 0) { 622 ierr = PCMGSetInterpolation(pc, level, prolong_restr_mat[level]); 623 CHKERRQ(ierr); 624 ierr = PCMGSetRestriction(pc, level, prolong_restr_mat[level]); 625 CHKERRQ(ierr); 626 } 627 } 628 629 // ------ PCMG coarse solve 630 KSP ksp_coarse; 631 PC pc_coarse; 632 633 // -------- Coarse KSP 634 ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr); 635 ierr = KSPSetType(ksp_coarse, KSPPREONLY); CHKERRQ(ierr); 636 ierr = KSPSetOperators(ksp_coarse, jacob_mat_coarse, jacob_mat_coarse); 637 CHKERRQ(ierr); 638 ierr = KSPSetOptionsPrefix(ksp_coarse, "coarse_"); CHKERRQ(ierr); 639 640 // -------- Coarse preconditioner 641 ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr); 642 ierr = PCSetType(pc_coarse, PCGAMG); CHKERRQ(ierr); 643 ierr = PCSetOptionsPrefix(pc_coarse, "coarse_"); CHKERRQ(ierr); 644 645 ierr = KSPSetFromOptions(ksp_coarse); CHKERRQ(ierr); 646 ierr = PCSetFromOptions(pc_coarse); CHKERRQ(ierr); 647 648 // ------ PCMG options 649 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 650 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 651 ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr); 652 } 653 ierr = KSPSetFromOptions(ksp); 654 ierr = PCSetFromOptions(pc); 655 } 656 { 657 // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 658 SNESLineSearch line_search; 659 660 ierr = SNESGetLineSearch(snes, &line_search); CHKERRQ(ierr); 661 ierr = SNESLineSearchSetType(line_search, SNESLINESEARCHCP); CHKERRQ(ierr); 662 } 663 664 ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 665 666 // Performance logging 667 ierr = PetscLogStagePop(); 668 669 // --------------------------------------------------------------------------- 670 // Set initial guess 671 // --------------------------------------------------------------------------- 672 ierr = PetscObjectSetName((PetscObject)U, ""); CHKERRQ(ierr); 673 ierr = VecSet(U, 0.0); CHKERRQ(ierr); 674 675 // View solution 676 if (app_ctx->view_soln) { 677 ierr = ViewSolution(comm, app_ctx, U, 0, 0.0); CHKERRQ(ierr); 678 } 679 680 // --------------------------------------------------------------------------- 681 // Solve SNES 682 // --------------------------------------------------------------------------- 683 PetscBool snes_monitor = PETSC_FALSE; 684 ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snes_monitor); 685 CHKERRQ(ierr); 686 687 // Performance logging 688 ierr = PetscLogStageRegister("SNES Solve Stage", &stage_snes_solve); 689 CHKERRQ(ierr); 690 ierr = PetscLogStagePush(stage_snes_solve); CHKERRQ(ierr); 691 692 // Timing 693 ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 694 start_time = MPI_Wtime(); 695 696 // Solve for each load increment 697 PetscInt increment; 698 for (increment = 1; increment <= app_ctx->num_increments; increment++) { 699 // -- Log increment count 700 if (snes_monitor) { 701 ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 702 CHKERRQ(ierr); 703 } 704 705 // -- Scale the problem 706 PetscScalar load_increment = 1.0*increment / app_ctx->num_increments, 707 scalingFactor = load_increment / 708 (increment == 1 ? 1 : res_ctx->load_increment); 709 res_ctx->load_increment = load_increment; 710 if (app_ctx->num_increments > 1 && app_ctx->forcing_choice != FORCE_NONE) { 711 ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 712 } 713 714 // -- Solve 715 ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 716 717 // -- View solution 718 if (app_ctx->view_soln) { 719 ierr = ViewSolution(comm, app_ctx, U, increment, load_increment); CHKERRQ(ierr); 720 } 721 722 // -- Update SNES iteration count 723 PetscInt its; 724 ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 725 snes_its += its; 726 ierr = SNESGetLinearSolveIterations(snes, &its); CHKERRQ(ierr); 727 ksp_its += its; 728 729 // -- Check for divergence 730 SNESConvergedReason reason; 731 ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 732 if (reason < 0) 733 break; 734 if (app_ctx->energy_viewer) { 735 // -- Log strain energy for current load increment 736 CeedScalar energy; 737 ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, 738 U, &energy); CHKERRQ(ierr); 739 740 if (!app_ctx->test_mode) { 741 // -- Output 742 ierr = PetscPrintf(comm, 743 " Strain Energy : %.12e\n", 744 energy); CHKERRQ(ierr); 745 } 746 ierr = PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", load_increment, 747 energy); CHKERRQ(ierr); 748 } 749 } 750 751 // Timing 752 elapsed_time = MPI_Wtime() - start_time; 753 754 // Performance logging 755 ierr = PetscLogStagePop(); 756 757 // --------------------------------------------------------------------------- 758 // Output summary 759 // --------------------------------------------------------------------------- 760 if (!app_ctx->test_mode) { 761 // -- SNES 762 SNESType snes_type; 763 SNESConvergedReason reason; 764 PetscReal rnorm; 765 ierr = SNESGetType(snes, &snes_type); CHKERRQ(ierr); 766 ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 767 ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 768 ierr = PetscPrintf(comm, 769 " SNES:\n" 770 " SNES Type : %s\n" 771 " SNES Convergence : %s\n" 772 " Number of Load Increments : %d\n" 773 " Completed Load Increments : %d\n" 774 " Total SNES Iterations : %D\n" 775 " Final rnorm : %e\n", 776 snes_type, SNESConvergedReasons[reason], 777 app_ctx->num_increments, increment - 1, 778 snes_its, (double)rnorm); CHKERRQ(ierr); 779 780 // -- KSP 781 KSP ksp; 782 KSPType ksp_type; 783 ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 784 ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 785 ierr = PetscPrintf(comm, 786 " Linear Solver:\n" 787 " KSP Type : %s\n" 788 " Total KSP Iterations : %D\n", 789 ksp_type, ksp_its); CHKERRQ(ierr); 790 791 // -- PC 792 PC pc; 793 PCType pc_type; 794 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 795 ierr = PCGetType(pc, &pc_type); CHKERRQ(ierr); 796 ierr = PetscPrintf(comm, 797 " PC Type : %s\n", 798 pc_type); CHKERRQ(ierr); 799 800 if (!strcmp(pc_type, PCMG)) { 801 PCMGType pcmg_type; 802 ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr); 803 ierr = PetscPrintf(comm, 804 " P-Multigrid:\n" 805 " PCMG Type : %s\n" 806 " PCMG Cycle Type : %s\n", 807 PCMGTypes[pcmg_type], 808 PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr); 809 810 // -- Coarse Solve 811 KSP ksp_coarse; 812 PC pc_coarse; 813 PCType pc_type; 814 815 ierr = PCMGGetCoarseSolve(pc, &ksp_coarse); CHKERRQ(ierr); 816 ierr = KSPGetType(ksp_coarse, &ksp_type); CHKERRQ(ierr); 817 ierr = KSPGetPC(ksp_coarse, &pc_coarse); CHKERRQ(ierr); 818 ierr = PCGetType(pc_coarse, &pc_type); CHKERRQ(ierr); 819 ierr = PetscPrintf(comm, 820 " Coarse Solve:\n" 821 " KSP Type : %s\n" 822 " PC Type : %s\n", 823 ksp_type, pc_type); CHKERRQ(ierr); 824 } 825 } 826 827 // --------------------------------------------------------------------------- 828 // Compute solve time 829 // --------------------------------------------------------------------------- 830 if (!app_ctx->test_mode) { 831 ierr = MPI_Allreduce(&elapsed_time, &min_time, 1, MPI_DOUBLE, MPI_MIN, comm); 832 CHKERRQ(ierr); 833 ierr = MPI_Allreduce(&elapsed_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, comm); 834 CHKERRQ(ierr); 835 ierr = PetscPrintf(comm, 836 " Performance:\n" 837 " SNES Solve Time : %g (%g) sec\n" 838 " DoFs/Sec in SNES : %g (%g) million\n", 839 max_time, min_time, 1e-6*U_g_size[fine_level]*ksp_its/max_time, 840 1e-6*U_g_size[fine_level]*ksp_its/min_time); CHKERRQ(ierr); 841 } 842 843 // --------------------------------------------------------------------------- 844 // Compute error 845 // --------------------------------------------------------------------------- 846 if (app_ctx->forcing_choice == FORCE_MMS) { 847 CeedScalar l2_error = 1., l2_U_norm = 1.; 848 const CeedScalar *true_array; 849 Vec error_vec, true_vec; 850 851 // -- Work vectors 852 ierr = VecDuplicate(U, &error_vec); CHKERRQ(ierr); 853 ierr = VecSet(error_vec, 0.0); CHKERRQ(ierr); 854 ierr = VecDuplicate(U, &true_vec); CHKERRQ(ierr); 855 ierr = VecSet(true_vec, 0.0); CHKERRQ(ierr); 856 857 // -- Assemble global true solution vector 858 CeedVectorGetArrayRead(ceed_data[fine_level]->true_soln, 859 CEED_MEM_HOST, &true_array); 860 ierr = VecPlaceArray(res_ctx->Y_loc, (PetscScalar *)true_array); 861 CHKERRQ(ierr); 862 ierr = DMLocalToGlobal(res_ctx->dm, res_ctx->Y_loc, INSERT_VALUES, true_vec); 863 CHKERRQ(ierr); 864 ierr = VecResetArray(res_ctx->Y_loc); CHKERRQ(ierr); 865 CeedVectorRestoreArrayRead(ceed_data[fine_level]->true_soln, &true_array); 866 867 // -- Compute L2 error 868 ierr = VecWAXPY(error_vec, -1.0, U, true_vec); CHKERRQ(ierr); 869 ierr = VecNorm(error_vec, NORM_2, &l2_error); CHKERRQ(ierr); 870 ierr = VecNorm(U, NORM_2, &l2_U_norm); CHKERRQ(ierr); 871 l2_error /= l2_U_norm; 872 873 // -- Output 874 if (!app_ctx->test_mode || l2_error > 0.05) { 875 ierr = PetscPrintf(comm, 876 " L2 Error : %e\n", 877 l2_error); CHKERRQ(ierr); 878 } 879 880 // -- Cleanup 881 ierr = VecDestroy(&error_vec); CHKERRQ(ierr); 882 ierr = VecDestroy(&true_vec); CHKERRQ(ierr); 883 } 884 885 // --------------------------------------------------------------------------- 886 // Compute energy 887 // --------------------------------------------------------------------------- 888 PetscReal energy; 889 ierr = ComputeStrainEnergy(dm_energy, res_ctx, ceed_data[fine_level]->op_energy, 890 U, &energy); CHKERRQ(ierr); 891 if (!app_ctx->test_mode) { 892 // -- Output 893 ierr = PetscPrintf(comm, 894 " Strain Energy : %.12e\n", 895 energy); CHKERRQ(ierr); 896 } 897 ierr = RegressionTests_solids(app_ctx, energy); CHKERRQ(ierr); 898 899 // --------------------------------------------------------------------------- 900 // Output diagnostic quantities 901 // --------------------------------------------------------------------------- 902 if (app_ctx->view_soln || app_ctx->view_final_soln) { 903 // -- Setup context 904 UserMult diagnostic_ctx; 905 ierr = PetscMalloc1(1, &diagnostic_ctx); CHKERRQ(ierr); 906 ierr = PetscMemcpy(diagnostic_ctx, res_ctx, sizeof(*res_ctx)); CHKERRQ(ierr); 907 diagnostic_ctx->dm = dm_diagnostic; 908 diagnostic_ctx->op = ceed_data[fine_level]->op_diagnostic; 909 910 // -- Compute and output 911 ierr = ViewDiagnosticQuantities(comm, level_dms[fine_level], diagnostic_ctx, 912 app_ctx, U, 913 ceed_data[fine_level]->elem_restr_diagnostic); 914 CHKERRQ(ierr); 915 916 // -- Cleanup 917 ierr = PetscFree(diagnostic_ctx); CHKERRQ(ierr); 918 } 919 920 // --------------------------------------------------------------------------- 921 // Free objects 922 // --------------------------------------------------------------------------- 923 // Data in arrays per level 924 for (PetscInt level = 0; level < num_levels; level++) { 925 // Vectors 926 ierr = VecDestroy(&U_g[level]); CHKERRQ(ierr); 927 ierr = VecDestroy(&U_loc[level]); CHKERRQ(ierr); 928 929 // Jacobian matrix and data 930 ierr = VecDestroy(&jacob_ctx[level]->Y_loc); CHKERRQ(ierr); 931 ierr = MatDestroy(&jacob_mat[level]); CHKERRQ(ierr); 932 ierr = PetscFree(jacob_ctx[level]); CHKERRQ(ierr); 933 934 // Prolongation/Restriction matrix and data 935 if (level > 0) { 936 ierr = PetscFree(prolong_restr_ctx[level]); CHKERRQ(ierr); 937 ierr = MatDestroy(&prolong_restr_mat[level]); CHKERRQ(ierr); 938 } 939 940 // DM 941 ierr = DMDestroy(&level_dms[level]); CHKERRQ(ierr); 942 943 // libCEED objects 944 ierr = CeedDataDestroy(level, ceed_data[level]); CHKERRQ(ierr); 945 } 946 947 ierr = PetscViewerDestroy(&app_ctx->energy_viewer); CHKERRQ(ierr); 948 949 // Arrays 950 ierr = PetscFree(U_g); CHKERRQ(ierr); 951 ierr = PetscFree(U_loc); CHKERRQ(ierr); 952 ierr = PetscFree(U_g_size); CHKERRQ(ierr); 953 ierr = PetscFree(U_l_size); CHKERRQ(ierr); 954 ierr = PetscFree(U_loc_size); CHKERRQ(ierr); 955 ierr = PetscFree(jacob_ctx); CHKERRQ(ierr); 956 ierr = PetscFree(jacob_mat); CHKERRQ(ierr); 957 ierr = PetscFree(prolong_restr_ctx); CHKERRQ(ierr); 958 ierr = PetscFree(prolong_restr_mat); CHKERRQ(ierr); 959 ierr = PetscFree(app_ctx->level_degrees); CHKERRQ(ierr); 960 ierr = PetscFree(ceed_data); CHKERRQ(ierr); 961 962 // libCEED objects 963 CeedVectorDestroy(&form_jacob_ctx->coo_values); 964 CeedQFunctionContextDestroy(&ctx_phys); 965 CeedQFunctionContextDestroy(&ctx_phys_smoother); 966 CeedDestroy(&ceed); 967 968 // PETSc objects 969 ierr = VecDestroy(&U); CHKERRQ(ierr); 970 ierr = VecDestroy(&R); CHKERRQ(ierr); 971 ierr = VecDestroy(&R_loc); CHKERRQ(ierr); 972 ierr = VecDestroy(&F); CHKERRQ(ierr); 973 ierr = VecDestroy(&F_loc); CHKERRQ(ierr); 974 ierr = VecDestroy(&neumann_bcs); CHKERRQ(ierr); 975 ierr = VecDestroy(&bcs_loc); CHKERRQ(ierr); 976 ierr = MatDestroy(&jacob_mat_coarse); CHKERRQ(ierr); 977 ierr = SNESDestroy(&snes); CHKERRQ(ierr); 978 ierr = DMDestroy(&dm_orig); CHKERRQ(ierr); 979 ierr = DMDestroy(&dm_energy); CHKERRQ(ierr); 980 ierr = DMDestroy(&dm_diagnostic); CHKERRQ(ierr); 981 ierr = PetscFree(level_dms); CHKERRQ(ierr); 982 983 // -- Function list 984 ierr = PetscFunctionListDestroy(&problem_functions->setupPhysics); 985 CHKERRQ(ierr); 986 ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedFineLevel); 987 CHKERRQ(ierr); 988 ierr = PetscFunctionListDestroy(&problem_functions->setupLibceedLevel); 989 CHKERRQ(ierr); 990 991 // Structs 992 ierr = PetscFree(res_ctx); CHKERRQ(ierr); 993 ierr = PetscFree(form_jacob_ctx); CHKERRQ(ierr); 994 ierr = PetscFree(jacob_coarse_ctx); CHKERRQ(ierr); 995 ierr = PetscFree(app_ctx); CHKERRQ(ierr); 996 ierr = PetscFree(problem_functions); CHKERRQ(ierr); 997 ierr = PetscFree(units); CHKERRQ(ierr); 998 999 return PetscFinalize(); 1000 } 1001