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