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 PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 72 { 73 PC_Deflation *def = (PC_Deflation*)pc->data; 74 75 PetscFunctionBegin; 76 if (type) def->spacetype = type; 77 if (size > 0) def->spacesize = size; 78 PetscFunctionReturn(0); 79 } 80 81 /*@ 82 PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 83 84 Logically Collective on PC 85 86 Input Parameters: 87 + pc - the preconditioner context 88 . type - deflation space type to compute (or PETSC_IGNORE) 89 - size - size of the space to compute (or PETSC_DEFAULT) 90 91 Notes: 92 For wavelet-based deflation, size represents number of levels. 93 The deflation space is computed in PCSetUP(). 94 95 Level: intermediate 96 97 .seealso: PCDEFLATION 98 @*/ 99 PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 100 { 101 PetscErrorCode ierr; 102 103 PetscFunctionBegin; 104 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 105 if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 106 if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 107 ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 111 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 112 { 113 PC_Deflation *def = (PC_Deflation*)pc->data; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 if (transpose) { 118 def->Wt = W; 119 def->W = NULL; 120 } else { 121 def->W = W; 122 } 123 ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 /* TODO create PCDeflationSetSpaceTranspose? */ 128 /*@ 129 PCDeflationSetSpace - Set deflation space matrix (or its transpose). 130 131 Logically Collective on PC 132 133 Input Parameters: 134 + pc - the preconditioner context 135 . W - deflation matrix 136 - tranpose - indicates that W is an explicit transpose of the deflation matrix 137 138 Level: intermediate 139 140 .seealso: PCDEFLATION 141 @*/ 142 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 143 { 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 148 PetscValidHeaderSpecific(W,MAT_CLASSID,2); 149 PetscValidLogicalCollectiveBool(pc,transpose,3); 150 ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 151 PetscFunctionReturn(0); 152 } 153 154 static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 155 { 156 PC_Deflation *def = (PC_Deflation*)pc->data; 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 161 def->WtA = mat; 162 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 163 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 164 PetscFunctionReturn(0); 165 } 166 167 /*@ 168 PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 169 170 Not Collective 171 172 Input Parameters: 173 + pc - preconditioner context 174 - mat - projection null space matrix 175 176 Level: developer 177 178 .seealso: PCDEFLATION 179 @*/ 180 PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 181 { 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 186 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 187 ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 188 PetscFunctionReturn(0); 189 } 190 191 static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 192 { 193 PC_Deflation *def = (PC_Deflation*)pc->data; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 198 def->WtAW = mat; 199 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 200 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 /*@ 205 PCDeflationSetCoarseMat - Set coarse problem Mat. 206 207 Not Collective 208 209 Input Parameters: 210 + pc - preconditioner context 211 - mat - coarse problem mat 212 213 Level: developer 214 215 .seealso: PCDEFLATION 216 @*/ 217 PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 218 { 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 223 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 224 ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 225 PetscFunctionReturn(0); 226 } 227 228 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 229 { 230 PC_Deflation *def = (PC_Deflation*)pc->data; 231 232 PetscFunctionBegin; 233 if (current) def->nestedlvl = current; 234 def->maxnestedlvl = max; 235 PetscFunctionReturn(0); 236 } 237 238 /*@ 239 PCDeflationSetMaxLvl - Set maximum level of deflation. 240 241 Logically Collective on PC 242 243 Input Parameters: 244 + pc - the preconditioner context 245 - max - maximum deflation level 246 247 Level: intermediate 248 249 .seealso: PCDEFLATION 250 @*/ 251 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 252 { 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 257 PetscValidLogicalCollectiveInt(pc,max,2); 258 ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 263 { 264 PC_Deflation *def = (PC_Deflation*)pc->data; 265 266 PetscFunctionBegin; 267 *apc = def->pc; 268 PetscFunctionReturn(0); 269 } 270 271 /*@ 272 PCDeflationGetPC - Returns a pointer to additional preconditioner. 273 274 Not Collective 275 276 Input Parameters: 277 . pc - the preconditioner context 278 279 Output Parameter: 280 . apc - additional preconditioner 281 282 Level: advanced 283 284 .seealso: PCDeflationSetPC(), PCDEFLATION 285 @*/ 286 PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 287 { 288 PetscErrorCode ierr; 289 290 PetscFunctionBegin; 291 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 292 PetscValidPointer(pc,2); 293 ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 294 PetscFunctionReturn(0); 295 } 296 297 static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 298 { 299 PC_Deflation *def = (PC_Deflation*)pc->data; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 304 def->pc = apc; 305 ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 306 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 307 PetscFunctionReturn(0); 308 } 309 310 /*@ 311 PCDeflationSetPC - Set additional preconditioner. 312 313 Collective on PC 314 315 Input Parameters: 316 + pc - the preconditioner context 317 - apc - additional preconditioner 318 319 Level: developer 320 321 .seealso: PCDeflationGetPC(), PCDEFLATION 322 @*/ 323 PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 324 { 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 329 PetscValidHeaderSpecific(apc,PC_CLASSID,2); 330 PetscCheckSameComm(pc,1,apc,2); 331 ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 332 PetscFunctionReturn(0); 333 } 334 335 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 336 { 337 PC_Deflation *def = (PC_Deflation*)pc->data; 338 339 PetscFunctionBegin; 340 *ksp = def->WtAWinv; 341 PetscFunctionReturn(0); 342 } 343 344 /*@ 345 PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 346 347 Not Collective 348 349 Input Parameters: 350 . pc - preconditioner context 351 352 Output Parameter: 353 . ksp - coarse problem KSP context 354 355 Level: developer 356 357 .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 358 @*/ 359 PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 360 { 361 PetscErrorCode ierr; 362 363 PetscFunctionBegin; 364 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 365 PetscValidPointer(ksp,2); 366 ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 367 PetscFunctionReturn(0); 368 } 369 370 static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 371 { 372 PC_Deflation *def = (PC_Deflation*)pc->data; 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 377 def->WtAWinv = ksp; 378 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 379 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 380 PetscFunctionReturn(0); 381 } 382 383 /*@ 384 PCDeflationSetCoarseKSP - Set coarse problem KSP. 385 386 Collective on PC 387 388 Input Parameters: 389 + pc - preconditioner context 390 - ksp - coarse problem KSP context 391 392 Level: developer 393 394 .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 395 @*/ 396 PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 397 { 398 PetscErrorCode ierr; 399 400 PetscFunctionBegin; 401 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 402 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 403 PetscCheckSameComm(pc,1,ksp,2); 404 ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 405 PetscFunctionReturn(0); 406 } 407 408 /* 409 x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 410 */ 411 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 412 { 413 PC_Deflation *def = (PC_Deflation*)pc->data; 414 Mat A; 415 Vec r,w1,w2; 416 PetscBool nonzero; 417 PetscErrorCode ierr; 418 419 PetscFunctionBegin; 420 w1 = def->workcoarse[0]; 421 w2 = def->workcoarse[1]; 422 r = def->work; 423 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 424 425 ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 426 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 427 if (nonzero) { 428 ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 429 ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 430 } else { 431 ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 432 } 433 434 if (def->Wt) { 435 ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 436 } else { 437 ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 438 } 439 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 440 ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 441 ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 442 PetscFunctionReturn(0); 443 } 444 445 /* 446 if (def->correct) { 447 z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 448 } else { 449 z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 450 } 451 */ 452 static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 453 { 454 PC_Deflation *def = (PC_Deflation*)pc->data; 455 Mat A; 456 Vec u,w1,w2; 457 PetscErrorCode ierr; 458 459 PetscFunctionBegin; 460 w1 = def->workcoarse[0]; 461 w2 = def->workcoarse[1]; 462 u = def->work; 463 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 464 465 ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 466 if (!def->init) { 467 ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 468 if (def->correct) { 469 if (def->Wt) { 470 ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 471 } else { 472 ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 473 } 474 ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 475 } 476 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 477 ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 478 ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 479 } 480 PetscFunctionReturn(0); 481 } 482 483 static PetscErrorCode PCSetUp_Deflation(PC pc) 484 { 485 PC_Deflation *def = (PC_Deflation*)pc->data; 486 KSP innerksp; 487 PC pcinner; 488 Mat Amat,nextDef=NULL,*mats; 489 PetscInt i,m,red,size,commsize; 490 PetscBool match,flgspd,transp=PETSC_FALSE; 491 MatCompositeType ctype; 492 MPI_Comm comm; 493 const char *prefix; 494 PetscErrorCode ierr; 495 496 PetscFunctionBegin; 497 if (pc->setupcalled) PetscFunctionReturn(0); 498 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 499 ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 500 501 /* compute a deflation space */ 502 if (def->W || def->Wt) { 503 def->spacetype = PC_DEFLATION_SPACE_USER; 504 } else { 505 ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 506 } 507 508 /* nested deflation */ 509 if (def->W) { 510 ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 511 if (match) { 512 ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 513 ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 514 } 515 } else { 516 ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 517 ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 518 if (match) { 519 ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 520 ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 521 } 522 transp = PETSC_TRUE; 523 } 524 if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 525 if (!transp) { 526 if (def->nestedlvl < def->maxnestedlvl) { 527 ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 528 for (i=0; i<size; i++) { 529 ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 530 } 531 size -= 1; 532 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 533 def->W = mats[size]; 534 ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 535 if (size > 1) { 536 ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 537 ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 538 } else { 539 nextDef = mats[0]; 540 ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 541 } 542 ierr = PetscFree(mats);CHKERRQ(ierr); 543 } else { 544 /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 545 ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 546 } 547 } else { 548 if (def->nestedlvl < def->maxnestedlvl) { 549 ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 550 for (i=0; i<size; i++) { 551 ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 552 } 553 size -= 1; 554 ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 555 def->Wt = mats[0]; 556 ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 557 if (size > 1) { 558 ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 559 ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 560 } else { 561 nextDef = mats[1]; 562 ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 563 } 564 ierr = PetscFree(mats);CHKERRQ(ierr); 565 } else { 566 /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 567 ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 568 } 569 } 570 } 571 572 if (transp) { 573 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 574 ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 575 } 576 577 578 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 579 580 /* assemble WtA */ 581 if (!def->WtA) { 582 if (def->Wt) { 583 ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 584 } else { 585 #if defined(PETSC_USE_COMPLEX) 586 ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 587 ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 588 #else 589 ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 590 #endif 591 } 592 } 593 /* setup coarse problem */ 594 if (!def->WtAWinv) { 595 ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 596 if (!def->WtAW) { 597 ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 598 /* TODO create MatInheritOption(Mat,MatOption) */ 599 ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 600 ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 601 #if defined(PETSC_USE_DEBUG) 602 /* Check columns of W are not in kernel of A */ 603 PetscReal *norms; 604 ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 605 ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 606 for (i=0; i<m; i++) { 607 if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 608 SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 609 } 610 } 611 ierr = PetscFree(norms);CHKERRQ(ierr); 612 #endif 613 } else { 614 ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 615 } 616 /* TODO use MATINV */ 617 ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 618 ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 619 ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 620 /* Setup KSP and PC */ 621 if (nextDef) { /* next level for multilevel deflation */ 622 innerksp = def->WtAWinv; 623 /* set default KSPtype */ 624 if (!def->ksptype) { 625 def->ksptype = KSPFGMRES; 626 if (flgspd) { /* SPD system */ 627 def->ksptype = KSPFCG; 628 } 629 } 630 ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 631 ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 632 ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 633 ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 634 /* inherit options */ 635 ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 636 ((PC_Deflation*)(pcinner->data))->correct = def->correct; 637 ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 638 } else { /* the last level */ 639 ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 640 ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 641 /* ugly hack to not have overwritten PCTELESCOPE */ 642 if (prefix) { 643 ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 644 } 645 ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 646 ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 647 /* Reduction factor choice */ 648 red = def->reductionfact; 649 if (red < 0) { 650 ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 651 red = ceil((float)commsize/ceil((float)m/commsize)); 652 ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 653 if (match) red = commsize; 654 ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 655 } 656 ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 657 ierr = PCSetUp(pcinner);CHKERRQ(ierr); 658 ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 659 if (innerksp) { 660 ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 661 /* TODO Cholesky if flgspd? */ 662 ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 663 //TODO remove explicit matSolverPackage 664 if (commsize == red) { 665 ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 666 } else { 667 ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 668 } 669 } 670 } 671 672 if (innerksp) { 673 /* TODO use def_[lvl]_ if lvl > 0? */ 674 if (prefix) { 675 ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 676 } 677 ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 678 ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 679 ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 680 } 681 } 682 ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 683 ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 684 685 if (!def->pc) { 686 ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 687 ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 688 ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 689 if (prefix) { 690 ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 691 } 692 ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 693 ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 694 ierr = PCSetUp(def->pc);CHKERRQ(ierr); 695 } 696 697 /* create work vecs */ 698 ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 699 ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 700 PetscFunctionReturn(0); 701 } 702 703 static PetscErrorCode PCReset_Deflation(PC pc) 704 { 705 PC_Deflation *def = (PC_Deflation*)pc->data; 706 PetscErrorCode ierr; 707 708 PetscFunctionBegin; 709 ierr = VecDestroy(&def->work);CHKERRQ(ierr); 710 ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 711 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 712 ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 713 ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 714 ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 715 ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 716 ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 717 PetscFunctionReturn(0); 718 } 719 720 /* 721 PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 722 that was created with PCCreate_Deflation(). 723 724 Input Parameter: 725 . pc - the preconditioner context 726 727 Application Interface Routine: PCDestroy() 728 */ 729 static PetscErrorCode PCDestroy_Deflation(PC pc) 730 { 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 735 ierr = PetscFree(pc->data);CHKERRQ(ierr); 736 PetscFunctionReturn(0); 737 } 738 739 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 740 { 741 PC_Deflation *def = (PC_Deflation*)pc->data; 742 PetscBool iascii; 743 PetscErrorCode ierr; 744 745 PetscFunctionBegin; 746 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 747 if (iascii) { 748 //if (cg->adaptiveconv) { 749 // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 750 //} 751 if (def->correct) { 752 ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 753 } 754 if (!def->nestedlvl) { 755 ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 756 ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 757 } 758 759 ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 760 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 761 ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 762 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 763 764 ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 765 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 766 ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 767 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 768 } 769 PetscFunctionReturn(0); 770 } 771 772 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 773 { 774 PC_Deflation *def = (PC_Deflation*)pc->data; 775 PetscErrorCode ierr; 776 777 PetscFunctionBegin; 778 ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 779 ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 780 ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 781 ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 782 //TODO add set function and fix manpages 783 ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 784 ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 785 ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr); 786 ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 787 ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 788 ierr = PetscOptionsTail();CHKERRQ(ierr); 789 PetscFunctionReturn(0); 790 } 791 792 /*MC 793 PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 794 or to a predefined value 795 796 Options Database Key: 797 + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 798 - -pc_jacobi_abs - use the absolute value of the diagonal entry 799 800 Level: beginner 801 802 Notes: 803 todo 804 805 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 806 PCDeflationSetType(), PCDeflationSetSpace() 807 M*/ 808 809 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 810 { 811 PC_Deflation *def; 812 PetscErrorCode ierr; 813 814 PetscFunctionBegin; 815 ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 816 pc->data = (void*)def; 817 818 def->init = PETSC_FALSE; 819 def->correct = PETSC_FALSE; 820 def->correctfact = 1.0; 821 def->reductionfact = -1; 822 def->spacetype = PC_DEFLATION_SPACE_HAAR; 823 def->spacesize = 1; 824 def->extendsp = PETSC_FALSE; 825 def->nestedlvl = 0; 826 def->maxnestedlvl = 0; 827 828 /* 829 Set the pointers for the functions that are provided above. 830 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 831 are called, they will automatically call these functions. Note we 832 choose not to provide a couple of these functions since they are 833 not needed. 834 */ 835 pc->ops->apply = PCApply_Deflation; 836 pc->ops->applytranspose = PCApply_Deflation; 837 pc->ops->presolve = PCPreSolve_Deflation; 838 pc->ops->postsolve = 0; 839 pc->ops->setup = PCSetUp_Deflation; 840 pc->ops->reset = PCReset_Deflation; 841 pc->ops->destroy = PCDestroy_Deflation; 842 pc->ops->setfromoptions = PCSetFromOptions_Deflation; 843 pc->ops->view = PCView_Deflation; 844 pc->ops->applyrichardson = 0; 845 pc->ops->applysymmetricleft = 0; 846 pc->ops->applysymmetricright = 0; 847 848 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 849 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 850 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 851 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 852 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 853 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 854 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 855 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 856 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 857 PetscFunctionReturn(0); 858 } 859 860