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