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