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