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 <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 52 53 const char *const PCDeflationTypes[] = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0}; 54 55 /* 56 Private context (data structure) for the deflation preconditioner. 57 */ 58 typedef struct { 59 PetscBool init; /* do only init step - error correction of direction is omitted */ 60 PetscBool pre; /* start with x0 being the solution in the deflation space */ 61 PetscBool correct; /* add CP (Qr) correction to descent direction */ 62 PetscBool truenorm; 63 PetscBool adaptiveconv; 64 PetscReal adaptiveconst; 65 PetscInt reductionfact; 66 Mat W,Wt,AW,WtAW; /* deflation space, coarse problem mats */ 67 KSP WtAWinv; /* deflation coarse problem */ 68 KSPType ksptype; 69 Vec work; 70 Vec *workcoarse; 71 72 PCDeflationSpaceType spacetype; 73 PetscInt spacesize; 74 PetscInt nestedlvl; 75 PetscInt maxnestedlvl; 76 PetscBool extendsp; 77 } PC_Deflation; 78 79 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 80 { 81 PC_Deflation *def = (PC_Deflation*)pc->data; 82 PetscErrorCode ierr; 83 84 PetscFunctionBegin; 85 if (transpose) { 86 def->Wt = W; 87 def->W = NULL; 88 } else { 89 def->W = W; 90 } 91 ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 95 /* TODO create PCDeflationSetSpaceTranspose? */ 96 /*@ 97 PCDeflationSetSpace - Set deflation space matrix (or its transpose). 98 99 Logically Collective on PC 100 101 Input Parameters: 102 + pc - the preconditioner context 103 . W - deflation matrix 104 - tranpose - indicates that W is an explicit transpose of the deflation matrix 105 106 Level: intermediate 107 108 .seealso: PCDEFLATION 109 @*/ 110 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 111 { 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 116 PetscValidHeaderSpecific(W,MAT_CLASSID,2); 117 PetscValidLogicalCollectiveBool(pc,transpose,3); 118 ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 123 { 124 PC_Deflation *def = (PC_Deflation*)pc->data; 125 126 PetscFunctionBegin; 127 def->nestedlvl = current; 128 def->maxnestedlvl = max; 129 PetscFunctionReturn(0); 130 } 131 132 /*@ 133 PCDeflationSetMaxLvl - Set maximum level of deflation. 134 135 Logically Collective on PC 136 137 Input Parameters: 138 + pc - the preconditioner context 139 . max - maximum deflation level 140 141 Level: intermediate 142 143 .seealso: PCDEFLATION 144 @*/ 145 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 146 { 147 PetscErrorCode ierr; 148 149 PetscFunctionBegin; 150 PetscValidLogicalCollectiveInt(pc,max,2); 151 ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 152 PetscFunctionReturn(0); 153 } 154 155 /* 156 TODO CP corection? 157 x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 158 */ 159 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 160 { 161 PC_Deflation *def = (PC_Deflation*)pc->data; 162 Mat A; 163 Vec r,w1,w2; 164 PetscBool nonzero; 165 PetscErrorCode ierr; 166 167 PetscFunctionBegin; 168 w1 = def->workcoarse[0]; 169 w2 = def->workcoarse[1]; 170 r = def->work; 171 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 172 173 ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 174 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 175 if (nonzero) { 176 ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 177 ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 178 } else { 179 ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 180 } 181 182 ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 183 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 184 ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 185 ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 186 PetscFunctionReturn(0); 187 } 188 189 /* 190 if (def->correct) { 191 z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*r 192 } else { 193 z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*r 194 } 195 */ 196 static PetscErrorCode PCApply_Deflation(PC pc,Vec r, Vec z) 197 { 198 PC_Deflation *def = (PC_Deflation*)pc->data; 199 Mat A; 200 Vec u,w1,w2; 201 PetscErrorCode ierr; 202 203 PetscFunctionBegin; 204 w1 = def->workcoarse[0]; 205 w2 = def->workcoarse[1]; 206 u = def->work; 207 ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 208 209 if (!def->AW) { 210 ierr = MatMult(A,r,u);CHKERRQ(ierr); /* u <- A*r */ 211 if (def->correct) ierr = VecAXPY(u,-1.0,r);CHKERRQ(ierr); /* u <- A*r -r */ 212 ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr); /* w1 <- W'*u */ 213 } else { 214 ierr = MatMultTranspose(def->AW,u,w1);CHKERRQ(ierr); /* u <- A*r */ 215 if (def->correct) { 216 ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*u */ 217 ierr = VecAXPY(w1,-1.0,w2);CHKERRQ(ierr); /* w1 <- w1 - w2 */ 218 } 219 } 220 ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 221 ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 222 ierr = VecWAXPY(z,-1.0,u,r);CHKERRQ(ierr); /* z <- r - u */ 223 PetscFunctionReturn(0); 224 } 225 226 static PetscErrorCode PCDeflationSetType_Deflation(PC pc,PCDeflationType type) 227 { 228 PC_Deflation *def = (PC_Deflation*)pc->data; 229 230 PetscFunctionBegin; 231 def->init = PETSC_FALSE; 232 def->pre = PETSC_FALSE; 233 if (type == PC_DEFLATION_POST) { 234 //pc->ops->postsolve = PCPostSolve_Deflation; 235 pc->ops->presolve = 0; 236 } else { 237 pc->ops->presolve = PCPreSolve_Deflation; 238 pc->ops->postsolve = 0; 239 if (type == PC_DEFLATION_INIT) { 240 def->init = PETSC_TRUE; 241 pc->ops->apply = 0; 242 } else { 243 def->pre = PETSC_TRUE; 244 } 245 } 246 PetscFunctionReturn(0); 247 } 248 249 /*@ 250 PCDeflationSetType - Causes the deflation preconditioner to use only a special 251 initial gues or pre/post solve solution update 252 253 Logically Collective on PC 254 255 Input Parameters: 256 + pc - the preconditioner context 257 - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 258 259 Options Database Key: 260 . -pc_deflation_type <pre,init,post> 261 262 Level: intermediate 263 264 Concepts: Deflation preconditioner 265 266 .seealso: PCDeflationGetType() 267 @*/ 268 PetscErrorCode PCDeflationSetType(PC pc,PCDeflationType type) 269 { 270 PetscErrorCode ierr; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 274 ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode PCDeflationGetType_Deflation(PC pc,PCDeflationType *type) 279 { 280 PC_Deflation *def = (PC_Deflation*)pc->data; 281 282 PetscFunctionBegin; 283 if (def->init) { 284 *type = PC_DEFLATION_INIT; 285 } else if (def->pre) { 286 *type = PC_DEFLATION_PRE; 287 } else { 288 *type = PC_DEFLATION_POST; 289 } 290 PetscFunctionReturn(0); 291 } 292 293 /*@ 294 PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner 295 296 Not Collective on PC 297 298 Input Parameter: 299 . pc - the preconditioner context 300 301 Output Parameter: 302 - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 303 304 Level: intermediate 305 306 Concepts: Deflation preconditioner 307 308 .seealso: PCDeflationSetType() 309 @*/ 310 PetscErrorCode PCDeflationGetType(PC pc,PCDeflationType *type) 311 { 312 PetscErrorCode ierr; 313 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 316 ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 static PetscErrorCode PCSetUp_Deflation(PC pc) 321 { 322 PC_Deflation *def = (PC_Deflation*)pc->data; 323 KSP innerksp; 324 PC pcinner; 325 Mat Amat,nextDef=NULL,*mats; 326 PetscInt i,m,red,size,commsize; 327 PetscBool match,flgspd,transp=PETSC_FALSE; 328 MatCompositeType ctype; 329 MPI_Comm comm; 330 const char *prefix; 331 PetscErrorCode ierr; 332 333 PetscFunctionBegin; 334 if (pc->setupcalled) PetscFunctionReturn(0); 335 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 336 ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 337 338 /* compute a deflation space */ 339 if (def->W || def->Wt) { 340 def->spacetype = PC_DEFLATION_SPACE_USER; 341 } else { 342 //ierr = KSPDCGComputeDeflationSpace(ksp);CHKERRQ(ierr); 343 } 344 345 /* nested deflation */ 346 if (def->W) { 347 ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 348 if (match) { 349 ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 350 ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 351 } 352 } else { 353 ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr); 354 ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 355 if (match) { 356 ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 357 ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 358 } 359 transp = PETSC_TRUE; 360 } 361 if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 362 ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 363 if (!transp) { 364 for (i=0; i<size; i++) { 365 ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 366 //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr); 367 } 368 if (def->nestedlvl < def->maxnestedlvl) { 369 size -= 1; 370 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 371 def->W = mats[size]; 372 ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 373 if (size > 1) { 374 ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 375 ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 376 } else { 377 nextDef = mats[0]; 378 ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 379 } 380 } else { 381 /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 382 ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 383 } 384 } 385 ierr = PetscFree(mats);CHKERRQ(ierr); 386 } 387 388 /* setup coarse problem */ 389 if (!def->WtAWinv) { 390 ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */ 391 if (!def->WtAW) { 392 /* TODO add implicit product version ? */ 393 if (!def->AW) { 394 ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 395 } else { 396 ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 397 } 398 /* TODO create MatInheritOption(Mat,MatOption) */ 399 ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 400 ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 401 #if defined(PETSC_USE_DEBUG) 402 /* Check WtAW is not sigular */ 403 PetscReal *norms; 404 ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 405 ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 406 for (i=0; i<m; i++) { 407 if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 408 SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 409 } 410 } 411 ierr = PetscFree(norms);CHKERRQ(ierr); 412 #endif 413 } else { 414 ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 415 } 416 /* TODO use MATINV */ 417 ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 418 ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 419 ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 420 ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 421 ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 422 /* Reduction factor choice */ 423 red = def->reductionfact; 424 if (red < 0) { 425 ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 426 red = ceil((float)commsize/ceil((float)m/commsize)); 427 ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 428 if (match) red = commsize; 429 ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 430 } 431 ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 432 ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 433 ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 434 /* Setup KSP and PC */ 435 if (nextDef) { /* next level for multilevel deflation */ 436 /* set default KSPtype */ 437 if (!def->ksptype) { 438 def->ksptype = KSPFGMRES; 439 if (flgspd) { /* SPD system */ 440 def->ksptype = KSPFCG; 441 } 442 } 443 ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */ 444 ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 445 ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 446 ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 447 /* inherit options TODO if not set */ 448 ((PC_Deflation*)(pcinner))->ksptype = def->ksptype; 449 ((PC_Deflation*)(pcinner))->correct = def->correct; 450 ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv; 451 ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst; 452 ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 453 } else { /* the last level */ 454 ierr = KSPSetType(innerksp,KSPPREONLY);CHKERRQ(ierr); 455 /* TODO Cholesky if flgspd? */ 456 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 457 //TODO remove explicit matSolverPackage 458 if (commsize == red) { 459 ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);CHKERRQ(ierr); 460 } else { 461 ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 462 } 463 } 464 /* TODO use def_[lvl]_ if lvl > 0? */ 465 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 466 if (prefix) { 467 ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 468 ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 469 } else { 470 ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 471 } 472 /* TODO: check if WtAWinv is KSP and move following from this if */ 473 ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 474 //if (def->adaptiveconv) { 475 // PetscReal *rnorm; 476 // PetscNew(&rnorm); 477 // ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr); 478 //} 479 ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 480 } 481 482 /* create work vecs */ 483 ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 484 ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 /* -------------------------------------------------------------------------- */ 488 static PetscErrorCode PCReset_Deflation(PC pc) 489 { 490 PC_Deflation *jac = (PC_Deflation*)pc->data; 491 PetscErrorCode ierr; 492 493 PetscFunctionBegin; 494 PetscFunctionReturn(0); 495 } 496 497 /* 498 PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 499 that was created with PCCreate_Deflation(). 500 501 Input Parameter: 502 . pc - the preconditioner context 503 504 Application Interface Routine: PCDestroy() 505 */ 506 static PetscErrorCode PCDestroy_Deflation(PC pc) 507 { 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 512 513 /* 514 Free the private data structure that was hanging off the PC 515 */ 516 ierr = PetscFree(pc->data);CHKERRQ(ierr); 517 PetscFunctionReturn(0); 518 } 519 520 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 521 { 522 PC_Deflation *jac = (PC_Deflation*)pc->data; 523 PetscErrorCode ierr; 524 PetscBool flg; 525 PCDeflationType deflt,type; 526 527 PetscFunctionBegin; 528 ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr); 529 ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 530 ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 531 if (flg) { 532 ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr); 533 } 534 ierr = PetscOptionsTail();CHKERRQ(ierr); 535 PetscFunctionReturn(0); 536 } 537 538 /*MC 539 PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 540 or to a predefined value 541 542 Options Database Key: 543 + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 544 - -pc_jacobi_abs - use the absolute value of the diagonal entry 545 546 Level: beginner 547 548 Notes: 549 todo 550 551 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 552 PCDeflationSetType(), PCDeflationSetSpace() 553 M*/ 554 555 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 556 { 557 PC_Deflation *def; 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 562 pc->data = (void*)def; 563 564 def->init = PETSC_FALSE; 565 def->pre = PETSC_TRUE; 566 def->correct = PETSC_FALSE; 567 def->truenorm = PETSC_TRUE; 568 def->reductionfact = -1; 569 def->spacetype = PC_DEFLATION_SPACE_HAAR; 570 def->spacesize = 1; 571 def->extendsp = PETSC_FALSE; 572 def->nestedlvl = 0; 573 def->maxnestedlvl = 0; 574 def->adaptiveconv = PETSC_FALSE; 575 def->adaptiveconst = 1.0; 576 577 /* 578 Set the pointers for the functions that are provided above. 579 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 580 are called, they will automatically call these functions. Note we 581 choose not to provide a couple of these functions since they are 582 not needed. 583 */ 584 pc->ops->apply = PCApply_Deflation; 585 pc->ops->applytranspose = PCApply_Deflation; 586 pc->ops->presolve = PCPreSolve_Deflation; 587 pc->ops->postsolve = 0; 588 pc->ops->setup = PCSetUp_Deflation; 589 pc->ops->reset = PCReset_Deflation; 590 pc->ops->destroy = PCDestroy_Deflation; 591 pc->ops->setfromoptions = PCSetFromOptions_Deflation; 592 pc->ops->view = 0; 593 pc->ops->applyrichardson = 0; 594 pc->ops->applysymmetricleft = 0; 595 pc->ops->applysymmetricright = 0; 596 597 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr); 598 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr); 599 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 600 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 601 PetscFunctionReturn(0); 602 } 603 604