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 linElas -forcing mms 31 // ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperSS -forcing none -ceed /cpu/self 32 // ./elasticity -mesh [.exo file] -degree 2 -E 1 -nu 0.3 -bc_zero 999 -bc_clamp 998 -problem hyperFS -forcing none -ceed /gpu/occa 33 // 34 // Sample meshes can be found at https://github.com/jeremylt/ceedSampleMeshes 35 // 36 //TESTARGS -ceed {ceed_resource} -test -degree 2 -nu 0.3 -E 1 -dm_plex_box_faces 3,3,3 37 38 /// @file 39 /// CEED elasticity example using PETSc with DMPlex 40 41 const char help[] = "Solve solid Problems with CEED and PETSc DMPlex\n"; 42 43 #include "elasticity.h" 44 45 int main(int argc, char **argv) { 46 PetscInt ierr; 47 MPI_Comm comm; 48 // Context structs 49 AppCtx appCtx; // Contains problem options 50 Physics phys; // Contains physical constants 51 Units units; // Contains units scaling 52 // PETSc objects 53 PetscLogStage stageDMSetup, stageLibceedSetup, 54 stageSnesSetup, stageSnesSolve; 55 DM dmOrig; // Distributed DM to clone 56 DM *levelDMs; 57 Vec U, *Ug, *Uloc; // U: solution, R: residual, F: forcing 58 Vec R, Rloc, F, Floc; // g: global, loc: local 59 SNES snes, snesCoarse; 60 Mat *jacobMat, jacobMatCoarse, *prolongRestrMat; 61 // PETSc data 62 UserMult resCtx, jacobCoarseCtx, *jacobCtx; 63 FormJacobCtx formJacobCtx; 64 UserMultProlongRestr *prolongRestrCtx; 65 PCMGCycleType pcmgCycleType = PC_MG_CYCLE_V; 66 // libCEED objects 67 Ceed ceed, ceedFine = NULL; 68 CeedData *ceedData; 69 CeedQFunction qfRestrict = NULL, qfProlong = NULL; 70 // Parameters 71 PetscInt ncompu = 3; // 3 DoFs in 3D 72 PetscInt numLevels = 1, fineLevel = 0; 73 PetscInt *Ugsz, *Ulsz, *Ulocsz; // sz: size 74 PetscInt snesIts = 0; 75 // Timing 76 double startTime, elapsedTime, minTime, maxTime; 77 78 ierr = PetscInitialize(&argc, &argv, NULL, help); 79 if (ierr) 80 return ierr; 81 82 // --------------------------------------------------------------------------- 83 // Process command line options 84 // --------------------------------------------------------------------------- 85 comm = PETSC_COMM_WORLD; 86 87 // -- Set mesh file, polynomial degree, problem type 88 ierr = PetscCalloc1(1, &appCtx); CHKERRQ(ierr); 89 ierr = ProcessCommandLineOptions(comm, appCtx); CHKERRQ(ierr); 90 numLevels = appCtx->numLevels; 91 fineLevel = numLevels - 1; 92 93 // -- Set Poison's ratio, Young's Modulus 94 ierr = PetscMalloc1(1, &phys); CHKERRQ(ierr); 95 ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 96 ierr = ProcessPhysics(comm, phys, units); CHKERRQ(ierr); 97 98 // --------------------------------------------------------------------------- 99 // Setup DM 100 // --------------------------------------------------------------------------- 101 // Performance logging 102 ierr = PetscLogStageRegister("DM and Vector Setup Stage", &stageDMSetup); 103 CHKERRQ(ierr); 104 ierr = PetscLogStagePush(stageDMSetup); CHKERRQ(ierr); 105 106 // -- Create distributed DM from mesh file 107 ierr = CreateDistributedDM(comm, appCtx, &dmOrig); CHKERRQ(ierr); 108 109 // -- Setup DM by polynomial degree 110 ierr = PetscMalloc1(numLevels, &levelDMs); CHKERRQ(ierr); 111 for (int level = 0; level < numLevels; level++) { 112 ierr = DMClone(dmOrig, &levelDMs[level]); CHKERRQ(ierr); 113 ierr = SetupDMByDegree(levelDMs[level], appCtx, appCtx->levelDegrees[level], 114 ncompu); CHKERRQ(ierr); 115 } 116 117 // --------------------------------------------------------------------------- 118 // Setup solution and work vectors 119 // --------------------------------------------------------------------------- 120 // Allocate arrays 121 ierr = PetscMalloc1(numLevels, &Ug); CHKERRQ(ierr); 122 ierr = PetscMalloc1(numLevels, &Uloc); CHKERRQ(ierr); 123 ierr = PetscMalloc1(numLevels, &Ugsz); CHKERRQ(ierr); 124 ierr = PetscMalloc1(numLevels, &Ulsz); CHKERRQ(ierr); 125 ierr = PetscMalloc1(numLevels, &Ulocsz); CHKERRQ(ierr); 126 127 // -- Setup solution vectors for each level 128 for (int level = 0; level < numLevels; level++) { 129 // -- Create global unknown vector U 130 ierr = DMCreateGlobalVector(levelDMs[level], &Ug[level]); CHKERRQ(ierr); 131 ierr = VecGetSize(Ug[level], &Ugsz[level]); CHKERRQ(ierr); 132 // Note: Local size for matShell 133 ierr = VecGetLocalSize(Ug[level], &Ulsz[level]); CHKERRQ(ierr); 134 135 // -- Create local unknown vector Uloc 136 ierr = DMCreateLocalVector(levelDMs[level], &Uloc[level]); CHKERRQ(ierr); 137 // Note: local size for libCEED 138 ierr = VecGetSize(Uloc[level], &Ulocsz[level]); CHKERRQ(ierr); 139 } 140 141 // -- Create residual and forcing vectors 142 ierr = VecDuplicate(Ug[fineLevel], &U); CHKERRQ(ierr); 143 ierr = VecDuplicate(Ug[fineLevel], &R); CHKERRQ(ierr); 144 ierr = VecDuplicate(Ug[fineLevel], &F); CHKERRQ(ierr); 145 ierr = VecDuplicate(Uloc[fineLevel], &Rloc); CHKERRQ(ierr); 146 ierr = VecDuplicate(Uloc[fineLevel], &Floc); CHKERRQ(ierr); 147 148 // Performance logging 149 ierr = PetscLogStagePop(); 150 151 // --------------------------------------------------------------------------- 152 // Set up libCEED 153 // --------------------------------------------------------------------------- 154 // Performance logging 155 ierr = PetscLogStageRegister("libCEED Setup Stage", &stageLibceedSetup); 156 CHKERRQ(ierr); 157 ierr = PetscLogStagePush(stageLibceedSetup); CHKERRQ(ierr); 158 159 // Initalize backend 160 CeedInit(appCtx->ceedResource, &ceed); 161 if (appCtx->degree > 4 && appCtx->ceedResourceFine[0]) 162 CeedInit(appCtx->ceedResourceFine, &ceedFine); 163 164 // -- Create libCEED local forcing vector 165 CeedVector forceCeed; 166 CeedScalar *f; 167 if (appCtx->forcingChoice != FORCE_NONE) { 168 ierr = VecGetArray(Floc, &f); CHKERRQ(ierr); 169 CeedVectorCreate(ceed, Ulocsz[fineLevel], &forceCeed); 170 CeedVectorSetArray(forceCeed, CEED_MEM_HOST, CEED_USE_POINTER, f); 171 } 172 173 // -- Restriction and prolongation QFunction 174 if (appCtx->multigridChoice != MULTIGRID_NONE) { 175 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 176 &qfRestrict); 177 CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 178 &qfProlong); 179 } 180 181 // -- Setup libCEED objects 182 ierr = PetscMalloc1(numLevels, &ceedData); CHKERRQ(ierr); 183 // ---- Setup residual evaluator and geometric information 184 ierr = PetscCalloc1(1, &ceedData[fineLevel]); CHKERRQ(ierr); 185 { 186 bool highOrder = (appCtx->levelDegrees[fineLevel] > 4 && ceedFine); 187 Ceed levelCeed = highOrder ? ceedFine : ceed; 188 ierr = SetupLibceedFineLevel(levelDMs[fineLevel], levelCeed, appCtx, 189 phys, ceedData, fineLevel, ncompu, 190 Ugsz[fineLevel], Ulocsz[fineLevel], forceCeed, 191 qfRestrict, qfProlong); 192 CHKERRQ(ierr); 193 } 194 // ---- Setup Jacobian evaluator and prolongation/restriction 195 for (int level = 0; level < numLevels; level++) { 196 if (level != fineLevel) { 197 ierr = PetscCalloc1(1, &ceedData[level]); CHKERRQ(ierr); 198 } 199 200 // Note: use high order ceed, if specified and degree > 3 201 bool highOrder = (appCtx->levelDegrees[level] > 4 && ceedFine); 202 Ceed levelCeed = highOrder ? ceedFine : ceed; 203 ierr = SetupLibceedLevel(levelDMs[level], levelCeed, appCtx, phys, 204 ceedData, level, ncompu, Ugsz[level], 205 Ulocsz[level], forceCeed, qfRestrict, 206 qfProlong); CHKERRQ(ierr); 207 } 208 209 // Performance logging 210 ierr = PetscLogStagePop(); 211 212 // --------------------------------------------------------------------------- 213 // Setup global forcing vector 214 // --------------------------------------------------------------------------- 215 ierr = VecZeroEntries(F); CHKERRQ(ierr); 216 217 if (appCtx->forcingChoice != FORCE_NONE) { 218 ierr = VecRestoreArray(Floc, &f); CHKERRQ(ierr); 219 ierr = DMLocalToGlobal(levelDMs[fineLevel], Floc, ADD_VALUES, F); 220 CHKERRQ(ierr); 221 CeedVectorDestroy(&forceCeed); 222 } 223 224 // --------------------------------------------------------------------------- 225 // Print problem summary 226 // --------------------------------------------------------------------------- 227 if (!appCtx->testMode) { 228 const char *usedresource; 229 CeedGetResource(ceed, &usedresource); 230 231 ierr = PetscPrintf(comm, 232 "\n-- Elastisticy Example - libCEED + PETSc --\n" 233 " libCEED:\n" 234 " libCEED Backend : %s\n", 235 usedresource); CHKERRQ(ierr); 236 237 if (ceedFine) { 238 ierr = PetscPrintf(comm, 239 " libCEED Backend - high order : %s\n", 240 appCtx->ceedResourceFine); CHKERRQ(ierr); 241 } 242 243 ierr = PetscPrintf(comm, 244 " Problem:\n" 245 " Problem Name : %s\n" 246 " Forcing Function : %s\n" 247 " Mesh:\n" 248 " File : %s\n" 249 " Number of 1D Basis Nodes (p) : %d\n" 250 " Number of 1D Quadrature Points (q) : %d\n" 251 " Global nodes : %D\n" 252 " Owned nodes : %D\n" 253 " DoF per node : %D\n" 254 " Multigrid:\n" 255 " Type : %s\n" 256 " Number of Levels : %d\n", 257 problemTypesForDisp[appCtx->problemChoice], 258 forcingTypesForDisp[appCtx->forcingChoice], 259 appCtx->meshFile[0] ? appCtx->meshFile : "Box Mesh", 260 appCtx->degree + 1, appCtx->degree + 1, 261 Ugsz[fineLevel]/ncompu, Ulsz[fineLevel]/ncompu, ncompu, 262 (appCtx->degree == 1 && 263 appCtx->multigridChoice != MULTIGRID_NONE) ? 264 "Algebraic multigrid" : 265 multigridTypesForDisp[appCtx->multigridChoice], 266 (appCtx->degree == 1 || 267 appCtx->multigridChoice == MULTIGRID_NONE) ? 268 0 : numLevels); CHKERRQ(ierr); 269 270 if (appCtx->multigridChoice != MULTIGRID_NONE) { 271 for (int i = 0; i < 2; i++) { 272 CeedInt level = i ? fineLevel : 0; 273 ierr = PetscPrintf(comm, 274 " Level %D (%s):\n" 275 " Number of 1D Basis Nodes (p) : %d\n" 276 " Global Nodes : %D\n" 277 " Owned Nodes : %D\n", 278 level, i ? "fine" : "coarse", 279 appCtx->levelDegrees[level] + 1, 280 Ugsz[level]/ncompu, Ulsz[level]/ncompu); 281 CHKERRQ(ierr); 282 } 283 } 284 } 285 286 // --------------------------------------------------------------------------- 287 // Setup SNES 288 // --------------------------------------------------------------------------- 289 // Performance logging 290 ierr = PetscLogStageRegister("SNES Setup Stage", &stageSnesSetup); 291 CHKERRQ(ierr); 292 ierr = PetscLogStagePush(stageSnesSetup); CHKERRQ(ierr); 293 294 // Create SNES 295 ierr = SNESCreate(comm, &snes); CHKERRQ(ierr); 296 ierr = SNESSetDM(snes, levelDMs[fineLevel]); CHKERRQ(ierr); 297 298 // -- Jacobian evaluators 299 ierr = PetscMalloc1(numLevels, &jacobCtx); CHKERRQ(ierr); 300 ierr = PetscMalloc1(numLevels, &jacobMat); CHKERRQ(ierr); 301 for (int level = 0; level < numLevels; level++) { 302 // -- Jacobian context for level 303 ierr = PetscMalloc1(1, &jacobCtx[level]); CHKERRQ(ierr); 304 ierr = SetupJacobianCtx(comm, appCtx, levelDMs[level], Ug[level], 305 Uloc[level], ceedData[level], ceed, 306 jacobCtx[level]); CHKERRQ(ierr); 307 308 // -- Form Action of Jacobian on delta_u 309 ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level], Ugsz[level], 310 Ugsz[level], jacobCtx[level], &jacobMat[level]); 311 CHKERRQ(ierr); 312 ierr = MatShellSetOperation(jacobMat[level], MATOP_MULT, 313 (void (*)(void))ApplyJacobian_Ceed); 314 CHKERRQ(ierr); 315 ierr = MatShellSetOperation(jacobMat[level], MATOP_GET_DIAGONAL, 316 (void(*)(void))GetDiag_Ceed); 317 318 } 319 // Note: FormJacobian updates the state count of the Jacobian diagonals 320 // and assembles the Jpre matrix, if needed 321 ierr = PetscMalloc1(1, &formJacobCtx); CHKERRQ(ierr); 322 formJacobCtx->jacobCtx = jacobCtx; 323 formJacobCtx->numLevels = numLevels; 324 formJacobCtx->jacobMat = jacobMat; 325 ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMat[fineLevel], 326 FormJacobian, formJacobCtx); CHKERRQ(ierr); 327 328 // -- Residual evaluation function 329 ierr = PetscMalloc1(1, &resCtx); CHKERRQ(ierr); 330 ierr = PetscMemcpy(resCtx, jacobCtx[fineLevel], 331 sizeof(*jacobCtx[fineLevel])); CHKERRQ(ierr); 332 resCtx->op = ceedData[fineLevel]->opApply; 333 ierr = SNESSetFunction(snes, R, FormResidual_Ceed, resCtx); CHKERRQ(ierr); 334 335 // -- Prolongation/Restriction evaluation 336 ierr = PetscMalloc1(numLevels, &prolongRestrCtx); CHKERRQ(ierr); 337 ierr = PetscMalloc1(numLevels, &prolongRestrMat); CHKERRQ(ierr); 338 for (int level = 1; level < numLevels; level++) { 339 // ---- Prolongation/restriction context for level 340 ierr = PetscMalloc1(1, &prolongRestrCtx[level]); CHKERRQ(ierr); 341 ierr = SetupProlongRestrictCtx(comm, levelDMs[level-1], levelDMs[level], 342 Ug[level], Uloc[level-1], Uloc[level], 343 ceedData[level-1], ceedData[level], ceed, 344 prolongRestrCtx[level]); CHKERRQ(ierr); 345 346 // ---- Form Action of Jacobian on delta_u 347 ierr = MatCreateShell(comm, Ulsz[level], Ulsz[level-1], Ugsz[level], 348 Ugsz[level-1], prolongRestrCtx[level], 349 &prolongRestrMat[level]); CHKERRQ(ierr); 350 // Note: In PCMG, restriction is the transpose of prolongation 351 ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT, 352 (void (*)(void))Prolong_Ceed); 353 ierr = MatShellSetOperation(prolongRestrMat[level], MATOP_MULT_TRANSPOSE, 354 (void (*)(void))Restrict_Ceed); 355 CHKERRQ(ierr); 356 } 357 358 // --------------------------------------------------------------------------- 359 // Setup dummy SNES for AMG coarse solve 360 // --------------------------------------------------------------------------- 361 ierr = SNESCreate(comm, &snesCoarse); CHKERRQ(ierr); 362 ierr = SNESSetDM(snesCoarse, levelDMs[0]); CHKERRQ(ierr); 363 ierr = SNESSetSolution(snesCoarse, Ug[0]); CHKERRQ(ierr); 364 365 // -- Jacobian matrix 366 ierr = DMSetMatType(levelDMs[0], MATAIJ); CHKERRQ(ierr); 367 ierr = DMCreateMatrix(levelDMs[0], &jacobMatCoarse); CHKERRQ(ierr); 368 ierr = SNESSetJacobian(snesCoarse, jacobMatCoarse, jacobMatCoarse, NULL, 369 NULL); CHKERRQ(ierr); 370 371 // -- Residual evaluation function 372 ierr = PetscMalloc1(1, &jacobCoarseCtx); CHKERRQ(ierr); 373 ierr = PetscMemcpy(jacobCoarseCtx, jacobCtx[0], sizeof(*jacobCtx[0])); 374 CHKERRQ(ierr); 375 ierr = SNESSetFunction(snesCoarse, Ug[0], ApplyJacobianCoarse_Ceed, 376 jacobCoarseCtx); CHKERRQ(ierr); 377 378 // -- Update formJacobCtx 379 formJacobCtx->Ucoarse = Ug[0]; 380 formJacobCtx->snesCoarse = snesCoarse; 381 formJacobCtx->jacobMatCoarse = jacobMatCoarse; 382 383 // --------------------------------------------------------------------------- 384 // Setup KSP 385 // --------------------------------------------------------------------------- 386 { 387 PC pc; 388 KSP ksp; 389 390 // -- KSP 391 ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 392 ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 393 ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 394 ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 395 PETSC_DEFAULT); CHKERRQ(ierr); 396 ierr = KSPSetOptionsPrefix(ksp, "outer_"); CHKERRQ(ierr); 397 398 // -- Preconditioning 399 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 400 ierr = PCSetDM(pc, levelDMs[fineLevel]); CHKERRQ(ierr); 401 ierr = PCSetOptionsPrefix(pc, "outer_"); CHKERRQ(ierr); 402 403 if (appCtx->multigridChoice == MULTIGRID_NONE) { 404 // ---- No Multigrid 405 ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 406 ierr = PCJacobiSetType(pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 407 } else if (appCtx->degree == 1) { 408 // ---- AMG for degree 1 409 ierr = SNESSetJacobian(snes, jacobMat[fineLevel], jacobMatCoarse, 410 FormJacobian, formJacobCtx); CHKERRQ(ierr); 411 ierr = KSPSetType(ksp, KSPPREONLY); CHKERRQ(ierr); 412 ierr = PCSetType(pc, PCGAMG); CHKERRQ(ierr); 413 } else { 414 // ---- PCMG 415 ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 416 417 // ------ PCMG levels 418 ierr = PCMGSetLevels(pc, numLevels, NULL); CHKERRQ(ierr); 419 for (int level = 0; level < numLevels; level++) { 420 // -------- Smoother 421 KSP kspSmoother, kspEst; 422 PC pcSmoother; 423 424 // ---------- Smoother KSP 425 ierr = PCMGGetSmoother(pc, level, &kspSmoother); CHKERRQ(ierr); 426 ierr = KSPSetDM(kspSmoother, levelDMs[level]); CHKERRQ(ierr); 427 ierr = KSPSetDMActive(kspSmoother, PETSC_FALSE); CHKERRQ(ierr); 428 429 // ---------- Chebyshev options 430 ierr = KSPSetType(kspSmoother, KSPCHEBYSHEV); CHKERRQ(ierr); 431 ierr = KSPChebyshevEstEigSet(kspSmoother, 0, 0.1, 0, 1.1); 432 CHKERRQ(ierr); 433 ierr = KSPChebyshevEstEigGetKSP(kspSmoother, &kspEst); CHKERRQ(ierr); 434 ierr = KSPSetType(kspEst, KSPCG); CHKERRQ(ierr); 435 ierr = KSPChebyshevEstEigSetUseNoisy(kspSmoother, PETSC_TRUE); 436 CHKERRQ(ierr); 437 ierr = KSPSetOperators(kspSmoother, jacobMat[level], jacobMat[level]); 438 CHKERRQ(ierr); 439 440 // ---------- Smoother preconditioner 441 ierr = KSPGetPC(kspSmoother, &pcSmoother); CHKERRQ(ierr); 442 ierr = PCSetType(pcSmoother, PCJACOBI); CHKERRQ(ierr); 443 ierr = PCJacobiSetType(pcSmoother, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 444 445 // -------- Work vector 446 if (level != fineLevel) { 447 ierr = PCMGSetX(pc, level, Ug[level]); CHKERRQ(ierr); 448 } 449 450 // -------- Level prolongation/restriction operator 451 if (level > 0) { 452 ierr = PCMGSetInterpolation(pc, level, prolongRestrMat[level]); 453 CHKERRQ(ierr); 454 ierr = PCMGSetRestriction(pc, level, prolongRestrMat[level]); 455 CHKERRQ(ierr); 456 } 457 } 458 459 // ------ PCMG coarse solve 460 KSP kspCoarse; 461 PC pcCoarse; 462 463 // -------- Coarse KSP 464 ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 465 ierr = KSPSetType(kspCoarse, KSPPREONLY); CHKERRQ(ierr); 466 ierr = KSPSetOperators(kspCoarse, jacobMatCoarse, jacobMatCoarse); 467 CHKERRQ(ierr); 468 ierr = KSPSetOptionsPrefix(kspCoarse, "coarse_"); CHKERRQ(ierr); 469 470 // -------- Coarse preconditioner 471 ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 472 ierr = PCSetType(pcCoarse, PCGAMG); CHKERRQ(ierr); 473 ierr = PCSetOptionsPrefix(pcCoarse, "coarse_"); CHKERRQ(ierr); 474 475 ierr = KSPSetFromOptions(kspCoarse); CHKERRQ(ierr); 476 ierr = PCSetFromOptions(pcCoarse); CHKERRQ(ierr); 477 478 // ------ PCMG options 479 ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 480 ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 481 ierr = PCMGSetCycleType(pc, pcmgCycleType); CHKERRQ(ierr); 482 } 483 ierr = KSPSetFromOptions(ksp); 484 ierr = PCSetFromOptions(pc); 485 } 486 { 487 // Default to critical-point (CP) line search (related to Wolfe's curvature condition) 488 SNESLineSearch linesearch; 489 490 ierr = SNESGetLineSearch(snes, &linesearch); CHKERRQ(ierr); 491 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP); CHKERRQ(ierr); 492 } 493 494 ierr = SNESSetFromOptions(snes); CHKERRQ(ierr); 495 496 // Performance logging 497 ierr = PetscLogStagePop(); 498 499 // --------------------------------------------------------------------------- 500 // Set initial guess 501 // --------------------------------------------------------------------------- 502 ierr = VecSet(U, 0.0); CHKERRQ(ierr); 503 504 // View solution 505 if (appCtx->viewSoln) { 506 ierr = ViewSolution(comm, U, 0, 0.0); CHKERRQ(ierr); 507 } 508 509 // --------------------------------------------------------------------------- 510 // Solve SNES 511 // --------------------------------------------------------------------------- 512 PetscBool snesMonitor = PETSC_FALSE; 513 ierr = PetscOptionsHasName(NULL, NULL, "-snes_monitor", &snesMonitor); 514 CHKERRQ(ierr); 515 516 // Performance logging 517 ierr = PetscLogStageRegister("SNES Solve Stage", &stageSnesSolve); 518 CHKERRQ(ierr); 519 ierr = PetscLogStagePush(stageSnesSolve); CHKERRQ(ierr); 520 521 // Timing 522 ierr = PetscBarrier((PetscObject)snes); CHKERRQ(ierr); 523 startTime = MPI_Wtime(); 524 525 // Solve for each load increment 526 PetscInt increment; 527 for (increment = 1; increment <= appCtx->numIncrements; increment++) { 528 // -- Log increment count 529 if (snesMonitor) { 530 ierr = PetscPrintf(comm, "%d Load Increment\n", increment - 1); 531 CHKERRQ(ierr); 532 } 533 534 // -- Scale the problem 535 PetscScalar loadIncrement = 1.0*increment / appCtx->numIncrements, 536 scalingFactor = loadIncrement / 537 (increment == 1 ? 1 : resCtx->loadIncrement); 538 resCtx->loadIncrement = loadIncrement; 539 if (appCtx->numIncrements > 1 && appCtx->forcingChoice != FORCE_NONE) { 540 ierr = VecScale(F, scalingFactor); CHKERRQ(ierr); 541 } 542 543 // -- Solve 544 ierr = SNESSolve(snes, F, U); CHKERRQ(ierr); 545 546 // -- View solution 547 if (appCtx->viewSoln) { 548 ierr = ViewSolution(comm, U, increment, loadIncrement); CHKERRQ(ierr); 549 } 550 551 // -- Update SNES iteration count 552 PetscInt its; 553 ierr = SNESGetIterationNumber(snes, &its); CHKERRQ(ierr); 554 snesIts += its; 555 556 // -- Check for divergence 557 SNESConvergedReason reason; 558 ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 559 if (reason < 0) 560 break; 561 } 562 563 // Timing 564 elapsedTime = MPI_Wtime() - startTime; 565 566 // Performance logging 567 ierr = PetscLogStagePop(); 568 569 // --------------------------------------------------------------------------- 570 // Output summary 571 // --------------------------------------------------------------------------- 572 if (!appCtx->testMode) { 573 // -- SNES 574 SNESType snesType; 575 SNESConvergedReason reason; 576 PetscReal rnorm; 577 ierr = SNESGetType(snes, &snesType); CHKERRQ(ierr); 578 ierr = SNESGetConvergedReason(snes, &reason); CHKERRQ(ierr); 579 ierr = SNESGetFunctionNorm(snes, &rnorm); CHKERRQ(ierr); 580 ierr = PetscPrintf(comm, 581 " SNES:\n" 582 " SNES Type : %s\n" 583 " SNES Convergence : %s\n" 584 " Number of Load Increments : %d\n" 585 " Completed Load Increments : %d\n" 586 " Total SNES Iterations : %D\n" 587 " Final rnorm : %e\n", 588 snesType, SNESConvergedReasons[reason], 589 appCtx->numIncrements, increment - 1, 590 snesIts, (double)rnorm); CHKERRQ(ierr); 591 592 // -- KSP 593 KSP ksp; 594 KSPType kspType; 595 ierr = SNESGetKSP(snes, &ksp); CHKERRQ(ierr); 596 ierr = KSPGetType(ksp, &kspType); CHKERRQ(ierr); 597 ierr = PetscPrintf(comm, 598 " Linear Solver:\n" 599 " KSP Type : %s\n", 600 kspType); CHKERRQ(ierr); 601 602 // -- PC 603 if (appCtx->multigridChoice != MULTIGRID_NONE && appCtx->degree > 1) { 604 PC pc; 605 PCMGType pcmgType; 606 ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 607 ierr = PCMGGetType(pc, &pcmgType); CHKERRQ(ierr); 608 ierr = PetscPrintf(comm, 609 " P-Multigrid:\n" 610 " PCMG Type : %s\n" 611 " PCMG Cycle Type : %s\n", 612 PCMGTypes[pcmgType], 613 PCMGCycleTypes[pcmgCycleType]); CHKERRQ(ierr); 614 615 // -- Coarse Solve 616 KSP kspCoarse; 617 PC pcCoarse; 618 PCType pcType; 619 620 ierr = PCMGGetCoarseSolve(pc, &kspCoarse); CHKERRQ(ierr); 621 ierr = KSPGetType(kspCoarse, &kspType); CHKERRQ(ierr); 622 ierr = KSPGetPC(kspCoarse, &pcCoarse); CHKERRQ(ierr); 623 ierr = PCGetType(pcCoarse, &pcType); CHKERRQ(ierr); 624 ierr = PetscPrintf(comm, 625 " Coarse Solve:\n" 626 " KSP Type : %s\n" 627 " PC Type : %s\n", 628 kspType, pcType); CHKERRQ(ierr); 629 } 630 } 631 632 // --------------------------------------------------------------------------- 633 // Compute solve time 634 // --------------------------------------------------------------------------- 635 if (!appCtx->testMode) { 636 ierr = MPI_Allreduce(&elapsedTime, &minTime, 1, MPI_DOUBLE, MPI_MIN, comm); 637 CHKERRQ(ierr); 638 ierr = MPI_Allreduce(&elapsedTime, &maxTime, 1, MPI_DOUBLE, MPI_MAX, comm); 639 CHKERRQ(ierr); 640 ierr = PetscPrintf(comm, 641 " Performance:\n" 642 " SNES Solve Time : %g (%g) sec\n", 643 maxTime, minTime); CHKERRQ(ierr); 644 } 645 646 // --------------------------------------------------------------------------- 647 // Compute error 648 // --------------------------------------------------------------------------- 649 if (appCtx->forcingChoice == FORCE_MMS) { 650 CeedScalar l2Error = 1., l2Unorm = 1.; 651 const CeedScalar *truearray; 652 Vec errorVec, trueVec; 653 654 // -- Work vectors 655 ierr = VecDuplicate(U, &errorVec); CHKERRQ(ierr); 656 ierr = VecSet(errorVec, 0.0); CHKERRQ(ierr); 657 ierr = VecDuplicate(U, &trueVec); CHKERRQ(ierr); 658 ierr = VecSet(trueVec, 0.0); CHKERRQ(ierr); 659 660 // -- Assemble global true solution vector 661 CeedVectorGetArrayRead(ceedData[fineLevel]->truesoln, CEED_MEM_HOST, 662 &truearray); 663 ierr = VecPlaceArray(resCtx->Yloc, truearray); CHKERRQ(ierr); 664 ierr = DMLocalToGlobal(resCtx->dm, resCtx->Yloc, INSERT_VALUES, trueVec); 665 CHKERRQ(ierr); 666 ierr = VecResetArray(resCtx->Yloc); CHKERRQ(ierr); 667 CeedVectorRestoreArrayRead(ceedData[fineLevel]->truesoln, &truearray); 668 669 // -- Compute L2 error 670 ierr = VecWAXPY(errorVec, -1.0, U, trueVec); CHKERRQ(ierr); 671 ierr = VecNorm(errorVec, NORM_2, &l2Error); CHKERRQ(ierr); 672 ierr = VecNorm(U, NORM_2, &l2Unorm); CHKERRQ(ierr); 673 l2Error /= l2Unorm; 674 675 // -- Output 676 if (!appCtx->testMode || l2Error > 0.05) { 677 ierr = PetscPrintf(comm, 678 " L2 Error : %e\n", 679 l2Error); CHKERRQ(ierr); 680 } 681 682 // -- Cleanup 683 ierr = VecDestroy(&errorVec); CHKERRQ(ierr); 684 ierr = VecDestroy(&trueVec); CHKERRQ(ierr); 685 } 686 687 // --------------------------------------------------------------------------- 688 // Free objects 689 // --------------------------------------------------------------------------- 690 // Data in arrays per level 691 for (int level = 0; level < numLevels; level++) { 692 // Vectors 693 ierr = VecDestroy(&Ug[level]); CHKERRQ(ierr); 694 ierr = VecDestroy(&Uloc[level]); CHKERRQ(ierr); 695 696 // Jacobian matrix and data 697 ierr = VecDestroy(&jacobCtx[level]->Yloc); CHKERRQ(ierr); 698 ierr = MatDestroy(&jacobMat[level]); CHKERRQ(ierr); 699 ierr = PetscFree(jacobCtx[level]); CHKERRQ(ierr); 700 701 // Prolongation/Restriction matrix and data 702 if (level > 0) { 703 ierr = VecDestroy(&prolongRestrCtx[level]->multVec); CHKERRQ(ierr); 704 ierr = PetscFree(prolongRestrCtx[level]); CHKERRQ(ierr); 705 ierr = MatDestroy(&prolongRestrMat[level]); CHKERRQ(ierr); 706 } 707 708 // DM 709 ierr = DMDestroy(&levelDMs[level]); CHKERRQ(ierr); 710 711 // libCEED objects 712 ierr = CeedDataDestroy(level, ceedData[level]); CHKERRQ(ierr); 713 } 714 715 // Arrays 716 ierr = PetscFree(Ug); CHKERRQ(ierr); 717 ierr = PetscFree(Uloc); CHKERRQ(ierr); 718 ierr = PetscFree(Ugsz); CHKERRQ(ierr); 719 ierr = PetscFree(Ulsz); CHKERRQ(ierr); 720 ierr = PetscFree(Ulocsz); CHKERRQ(ierr); 721 ierr = PetscFree(jacobCtx); CHKERRQ(ierr); 722 ierr = PetscFree(jacobMat); CHKERRQ(ierr); 723 ierr = PetscFree(prolongRestrCtx); CHKERRQ(ierr); 724 ierr = PetscFree(prolongRestrMat); CHKERRQ(ierr); 725 ierr = PetscFree(appCtx->levelDegrees); CHKERRQ(ierr); 726 ierr = PetscFree(ceedData); CHKERRQ(ierr); 727 728 // libCEED objects 729 CeedQFunctionDestroy(&qfRestrict); 730 CeedQFunctionDestroy(&qfProlong); 731 CeedDestroy(&ceed); 732 CeedDestroy(&ceedFine); 733 734 // PETSc objects 735 ierr = VecDestroy(&U); CHKERRQ(ierr); 736 ierr = VecDestroy(&R); CHKERRQ(ierr); 737 ierr = VecDestroy(&Rloc); CHKERRQ(ierr); 738 ierr = VecDestroy(&F); CHKERRQ(ierr); 739 ierr = VecDestroy(&Floc); CHKERRQ(ierr); 740 ierr = MatDestroy(&jacobMatCoarse); CHKERRQ(ierr); 741 ierr = SNESDestroy(&snes); CHKERRQ(ierr); 742 ierr = SNESDestroy(&snesCoarse); CHKERRQ(ierr); 743 ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 744 ierr = PetscFree(levelDMs); CHKERRQ(ierr); 745 746 // Structs 747 ierr = PetscFree(resCtx); CHKERRQ(ierr); 748 ierr = PetscFree(formJacobCtx); CHKERRQ(ierr); 749 ierr = PetscFree(jacobCoarseCtx); CHKERRQ(ierr); 750 ierr = PetscFree(appCtx); CHKERRQ(ierr); 751 ierr = PetscFree(phys); CHKERRQ(ierr); 752 ierr = PetscFree(units); CHKERRQ(ierr); 753 754 return PetscFinalize(); 755 } 756