1 2 /* -------------------------------------------------------------------- 3 4 This file implements a Deflation preconditioner in PETSc as part of PC. 5 You can use this as a starting point for implementing your own 6 preconditioner that is not provided with PETSc. (You might also consider 7 just using PCSHELL) 8 9 The following basic routines are required for each preconditioner. 10 PCCreate_XXX() - Creates a preconditioner context 11 PCSetFromOptions_XXX() - Sets runtime options 12 PCApply_XXX() - Applies the preconditioner 13 PCDestroy_XXX() - Destroys the preconditioner context 14 where the suffix "_XXX" denotes a particular implementation, in 15 this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation). 16 These routines are actually called via the common user interface 17 routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 18 so the application code interface remains identical for all 19 preconditioners. 20 21 Another key routine is: 22 PCSetUp_XXX() - Prepares for the use of a preconditioner 23 by setting data structures and options. The interface routine PCSetUp() 24 is not usually called directly by the user, but instead is called by 25 PCApply() if necessary. 26 27 Additional basic routines are: 28 PCView_XXX() - Prints details of runtime options that 29 have actually been used. 30 These are called by application codes via the interface routines 31 PCView(). 32 33 The various types of solvers (preconditioners, Krylov subspace methods, 34 nonlinear solvers, timesteppers) are all organized similarly, so the 35 above description applies to these categories also. One exception is 36 that the analogues of PCApply() for these components are KSPSolve(), 37 SNESSolve(), and TSSolve(). 38 39 Additional optional functionality unique to preconditioners is left and 40 right symmetric preconditioner application via PCApplySymmetricLeft() 41 and PCApplySymmetricRight(). The Deflation implementation is 42 PCApplySymmetricLeftOrRight_Deflation(). 43 44 -------------------------------------------------------------------- */ 45 46 /* 47 Include files needed for the Deflation preconditioner: 48 pcimpl.h - private include file intended for use by all preconditioners 49 */ 50 51 #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 52 53 const char *const PCDeflationSpaceTypes[] = { 54 "haar", 55 "jacket-haar", 56 "db2", 57 "db4", 58 "db8", 59 "db16", 60 "biorth22", 61 "meyer", 62 "aggregation", 63 "slepc", 64 "slepc-cheap", 65 "user", 66 "DdefSpaceType", 67 "Ddef_SPACE_", 68 0 69 }; 70 71 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 72 { 73 PC_Deflation *def = (PC_Deflation*)pc->data; 74 75 PetscFunctionBegin; 76 if (current) def->nestedlvl = current; 77 def->maxnestedlvl = max; 78 PetscFunctionReturn(0); 79 } 80 81 /*@ 82 PCDeflationSetMaxLvl - Set maximum level of deflation. 83 84 Logically Collective on PC 85 86 Input Parameters: 87 + pc - the preconditioner context 88 - max - maximum deflation level 89 90 Level: intermediate 91 92 .seealso: PCDEFLATION 93 @*/ 94 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 95 { 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 100 PetscValidLogicalCollectiveInt(pc,max,2); 101 ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red) 106 { 107 PC_Deflation *def = (PC_Deflation*)pc->data; 108 109 PetscFunctionBegin; 110 def->reductionfact = red; 111 PetscFunctionReturn(0); 112 } 113 114 /*@ 115 PCDeflationSetReductionFactor - Set reduction factor for PCTELESCOPE coarse problem solver. 116 117 Logically Collective on PC 118 119 Input Parameters: 120 + pc - the preconditioner context 121 - red - reduction factor (or PETSC_DETERMINE) 122 123 Level: intermediate 124 125 .seealso: PCDEFLATION 126 @*/ 127 PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red) 128 { 129 PetscErrorCode ierr; 130 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 133 PetscValidLogicalCollectiveInt(pc,red,2); 134 ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr); 135 PetscFunctionReturn(0); 136 } 137 138 static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact) 139 { 140 PC_Deflation *def = (PC_Deflation*)pc->data; 141 142 PetscFunctionBegin; 143 /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 144 def->correct = PETSC_TRUE; 145 def->correctfact = fact; 146 if (fact == PETSC_DEFAULT) { 147 def->correctfact = 1.0; 148 } else if (def->correct == 0.0) { 149 def->correct = PETSC_FALSE; 150 } 151 PetscFunctionReturn(0); 152 } 153 154 /*@ 155 PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 156 The Preconditioner becomes P*M^{-1} + factor*Q. 157 158 Logically Collective on PC 159 160 Input Parameters: 161 + pc - the preconditioner context 162 - fact - correction factor (set 0.0 to disable, PETSC_DEFAULT = 1.0) 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 on PC 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 Notes: 200 For wavelet-based deflation, size represents number of levels. 201 The deflation space is computed in PCSetUP(). 202 203 Level: intermediate 204 205 .seealso: PCDEFLATION 206 @*/ 207 PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 208 { 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 213 if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 214 if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 215 ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 220 { 221 PC_Deflation *def = (PC_Deflation*)pc->data; 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 if (transpose) { 226 def->Wt = W; 227 def->W = NULL; 228 } else { 229 def->W = W; 230 } 231 ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 232 PetscFunctionReturn(0); 233 } 234 235 /* TODO create PCDeflationSetSpaceTranspose? */ 236 /*@ 237 PCDeflationSetSpace - Set deflation space matrix (or its transpose). 238 239 Logically Collective on PC 240 241 Input Parameters: 242 + pc - the preconditioner context 243 . W - deflation matrix 244 - tranpose - indicates that W is an explicit transpose of the deflation matrix 245 246 Level: intermediate 247 248 .seealso: PCDEFLATION 249 @*/ 250 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 251 { 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 256 PetscValidHeaderSpecific(W,MAT_CLASSID,2); 257 PetscValidLogicalCollectiveBool(pc,transpose,3); 258 ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 263 { 264 PC_Deflation *def = (PC_Deflation*)pc->data; 265 PetscErrorCode ierr; 266 267 PetscFunctionBegin; 268 ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 269 def->WtA = mat; 270 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 271 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 272 PetscFunctionReturn(0); 273 } 274 275 /*@ 276 PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 277 278 Not 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 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 294 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 295 ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 296 PetscFunctionReturn(0); 297 } 298 299 static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 300 { 301 PC_Deflation *def = (PC_Deflation*)pc->data; 302 PetscErrorCode ierr; 303 304 PetscFunctionBegin; 305 ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 306 def->WtAW = mat; 307 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 308 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 /*@ 313 PCDeflationSetCoarseMat - Set coarse problem Mat. 314 315 Not Collective 316 317 Input Parameters: 318 + pc - preconditioner context 319 - mat - coarse problem mat 320 321 Level: developer 322 323 .seealso: PCDEFLATION 324 @*/ 325 PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 326 { 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 331 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 332 ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 337 { 338 PC_Deflation *def = (PC_Deflation*)pc->data; 339 340 PetscFunctionBegin; 341 *ksp = def->WtAWinv; 342 PetscFunctionReturn(0); 343 } 344 345 /*@ 346 PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 347 348 Not Collective 349 350 Input Parameters: 351 . pc - preconditioner context 352 353 Output Parameter: 354 . ksp - coarse problem KSP context 355 356 Level: developer 357 358 .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 359 @*/ 360 PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 361 { 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 366 PetscValidPointer(ksp,2); 367 ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 368 PetscFunctionReturn(0); 369 } 370 371 static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 372 { 373 PC_Deflation *def = (PC_Deflation*)pc->data; 374 PetscErrorCode ierr; 375 376 PetscFunctionBegin; 377 ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 378 def->WtAWinv = ksp; 379 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 380 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 381 PetscFunctionReturn(0); 382 } 383 384 /*@ 385 PCDeflationSetCoarseKSP - Set coarse problem KSP. 386 387 Collective on PC 388 389 Input Parameters: 390 + pc - preconditioner context 391 - ksp - coarse problem KSP context 392 393 Level: developer 394 395 .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 396 @*/ 397 PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 398 { 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 403 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 404 PetscCheckSameComm(pc,1,ksp,2); 405 ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 406 PetscFunctionReturn(0); 407 } 408 409 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 410 { 411 PC_Deflation *def = (PC_Deflation*)pc->data; 412 413 PetscFunctionBegin; 414 *apc = def->pc; 415 PetscFunctionReturn(0); 416 } 417 418 /*@ 419 PCDeflationGetPC - Returns a pointer to additional preconditioner. 420 421 Not Collective 422 423 Input Parameters: 424 . pc - the preconditioner context 425 426 Output Parameter: 427 . apc - additional preconditioner 428 429 Level: advanced 430 431 .seealso: PCDeflationSetPC(), PCDEFLATION 432 @*/ 433 PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 434 { 435 PetscErrorCode ierr; 436 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 439 PetscValidPointer(pc,2); 440 ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 445 { 446 PC_Deflation *def = (PC_Deflation*)pc->data; 447 PetscErrorCode ierr; 448 449 PetscFunctionBegin; 450 ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 451 def->pc = apc; 452 ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 453 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 454 PetscFunctionReturn(0); 455 } 456 457 /*@ 458 PCDeflationSetPC - Set additional preconditioner. 459 460 Collective on PC 461 462 Input Parameters: 463 + pc - the preconditioner context 464 - apc - additional preconditioner 465 466 Level: developer 467 468 .seealso: PCDeflationGetPC(), PCDEFLATION 469 @*/ 470 PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 471 { 472 PetscErrorCode ierr; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 476 PetscValidHeaderSpecific(apc,PC_CLASSID,2); 477 PetscCheckSameComm(pc,1,apc,2); 478 ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 479 PetscFunctionReturn(0); 480 } 481 482 /* 483 x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 484 */ 485 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 486 { 487 PC_Deflation *def = (PC_Deflation*)pc->data; 488 Mat A; 489 Vec r,w1,w2; 490 PetscBool nonzero; 491 PetscErrorCode ierr; 492 493 PetscFunctionBegin; 494 w1 = def->workcoarse[0]; 495 w2 = def->workcoarse[1]; 496 r = def->work; 497 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 498 499 ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 500 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 501 if (nonzero) { 502 ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 503 ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 504 } else { 505 ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 506 } 507 508 if (def->Wt) { 509 ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 510 } else { 511 ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 512 } 513 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 514 ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 515 ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 516 PetscFunctionReturn(0); 517 } 518 519 /* 520 if (def->correct) { 521 z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 522 } else { 523 z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 524 } 525 */ 526 static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 527 { 528 PC_Deflation *def = (PC_Deflation*)pc->data; 529 Mat A; 530 Vec u,w1,w2; 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 w1 = def->workcoarse[0]; 535 w2 = def->workcoarse[1]; 536 u = def->work; 537 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 538 539 ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 540 if (!def->init) { 541 ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 542 if (def->correct) { 543 if (def->Wt) { 544 ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 545 } else { 546 ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 547 } 548 ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 549 } 550 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 551 ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 552 ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 553 } 554 PetscFunctionReturn(0); 555 } 556 557 static PetscErrorCode PCSetUp_Deflation(PC pc) 558 { 559 PC_Deflation *def = (PC_Deflation*)pc->data; 560 KSP innerksp; 561 PC pcinner; 562 Mat Amat,nextDef=NULL,*mats; 563 PetscInt i,m,red,size,commsize; 564 PetscBool match,flgspd,transp=PETSC_FALSE; 565 MatCompositeType ctype; 566 MPI_Comm comm; 567 const char *prefix; 568 PetscErrorCode ierr; 569 570 PetscFunctionBegin; 571 if (pc->setupcalled) PetscFunctionReturn(0); 572 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 573 ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 574 575 /* compute a deflation space */ 576 if (def->W || def->Wt) { 577 def->spacetype = PC_DEFLATION_SPACE_USER; 578 } else { 579 ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 580 } 581 582 /* nested deflation */ 583 if (def->W) { 584 ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 585 if (match) { 586 ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 587 ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 588 } 589 } else { 590 ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 591 ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 592 if (match) { 593 ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 594 ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 595 } 596 transp = PETSC_TRUE; 597 } 598 if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 599 if (!transp) { 600 if (def->nestedlvl < def->maxnestedlvl) { 601 ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 602 for (i=0; i<size; i++) { 603 ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 604 } 605 size -= 1; 606 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 607 def->W = mats[size]; 608 ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 609 if (size > 1) { 610 ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 611 ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 612 } else { 613 nextDef = mats[0]; 614 ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 615 } 616 ierr = PetscFree(mats);CHKERRQ(ierr); 617 } else { 618 /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 619 ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 620 } 621 } else { 622 if (def->nestedlvl < def->maxnestedlvl) { 623 ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 624 for (i=0; i<size; i++) { 625 ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 626 } 627 size -= 1; 628 ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 629 def->Wt = mats[0]; 630 ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 631 if (size > 1) { 632 ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 633 ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 634 } else { 635 nextDef = mats[1]; 636 ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 637 } 638 ierr = PetscFree(mats);CHKERRQ(ierr); 639 } else { 640 /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 641 ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 642 } 643 } 644 } 645 646 if (transp) { 647 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 648 ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 649 } 650 651 652 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 653 654 /* assemble WtA */ 655 if (!def->WtA) { 656 if (def->Wt) { 657 ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 658 } else { 659 #if defined(PETSC_USE_COMPLEX) 660 ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 661 ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 662 #else 663 ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 664 #endif 665 } 666 } 667 /* setup coarse problem */ 668 if (!def->WtAWinv) { 669 ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 670 if (!def->WtAW) { 671 ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 672 /* TODO create MatInheritOption(Mat,MatOption) */ 673 ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 674 ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 675 #if defined(PETSC_USE_DEBUG) 676 /* Check columns of W are not in kernel of A */ 677 PetscReal *norms; 678 ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 679 ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 680 for (i=0; i<m; i++) { 681 if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 682 SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 683 } 684 } 685 ierr = PetscFree(norms);CHKERRQ(ierr); 686 #endif 687 } else { 688 ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 689 } 690 /* TODO use MATINV */ 691 ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 692 ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 693 ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 694 /* Setup KSP and PC */ 695 if (nextDef) { /* next level for multilevel deflation */ 696 innerksp = def->WtAWinv; 697 /* set default KSPtype */ 698 if (!def->ksptype) { 699 def->ksptype = KSPFGMRES; 700 if (flgspd) { /* SPD system */ 701 def->ksptype = KSPFCG; 702 } 703 } 704 ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 705 ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 706 ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 707 ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 708 /* inherit options */ 709 ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 710 ((PC_Deflation*)(pcinner->data))->correct = def->correct; 711 ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 712 } else { /* the last level */ 713 ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 714 ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 715 /* ugly hack to not have overwritten PCTELESCOPE */ 716 if (prefix) { 717 ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 718 } 719 ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 720 ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 721 /* Reduction factor choice */ 722 red = def->reductionfact; 723 if (red < 0) { 724 ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 725 red = ceil((float)commsize/ceil((float)m/commsize)); 726 ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 727 if (match) red = commsize; 728 ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 729 } 730 ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 731 ierr = PCSetUp(pcinner);CHKERRQ(ierr); 732 ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 733 if (innerksp) { 734 ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 735 /* TODO Cholesky if flgspd? */ 736 ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 737 //TODO remove explicit matSolverPackage 738 if (commsize == red) { 739 ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 740 } else { 741 ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 742 } 743 } 744 } 745 746 if (innerksp) { 747 /* TODO use def_[lvl]_ if lvl > 0? */ 748 if (prefix) { 749 ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 750 } 751 ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 752 ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 753 ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 754 } 755 } 756 ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 757 ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 758 759 if (!def->pc) { 760 ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 761 ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 762 ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 763 if (prefix) { 764 ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 765 } 766 ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 767 ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 768 ierr = PCSetUp(def->pc);CHKERRQ(ierr); 769 } 770 771 /* create work vecs */ 772 ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 773 ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 774 PetscFunctionReturn(0); 775 } 776 777 static PetscErrorCode PCReset_Deflation(PC pc) 778 { 779 PC_Deflation *def = (PC_Deflation*)pc->data; 780 PetscErrorCode ierr; 781 782 PetscFunctionBegin; 783 ierr = VecDestroy(&def->work);CHKERRQ(ierr); 784 ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 785 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 786 ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 787 ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 788 ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 789 ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 790 ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 791 PetscFunctionReturn(0); 792 } 793 794 /* 795 PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 796 that was created with PCCreate_Deflation(). 797 798 Input Parameter: 799 . pc - the preconditioner context 800 801 Application Interface Routine: PCDestroy() 802 */ 803 static PetscErrorCode PCDestroy_Deflation(PC pc) 804 { 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 809 ierr = PetscFree(pc->data);CHKERRQ(ierr); 810 PetscFunctionReturn(0); 811 } 812 813 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 814 { 815 PC_Deflation *def = (PC_Deflation*)pc->data; 816 PetscBool iascii; 817 PetscErrorCode ierr; 818 819 PetscFunctionBegin; 820 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 821 if (iascii) { 822 //if (cg->adaptiveconv) { 823 // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 824 //} 825 if (def->correct) { 826 ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 827 } 828 if (!def->nestedlvl) { 829 ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 830 ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 831 } 832 833 ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 834 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 835 ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 836 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 837 838 ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 839 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 840 ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 841 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 842 } 843 PetscFunctionReturn(0); 844 } 845 846 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 847 { 848 PC_Deflation *def = (PC_Deflation*)pc->data; 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 853 ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 854 ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 855 ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr); 856 ierr = PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr); 857 ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 858 ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 859 ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 860 //TODO add set function and fix manpages 861 ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 862 ierr = PetscOptionsTail();CHKERRQ(ierr); 863 PetscFunctionReturn(0); 864 } 865 866 /*MC 867 PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 868 or to a predefined value 869 870 Options Database Key: 871 + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 872 - -pc_jacobi_abs - use the absolute value of the diagonal entry 873 874 Level: beginner 875 876 Notes: 877 todo 878 879 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 880 PCDeflationSetType(), PCDeflationSetSpace() 881 M*/ 882 883 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 884 { 885 PC_Deflation *def; 886 PetscErrorCode ierr; 887 888 PetscFunctionBegin; 889 ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 890 pc->data = (void*)def; 891 892 def->init = PETSC_FALSE; 893 def->correct = PETSC_FALSE; 894 def->correctfact = 1.0; 895 def->reductionfact = -1; 896 def->spacetype = PC_DEFLATION_SPACE_HAAR; 897 def->spacesize = 1; 898 def->extendsp = PETSC_FALSE; 899 def->nestedlvl = 0; 900 def->maxnestedlvl = 0; 901 902 /* 903 Set the pointers for the functions that are provided above. 904 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 905 are called, they will automatically call these functions. Note we 906 choose not to provide a couple of these functions since they are 907 not needed. 908 */ 909 pc->ops->apply = PCApply_Deflation; 910 pc->ops->applytranspose = PCApply_Deflation; 911 pc->ops->presolve = PCPreSolve_Deflation; 912 pc->ops->postsolve = 0; 913 pc->ops->setup = PCSetUp_Deflation; 914 pc->ops->reset = PCReset_Deflation; 915 pc->ops->destroy = PCDestroy_Deflation; 916 pc->ops->setfromoptions = PCSetFromOptions_Deflation; 917 pc->ops->view = PCView_Deflation; 918 pc->ops->applyrichardson = 0; 919 pc->ops->applysymmetricleft = 0; 920 pc->ops->applysymmetricright = 0; 921 922 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 923 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr); 924 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr); 925 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 926 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 927 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 928 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 929 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 930 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 931 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 932 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 933 PetscFunctionReturn(0); 934 } 935 936