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