1 #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/ /* includes for fortran wrappers */ 2 3 const char *const PCDeflationSpaceTypes[] = {"haar", "db2", "db4", "db8", "db16", "biorth22", "meyer", "aggregation", "user", "PCDeflationSpaceType", "PC_DEFLATION_SPACE_", NULL}; 4 5 static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc, PetscBool flg) { 6 PC_Deflation *def = (PC_Deflation *)pc->data; 7 8 PetscFunctionBegin; 9 def->init = flg; 10 PetscFunctionReturn(0); 11 } 12 13 /*@ 14 PCDeflationSetInitOnly - Do only initialization step. 15 Sets initial guess to the solution on the deflation space but does not apply 16 the deflation preconditioner. The additional preconditioner is still applied. 17 18 Logically Collective 19 20 Input Parameters: 21 + pc - the preconditioner context 22 - flg - default PETSC_FALSE 23 24 Options Database Keys: 25 . -pc_deflation_init_only <false> - if true computes only the special guess 26 27 Level: intermediate 28 29 .seealso: `PCDEFLATION` 30 @*/ 31 PetscErrorCode PCDeflationSetInitOnly(PC pc, PetscBool flg) { 32 PetscFunctionBegin; 33 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 34 PetscValidLogicalCollectiveBool(pc, flg, 2); 35 PetscTryMethod(pc, "PCDeflationSetInitOnly_C", (PC, PetscBool), (pc, flg)); 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc, PetscInt current, PetscInt max) { 40 PC_Deflation *def = (PC_Deflation *)pc->data; 41 42 PetscFunctionBegin; 43 if (current) def->lvl = current; 44 def->maxlvl = max; 45 PetscFunctionReturn(0); 46 } 47 48 /*@ 49 PCDeflationSetLevels - Set the maximum level of deflation nesting. 50 51 Logically Collective 52 53 Input Parameters: 54 + pc - the preconditioner context 55 - max - maximum deflation level 56 57 Options Database Keys: 58 . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 59 60 Level: intermediate 61 62 .seealso: `PCDeflationSetSpaceToCompute()`, `PCDeflationSetSpace()`, `PCDEFLATION` 63 @*/ 64 PetscErrorCode PCDeflationSetLevels(PC pc, PetscInt max) { 65 PetscFunctionBegin; 66 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 67 PetscValidLogicalCollectiveInt(pc, max, 2); 68 PetscTryMethod(pc, "PCDeflationSetLevels_C", (PC, PetscInt, PetscInt), (pc, 0, max)); 69 PetscFunctionReturn(0); 70 } 71 72 static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc, PetscInt red) { 73 PC_Deflation *def = (PC_Deflation *)pc->data; 74 75 PetscFunctionBegin; 76 def->reductionfact = red; 77 PetscFunctionReturn(0); 78 } 79 80 /*@ 81 PCDeflationSetReductionFactor - Set reduction factor for the bottom PCTELESCOPE coarse problem solver. 82 83 Logically Collective 84 85 Input Parameters: 86 + pc - the preconditioner context 87 - red - reduction factor (or PETSC_DETERMINE) 88 89 Options Database Keys: 90 . -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for PCTELESCOPE 91 92 Notes: 93 Default is computed based on the size of the coarse problem. 94 95 Level: intermediate 96 97 .seealso: `PCTELESCOPE`, `PCDEFLATION` 98 @*/ 99 PetscErrorCode PCDeflationSetReductionFactor(PC pc, PetscInt red) { 100 PetscFunctionBegin; 101 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 102 PetscValidLogicalCollectiveInt(pc, red, 2); 103 PetscTryMethod(pc, "PCDeflationSetReductionFactor_C", (PC, PetscInt), (pc, red)); 104 PetscFunctionReturn(0); 105 } 106 107 static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc, PetscScalar fact) { 108 PC_Deflation *def = (PC_Deflation *)pc->data; 109 110 PetscFunctionBegin; 111 /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 112 def->correct = PETSC_TRUE; 113 def->correctfact = fact; 114 if (def->correct == 0.0) { def->correct = PETSC_FALSE; } 115 PetscFunctionReturn(0); 116 } 117 118 /*@ 119 PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 120 The Preconditioner becomes P*M^{-1} + fact*Q. 121 122 Logically Collective 123 124 Input Parameters: 125 + pc - the preconditioner context 126 - fact - correction factor 127 128 Options Database Keys: 129 + -pc_deflation_correction <false> - if true apply coarse problem correction 130 - -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 131 132 Notes: 133 Any non-zero fact enables the coarse problem correction. 134 135 Level: intermediate 136 137 .seealso: `PCDEFLATION` 138 @*/ 139 PetscErrorCode PCDeflationSetCorrectionFactor(PC pc, PetscScalar fact) { 140 PetscFunctionBegin; 141 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 142 PetscValidLogicalCollectiveScalar(pc, fact, 2); 143 PetscTryMethod(pc, "PCDeflationSetCorrectionFactor_C", (PC, PetscScalar), (pc, fact)); 144 PetscFunctionReturn(0); 145 } 146 147 static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc, PCDeflationSpaceType type, PetscInt size) { 148 PC_Deflation *def = (PC_Deflation *)pc->data; 149 150 PetscFunctionBegin; 151 if (type) def->spacetype = type; 152 if (size > 0) def->spacesize = size; 153 PetscFunctionReturn(0); 154 } 155 156 /*@ 157 PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 158 159 Logically Collective 160 161 Input Parameters: 162 + pc - the preconditioner context 163 . type - deflation space type to compute (or PETSC_IGNORE) 164 - size - size of the space to compute (or PETSC_DEFAULT) 165 166 Options Database Keys: 167 + -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 168 - -pc_deflation_compute_space_size <1> - size of the deflation space 169 170 Notes: 171 For wavelet-based deflation, size represents number of levels. 172 173 The deflation space is computed in PCSetUp(). 174 175 Level: intermediate 176 177 .seealso: `PCDeflationSetLevels()`, `PCDEFLATION` 178 @*/ 179 PetscErrorCode PCDeflationSetSpaceToCompute(PC pc, PCDeflationSpaceType type, PetscInt size) { 180 PetscFunctionBegin; 181 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 182 if (type) PetscValidLogicalCollectiveEnum(pc, type, 2); 183 if (size > 0) PetscValidLogicalCollectiveInt(pc, size, 3); 184 PetscTryMethod(pc, "PCDeflationSetSpaceToCompute_C", (PC, PCDeflationSpaceType, PetscInt), (pc, type, size)); 185 PetscFunctionReturn(0); 186 } 187 188 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc, Mat W, PetscBool transpose) { 189 PC_Deflation *def = (PC_Deflation *)pc->data; 190 191 PetscFunctionBegin; 192 /* possibly allows W' = Wt (which is valid but not tested) */ 193 PetscCall(PetscObjectReference((PetscObject)W)); 194 if (transpose) { 195 PetscCall(MatDestroy(&def->Wt)); 196 def->Wt = W; 197 } else { 198 PetscCall(MatDestroy(&def->W)); 199 def->W = W; 200 } 201 PetscFunctionReturn(0); 202 } 203 204 /*@ 205 PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose). 206 207 Logically Collective 208 209 Input Parameters: 210 + pc - the preconditioner context 211 . W - deflation matrix 212 - transpose - indicates that W is an explicit transpose of the deflation matrix 213 214 Notes: 215 Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel 216 deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and 217 the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with 218 W1 as the deflation matrix. This repeats until the maximum level set by 219 PCDeflationSetLevels() is reached or there are no more matrices available. 220 If there are matrices left after reaching the maximum level, 221 they are merged into a deflation matrix ...*W{n-1}*Wn. 222 223 Level: intermediate 224 225 .seealso: `PCDeflationSetLevels()`, `PCDEFLATION` 226 @*/ 227 PetscErrorCode PCDeflationSetSpace(PC pc, Mat W, PetscBool transpose) { 228 PetscFunctionBegin; 229 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 230 PetscValidHeaderSpecific(W, MAT_CLASSID, 2); 231 PetscValidLogicalCollectiveBool(pc, transpose, 3); 232 PetscTryMethod(pc, "PCDeflationSetSpace_C", (PC, Mat, PetscBool), (pc, W, transpose)); 233 PetscFunctionReturn(0); 234 } 235 236 static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc, Mat mat) { 237 PC_Deflation *def = (PC_Deflation *)pc->data; 238 239 PetscFunctionBegin; 240 PetscCall(PetscObjectReference((PetscObject)mat)); 241 PetscCall(MatDestroy(&def->WtA)); 242 def->WtA = mat; 243 PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)def->WtA)); 244 PetscFunctionReturn(0); 245 } 246 247 /*@ 248 PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A). 249 250 Collective 251 252 Input Parameters: 253 + pc - preconditioner context 254 - mat - projection null space matrix 255 256 Level: developer 257 258 .seealso: `PCDEFLATION` 259 @*/ 260 PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc, Mat mat) { 261 PetscFunctionBegin; 262 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 263 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 264 PetscTryMethod(pc, "PCDeflationSetProjectionNullSpaceMat_C", (PC, Mat), (pc, mat)); 265 PetscFunctionReturn(0); 266 } 267 268 static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc, Mat mat) { 269 PC_Deflation *def = (PC_Deflation *)pc->data; 270 271 PetscFunctionBegin; 272 PetscCall(PetscObjectReference((PetscObject)mat)); 273 PetscCall(MatDestroy(&def->WtAW)); 274 def->WtAW = mat; 275 PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)def->WtAW)); 276 PetscFunctionReturn(0); 277 } 278 279 /*@ 280 PCDeflationSetCoarseMat - Set the coarse problem Mat. 281 282 Collective 283 284 Input Parameters: 285 + pc - preconditioner context 286 - mat - coarse problem mat 287 288 Level: developer 289 290 .seealso: `PCDEFLATION` 291 @*/ 292 PetscErrorCode PCDeflationSetCoarseMat(PC pc, Mat mat) { 293 PetscFunctionBegin; 294 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 295 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 296 PetscTryMethod(pc, "PCDeflationSetCoarseMat_C", (PC, Mat), (pc, mat)); 297 PetscFunctionReturn(0); 298 } 299 300 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc, KSP *ksp) { 301 PC_Deflation *def = (PC_Deflation *)pc->data; 302 303 PetscFunctionBegin; 304 *ksp = def->WtAWinv; 305 PetscFunctionReturn(0); 306 } 307 308 /*@ 309 PCDeflationGetCoarseKSP - Returns the coarse problem KSP. 310 311 Not Collective 312 313 Input Parameters: 314 . pc - preconditioner context 315 316 Output Parameters: 317 . ksp - coarse problem KSP context 318 319 Level: advanced 320 321 .seealso: `PCDEFLATION` 322 @*/ 323 PetscErrorCode PCDeflationGetCoarseKSP(PC pc, KSP *ksp) { 324 PetscFunctionBegin; 325 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 326 PetscValidPointer(ksp, 2); 327 PetscTryMethod(pc, "PCDeflationGetCoarseKSP_C", (PC, KSP *), (pc, ksp)); 328 PetscFunctionReturn(0); 329 } 330 331 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc, PC *apc) { 332 PC_Deflation *def = (PC_Deflation *)pc->data; 333 334 PetscFunctionBegin; 335 *apc = def->pc; 336 PetscFunctionReturn(0); 337 } 338 339 /*@ 340 PCDeflationGetPC - Returns the additional preconditioner M^{-1}. 341 342 Not Collective 343 344 Input Parameters: 345 . pc - the preconditioner context 346 347 Output Parameters: 348 . apc - additional preconditioner 349 350 Level: advanced 351 352 .seealso: `PCDEFLATION` 353 @*/ 354 PetscErrorCode PCDeflationGetPC(PC pc, PC *apc) { 355 PetscFunctionBegin; 356 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 357 PetscValidPointer(pc, 1); 358 PetscTryMethod(pc, "PCDeflationGetPC_C", (PC, PC *), (pc, apc)); 359 PetscFunctionReturn(0); 360 } 361 362 /* 363 x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 364 */ 365 static PetscErrorCode PCPreSolve_Deflation(PC pc, KSP ksp, Vec b, Vec x) { 366 PC_Deflation *def = (PC_Deflation *)pc->data; 367 Mat A; 368 Vec r, w1, w2; 369 PetscBool nonzero; 370 371 PetscFunctionBegin; 372 w1 = def->workcoarse[0]; 373 w2 = def->workcoarse[1]; 374 r = def->work; 375 PetscCall(PCGetOperators(pc, NULL, &A)); 376 377 PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero)); 378 PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE)); 379 if (nonzero) { 380 PetscCall(MatMult(A, x, r)); /* r <- b - Ax */ 381 PetscCall(VecAYPX(r, -1.0, b)); 382 } else { 383 PetscCall(VecCopy(b, r)); /* r <- b (x is 0) */ 384 } 385 386 if (def->Wt) { 387 PetscCall(MatMult(def->Wt, r, w1)); /* w1 <- W'*r */ 388 } else { 389 PetscCall(MatMultHermitianTranspose(def->W, r, w1)); /* w1 <- W'*r */ 390 } 391 PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /* w2 <- (W'*A*W)^{-1}*w1 */ 392 PetscCall(MatMult(def->W, w2, r)); /* r <- W*w2 */ 393 PetscCall(VecAYPX(x, 1.0, r)); 394 PetscFunctionReturn(0); 395 } 396 397 /* 398 if (def->correct) { 399 z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 400 } else { 401 z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 402 } 403 */ 404 static PetscErrorCode PCApply_Deflation(PC pc, Vec r, Vec z) { 405 PC_Deflation *def = (PC_Deflation *)pc->data; 406 Mat A; 407 Vec u, w1, w2; 408 409 PetscFunctionBegin; 410 w1 = def->workcoarse[0]; 411 w2 = def->workcoarse[1]; 412 u = def->work; 413 PetscCall(PCGetOperators(pc, NULL, &A)); 414 415 PetscCall(PCApply(def->pc, r, z)); /* z <- M^{-1}*r */ 416 if (!def->init) { 417 PetscCall(MatMult(def->WtA, z, w1)); /* w1 <- W'*A*z */ 418 if (def->correct) { 419 if (def->Wt) { 420 PetscCall(MatMult(def->Wt, r, w2)); /* w2 <- W'*r */ 421 } else { 422 PetscCall(MatMultHermitianTranspose(def->W, r, w2)); /* w2 <- W'*r */ 423 } 424 PetscCall(VecAXPY(w1, -1.0 * def->correctfact, w2)); /* w1 <- w1 - l*w2 */ 425 } 426 PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /* w2 <- (W'*A*W)^{-1}*w1 */ 427 PetscCall(MatMult(def->W, w2, u)); /* u <- W*w2 */ 428 PetscCall(VecAXPY(z, -1.0, u)); /* z <- z - u */ 429 } 430 PetscFunctionReturn(0); 431 } 432 433 static PetscErrorCode PCSetUp_Deflation(PC pc) { 434 PC_Deflation *def = (PC_Deflation *)pc->data; 435 KSP innerksp; 436 PC pcinner; 437 Mat Amat, nextDef = NULL, *mats; 438 PetscInt i, m, red, size; 439 PetscMPIInt commsize; 440 PetscBool match, flgspd, isset, transp = PETSC_FALSE; 441 MatCompositeType ctype; 442 MPI_Comm comm; 443 char prefix[128] = ""; 444 445 PetscFunctionBegin; 446 if (pc->setupcalled) PetscFunctionReturn(0); 447 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 448 PetscCall(PCGetOperators(pc, NULL, &Amat)); 449 if (!def->lvl && !def->prefix) { PetscCall(PCGetOptionsPrefix(pc, &def->prefix)); } 450 if (def->lvl) { PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%d_", (int)def->lvl)); } 451 452 /* compute a deflation space */ 453 if (def->W || def->Wt) { 454 def->spacetype = PC_DEFLATION_SPACE_USER; 455 } else { 456 PetscCall(PCDeflationComputeSpace(pc)); 457 } 458 459 /* nested deflation */ 460 if (def->W) { 461 PetscCall(PetscObjectTypeCompare((PetscObject)def->W, MATCOMPOSITE, &match)); 462 if (match) { 463 PetscCall(MatCompositeGetType(def->W, &ctype)); 464 PetscCall(MatCompositeGetNumberMat(def->W, &size)); 465 } 466 } else { 467 PetscCall(MatCreateHermitianTranspose(def->Wt, &def->W)); 468 PetscCall(PetscObjectTypeCompare((PetscObject)def->Wt, MATCOMPOSITE, &match)); 469 if (match) { 470 PetscCall(MatCompositeGetType(def->Wt, &ctype)); 471 PetscCall(MatCompositeGetNumberMat(def->Wt, &size)); 472 } 473 transp = PETSC_TRUE; 474 } 475 if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 476 if (!transp) { 477 if (def->lvl < def->maxlvl) { 478 PetscCall(PetscMalloc1(size, &mats)); 479 for (i = 0; i < size; i++) { PetscCall(MatCompositeGetMat(def->W, i, &mats[i])); } 480 size -= 1; 481 PetscCall(MatDestroy(&def->W)); 482 def->W = mats[size]; 483 PetscCall(PetscObjectReference((PetscObject)mats[size])); 484 if (size > 1) { 485 PetscCall(MatCreateComposite(comm, size, mats, &nextDef)); 486 PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE)); 487 } else { 488 nextDef = mats[0]; 489 PetscCall(PetscObjectReference((PetscObject)mats[0])); 490 } 491 PetscCall(PetscFree(mats)); 492 } else { 493 /* TODO test merge side performance */ 494 /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */ 495 PetscCall(MatCompositeMerge(def->W)); 496 } 497 } else { 498 if (def->lvl < def->maxlvl) { 499 PetscCall(PetscMalloc1(size, &mats)); 500 for (i = 0; i < size; i++) { PetscCall(MatCompositeGetMat(def->Wt, i, &mats[i])); } 501 size -= 1; 502 PetscCall(MatDestroy(&def->Wt)); 503 def->Wt = mats[0]; 504 PetscCall(PetscObjectReference((PetscObject)mats[0])); 505 if (size > 1) { 506 PetscCall(MatCreateComposite(comm, size, &mats[1], &nextDef)); 507 PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE)); 508 } else { 509 nextDef = mats[1]; 510 PetscCall(PetscObjectReference((PetscObject)mats[1])); 511 } 512 PetscCall(PetscFree(mats)); 513 } else { 514 /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */ 515 PetscCall(MatCompositeMerge(def->Wt)); 516 } 517 } 518 } 519 520 if (transp) { 521 PetscCall(MatDestroy(&def->W)); 522 PetscCall(MatHermitianTranspose(def->Wt, MAT_INITIAL_MATRIX, &def->W)); 523 } 524 525 /* assemble WtA */ 526 if (!def->WtA) { 527 if (def->Wt) { 528 PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA)); 529 } else { 530 #if defined(PETSC_USE_COMPLEX) 531 PetscCall(MatHermitianTranspose(def->W, MAT_INITIAL_MATRIX, &def->Wt)); 532 PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA)); 533 #else 534 PetscCall(MatTransposeMatMult(def->W, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA)); 535 #endif 536 } 537 } 538 /* setup coarse problem */ 539 if (!def->WtAWinv) { 540 PetscCall(MatGetSize(def->W, NULL, &m)); 541 if (!def->WtAW) { 542 PetscCall(MatMatMult(def->WtA, def->W, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtAW)); 543 /* TODO create MatInheritOption(Mat,MatOption) */ 544 PetscCall(MatIsSPDKnown(Amat, &isset, &flgspd)); 545 if (isset) PetscCall(MatSetOption(def->WtAW, MAT_SPD, flgspd)); 546 if (PetscDefined(USE_DEBUG)) { 547 /* Check columns of W are not in kernel of A */ 548 PetscReal *norms; 549 550 PetscCall(PetscMalloc1(m, &norms)); 551 PetscCall(MatGetColumnNorms(def->WtAW, NORM_INFINITY, norms)); 552 for (i = 0; i < m; i++) { PetscCheck(norms[i] > 100 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)def->WtAW), PETSC_ERR_SUP, "Column %" PetscInt_FMT " of W is in kernel of A.", i); } 553 PetscCall(PetscFree(norms)); 554 } 555 } else PetscCall(MatIsSPDKnown(def->WtAW, &isset, &flgspd)); 556 557 /* TODO use MATINV ? */ 558 PetscCall(KSPCreate(comm, &def->WtAWinv)); 559 PetscCall(KSPSetOperators(def->WtAWinv, def->WtAW, def->WtAW)); 560 PetscCall(KSPGetPC(def->WtAWinv, &pcinner)); 561 /* Setup KSP and PC */ 562 if (nextDef) { /* next level for multilevel deflation */ 563 innerksp = def->WtAWinv; 564 /* set default KSPtype */ 565 if (!def->ksptype) { 566 def->ksptype = KSPFGMRES; 567 if (isset && flgspd) { /* SPD system */ 568 def->ksptype = KSPFCG; 569 } 570 } 571 PetscCall(KSPSetType(innerksp, def->ksptype)); /* TODO iherit from KSP + tolerances */ 572 PetscCall(PCSetType(pcinner, PCDEFLATION)); /* TODO create coarse preconditinoner M_c = WtMW ? */ 573 PetscCall(PCDeflationSetSpace(pcinner, nextDef, transp)); 574 PetscCall(PCDeflationSetLevels_Deflation(pcinner, def->lvl + 1, def->maxlvl)); 575 /* inherit options */ 576 if (def->prefix) ((PC_Deflation *)(pcinner->data))->prefix = def->prefix; 577 ((PC_Deflation *)(pcinner->data))->init = def->init; 578 ((PC_Deflation *)(pcinner->data))->ksptype = def->ksptype; 579 ((PC_Deflation *)(pcinner->data))->correct = def->correct; 580 ((PC_Deflation *)(pcinner->data))->correctfact = def->correctfact; 581 ((PC_Deflation *)(pcinner->data))->reductionfact = def->reductionfact; 582 PetscCall(MatDestroy(&nextDef)); 583 } else { /* the last level */ 584 PetscCall(KSPSetType(def->WtAWinv, KSPPREONLY)); 585 PetscCall(PCSetType(pcinner, PCTELESCOPE)); 586 /* do not overwrite PCTELESCOPE */ 587 if (def->prefix) PetscCall(KSPSetOptionsPrefix(def->WtAWinv, def->prefix)); 588 PetscCall(KSPAppendOptionsPrefix(def->WtAWinv, "deflation_tel_")); 589 PetscCall(PCSetFromOptions(pcinner)); 590 PetscCall(PetscObjectTypeCompare((PetscObject)pcinner, PCTELESCOPE, &match)); 591 PetscCheck(match, comm, PETSC_ERR_SUP, "User can not owerwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead."); 592 /* Reduction factor choice */ 593 red = def->reductionfact; 594 if (red < 0) { 595 PetscCallMPI(MPI_Comm_size(comm, &commsize)); 596 red = PetscCeilInt(commsize, PetscCeilInt(m, commsize)); 597 PetscCall(PetscObjectTypeCompareAny((PetscObject)(def->WtAW), &match, MATSEQDENSE, MATMPIDENSE, MATDENSE, "")); 598 if (match) red = commsize; 599 PetscCall(PetscInfo(pc, "Auto choosing reduction factor %" PetscInt_FMT "\n", red)); 600 } 601 PetscCall(PCTelescopeSetReductionFactor(pcinner, red)); 602 PetscCall(PCSetUp(pcinner)); 603 PetscCall(PCTelescopeGetKSP(pcinner, &innerksp)); 604 if (innerksp) { 605 PetscCall(KSPGetPC(innerksp, &pcinner)); 606 PetscCall(PCSetType(pcinner, PCLU)); 607 #if defined(PETSC_HAVE_SUPERLU) 608 PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU, MAT_FACTOR_LU, &match)); 609 if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU)); 610 #endif 611 #if defined(PETSC_HAVE_SUPERLU_DIST) 612 PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &match)); 613 if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU_DIST)); 614 #endif 615 } 616 } 617 618 if (innerksp) { 619 if (def->prefix) { 620 PetscCall(KSPSetOptionsPrefix(innerksp, def->prefix)); 621 PetscCall(KSPAppendOptionsPrefix(innerksp, "deflation_")); 622 } else { 623 PetscCall(KSPSetOptionsPrefix(innerksp, "deflation_")); 624 } 625 PetscCall(KSPAppendOptionsPrefix(innerksp, prefix)); 626 PetscCall(KSPSetFromOptions(innerksp)); 627 PetscCall(KSPSetUp(innerksp)); 628 } 629 } 630 PetscCall(KSPSetFromOptions(def->WtAWinv)); 631 PetscCall(KSPSetUp(def->WtAWinv)); 632 633 /* create preconditioner */ 634 if (!def->pc) { 635 PetscCall(PCCreate(comm, &def->pc)); 636 PetscCall(PCSetOperators(def->pc, Amat, Amat)); 637 PetscCall(PCSetType(def->pc, PCNONE)); 638 if (def->prefix) PetscCall(PCSetOptionsPrefix(def->pc, def->prefix)); 639 PetscCall(PCAppendOptionsPrefix(def->pc, "deflation_")); 640 PetscCall(PCAppendOptionsPrefix(def->pc, prefix)); 641 PetscCall(PCAppendOptionsPrefix(def->pc, "pc_")); 642 PetscCall(PCSetFromOptions(def->pc)); 643 PetscCall(PCSetUp(def->pc)); 644 } 645 646 /* create work vecs */ 647 PetscCall(MatCreateVecs(Amat, NULL, &def->work)); 648 PetscCall(KSPCreateVecs(def->WtAWinv, 2, &def->workcoarse, 0, NULL)); 649 PetscFunctionReturn(0); 650 } 651 652 static PetscErrorCode PCReset_Deflation(PC pc) { 653 PC_Deflation *def = (PC_Deflation *)pc->data; 654 655 PetscFunctionBegin; 656 PetscCall(VecDestroy(&def->work)); 657 PetscCall(VecDestroyVecs(2, &def->workcoarse)); 658 PetscCall(MatDestroy(&def->W)); 659 PetscCall(MatDestroy(&def->Wt)); 660 PetscCall(MatDestroy(&def->WtA)); 661 PetscCall(MatDestroy(&def->WtAW)); 662 PetscCall(KSPDestroy(&def->WtAWinv)); 663 PetscCall(PCDestroy(&def->pc)); 664 PetscFunctionReturn(0); 665 } 666 667 static PetscErrorCode PCDestroy_Deflation(PC pc) { 668 PetscFunctionBegin; 669 PetscCall(PCReset_Deflation(pc)); 670 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", NULL)); 671 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", NULL)); 672 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", NULL)); 673 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", NULL)); 674 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", NULL)); 675 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", NULL)); 676 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", NULL)); 677 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", NULL)); 678 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", NULL)); 679 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", NULL)); 680 PetscCall(PetscFree(pc->data)); 681 PetscFunctionReturn(0); 682 } 683 684 static PetscErrorCode PCView_Deflation(PC pc, PetscViewer viewer) { 685 PC_Deflation *def = (PC_Deflation *)pc->data; 686 PetscInt its; 687 PetscBool iascii; 688 689 PetscFunctionBegin; 690 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 691 if (iascii) { 692 if (def->correct) { PetscCall(PetscViewerASCIIPrintf(viewer, "using CP correction, factor = %g+%gi\n", (double)PetscRealPart(def->correctfact), (double)PetscImaginaryPart(def->correctfact))); } 693 if (!def->lvl) { PetscCall(PetscViewerASCIIPrintf(viewer, "deflation space type: %s\n", PCDeflationSpaceTypes[def->spacetype])); } 694 695 PetscCall(PetscViewerASCIIPrintf(viewer, "--- Additional PC:\n")); 696 PetscCall(PetscViewerASCIIPushTab(viewer)); 697 PetscCall(PCView(def->pc, viewer)); 698 PetscCall(PetscViewerASCIIPopTab(viewer)); 699 700 PetscCall(PetscViewerASCIIPrintf(viewer, "--- Coarse problem solver:\n")); 701 PetscCall(PetscViewerASCIIPushTab(viewer)); 702 PetscCall(KSPGetTotalIterations(def->WtAWinv, &its)); 703 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations: %" PetscInt_FMT "\n", its)); 704 PetscCall(KSPView(def->WtAWinv, viewer)); 705 PetscCall(PetscViewerASCIIPopTab(viewer)); 706 } 707 PetscFunctionReturn(0); 708 } 709 710 static PetscErrorCode PCSetFromOptions_Deflation(PC pc, PetscOptionItems *PetscOptionsObject) { 711 PC_Deflation *def = (PC_Deflation *)pc->data; 712 713 PetscFunctionBegin; 714 PetscOptionsHeadBegin(PetscOptionsObject, "Deflation options"); 715 PetscCall(PetscOptionsBool("-pc_deflation_init_only", "Use only initialization step - Initdef", "PCDeflationSetInitOnly", def->init, &def->init, NULL)); 716 PetscCall(PetscOptionsInt("-pc_deflation_levels", "Maximum of deflation levels", "PCDeflationSetLevels", def->maxlvl, &def->maxlvl, NULL)); 717 PetscCall(PetscOptionsInt("-pc_deflation_reduction_factor", "Reduction factor for coarse problem solution using PCTELESCOPE", "PCDeflationSetReductionFactor", def->reductionfact, &def->reductionfact, NULL)); 718 PetscCall(PetscOptionsBool("-pc_deflation_correction", "Add coarse problem correction Q to P", "PCDeflationSetCorrectionFactor", def->correct, &def->correct, NULL)); 719 PetscCall(PetscOptionsScalar("-pc_deflation_correction_factor", "Set multiple of Q to use as coarse problem correction", "PCDeflationSetCorrectionFactor", def->correctfact, &def->correctfact, NULL)); 720 PetscCall(PetscOptionsEnum("-pc_deflation_compute_space", "Compute deflation space", "PCDeflationSetSpace", PCDeflationSpaceTypes, (PetscEnum)def->spacetype, (PetscEnum *)&def->spacetype, NULL)); 721 PetscCall(PetscOptionsInt("-pc_deflation_compute_space_size", "Set size of the deflation space to compute", "PCDeflationSetSpace", def->spacesize, &def->spacesize, NULL)); 722 PetscCall(PetscOptionsBool("-pc_deflation_space_extend", "Extend deflation space instead of truncating (wavelets)", "PCDeflation", def->extendsp, &def->extendsp, NULL)); 723 PetscOptionsHeadEnd(); 724 PetscFunctionReturn(0); 725 } 726 727 /*MC 728 PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value. 729 730 Options Database Keys: 731 + -pc_deflation_init_only <false> - if true computes only the special guess 732 . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 733 . -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem) 734 . -pc_deflation_correction <false> - if true apply coarse problem correction 735 . -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 736 . -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 737 - -pc_deflation_compute_space_size <1> - size of the deflation space (corresponds to number of levels for wavelet-based deflation) 738 739 Notes: 740 Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2]) 741 preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat. 742 743 The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space. 744 If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the 745 preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner 746 application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues 747 to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some 748 eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates 749 of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper. 750 751 The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can 752 be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size. 753 User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix 754 is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the 755 deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned 756 by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum 757 level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged 758 (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by 759 PCDeflationSetLevels() or by -pc_deflation_levels. 760 761 The coarse problem KSP can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1] 762 from the second level onward. You can also use 763 PCDeflationGetCoarseKSP() to control it from code. The bottom level KSP defaults to 764 KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE. 765 For convenience, the reduction factor can be set by PCDeflationSetReductionFactor() 766 or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size. 767 768 The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for 769 coarse problem KSP apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use 770 PCDeflationGetPC() to control the additional preconditioner from code. It defaults to PCNONE. 771 772 The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can 773 be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can 774 significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We 775 recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create 776 an isolated eigenvalue. 777 778 The options are automatically inherited from the previous deflation level. 779 780 The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also 781 recommend limiting the number of iterations for the coarse problems. 782 783 See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients. 784 Section 4 describes some possible choices for the deflation space. 785 786 Developer Notes: 787 Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech 788 Academy of Sciences and VSB - TU Ostrava. 789 790 Developed from PERMON code used in [4] while on a research stay with 791 Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin. 792 793 References: 794 + * - A. Nicolaides. "Deflation of conjugate gradients with applications to boundary value problems", SIAM J. Numer. Anal. 24.2, 1987. 795 . * - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988. 796 . * - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008. 797 - * - J. Kruzik "Implementation of the Deflated Variants of the Conjugate Gradient Method", Master's thesis, VSB-TUO, 2018 - http://dspace5.vsb.cz/bitstream/handle/10084/130303/KRU0097_USP_N2658_2612T078_2018.pdf 798 799 Level: intermediate 800 801 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 802 `PCDeflationSetInitOnly()`, `PCDeflationSetLevels()`, `PCDeflationSetReductionFactor()`, 803 `PCDeflationSetCorrectionFactor()`, `PCDeflationSetSpaceToCompute()`, 804 `PCDeflationSetSpace()`, `PCDeflationSpaceType`, `PCDeflationSetProjectionNullSpaceMat()`, 805 `PCDeflationSetCoarseMat()`, `PCDeflationGetCoarseKSP()`, `PCDeflationGetPC()` 806 M*/ 807 808 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) { 809 PC_Deflation *def; 810 811 PetscFunctionBegin; 812 PetscCall(PetscNewLog(pc, &def)); 813 pc->data = (void *)def; 814 815 def->init = PETSC_FALSE; 816 def->correct = PETSC_FALSE; 817 def->correctfact = 1.0; 818 def->reductionfact = -1; 819 def->spacetype = PC_DEFLATION_SPACE_HAAR; 820 def->spacesize = 1; 821 def->extendsp = PETSC_FALSE; 822 def->lvl = 0; 823 def->maxlvl = 0; 824 def->W = NULL; 825 def->Wt = NULL; 826 827 pc->ops->apply = PCApply_Deflation; 828 pc->ops->presolve = PCPreSolve_Deflation; 829 pc->ops->setup = PCSetUp_Deflation; 830 pc->ops->reset = PCReset_Deflation; 831 pc->ops->destroy = PCDestroy_Deflation; 832 pc->ops->setfromoptions = PCSetFromOptions_Deflation; 833 pc->ops->view = PCView_Deflation; 834 835 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", PCDeflationSetInitOnly_Deflation)); 836 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", PCDeflationSetLevels_Deflation)); 837 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", PCDeflationSetReductionFactor_Deflation)); 838 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", PCDeflationSetCorrectionFactor_Deflation)); 839 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", PCDeflationSetSpaceToCompute_Deflation)); 840 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", PCDeflationSetSpace_Deflation)); 841 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", PCDeflationSetProjectionNullSpaceMat_Deflation)); 842 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", PCDeflationSetCoarseMat_Deflation)); 843 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", PCDeflationGetCoarseKSP_Deflation)); 844 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", PCDeflationGetPC_Deflation)); 845 PetscFunctionReturn(0); 846 } 847