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