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