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