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