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