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 = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 442 /* Setup KSP and PC */ 443 if (nextDef) { /* next level for multilevel deflation */ 444 innerksp = def->WtAWinv; 445 /* set default KSPtype */ 446 if (!def->ksptype) { 447 def->ksptype = KSPFGMRES; 448 if (flgspd) { /* SPD system */ 449 def->ksptype = KSPFCG; 450 } 451 } 452 ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */ 453 ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 454 ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 455 ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 456 /* inherit options TODO if not set */ 457 ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 458 ((PC_Deflation*)(pcinner->data))->correct = def->correct; 459 ((PC_Deflation*)(pcinner->data))->adaptiveconv = def->adaptiveconv; 460 ((PC_Deflation*)(pcinner->data))->adaptiveconst = def->adaptiveconst; 461 ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 462 } else { /* the last level */ 463 ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 464 ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 465 /* ugly hack to not have overwritten PCTELESCOPE */ 466 if (prefix) { 467 ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 468 ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 469 } else { 470 ierr = KSPSetOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 471 } 472 ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 473 /* Reduction factor choice */ 474 red = def->reductionfact; 475 if (red < 0) { 476 ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 477 red = ceil((float)commsize/ceil((float)m/commsize)); 478 ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 479 if (match) red = commsize; 480 ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 481 } 482 ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 483 ierr = PCSetUp(pcinner);CHKERRQ(ierr); 484 ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 485 if (innerksp) { 486 ierr = KSPGetPC(innerksp,&pcinner);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 } 497 498 if (innerksp) { 499 /* TODO use def_[lvl]_ if lvl > 0? */ 500 if (prefix) { 501 ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 502 ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 503 } else { 504 ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 505 } 506 ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 507 ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 508 } 509 /* TODO: check if WtAWinv is KSP and move following from this if */ 510 //if (def->adaptiveconv) { 511 // PetscReal *rnorm; 512 // PetscNew(&rnorm); 513 // ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr); 514 //} 515 } 516 ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 517 ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 518 519 /* create work vecs */ 520 ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 521 ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 522 PetscFunctionReturn(0); 523 } 524 525 static PetscErrorCode PCReset_Deflation(PC pc) 526 { 527 PC_Deflation *def = (PC_Deflation*)pc->data; 528 PetscErrorCode ierr; 529 530 PetscFunctionBegin; 531 ierr = VecDestroy(&def->work);CHKERRQ(ierr); 532 ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 533 ierr = MatDestroy(&def->W);CHKERRQ(ierr); 534 ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 535 ierr = MatDestroy(&def->AW);CHKERRQ(ierr); 536 ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 537 ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 538 PetscFunctionReturn(0); 539 } 540 541 /* 542 PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 543 that was created with PCCreate_Deflation(). 544 545 Input Parameter: 546 . pc - the preconditioner context 547 548 Application Interface Routine: PCDestroy() 549 */ 550 static PetscErrorCode PCDestroy_Deflation(PC pc) 551 { 552 PetscErrorCode ierr; 553 554 PetscFunctionBegin; 555 ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 556 ierr = PetscFree(pc->data);CHKERRQ(ierr); 557 PetscFunctionReturn(0); 558 } 559 560 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 561 { 562 PC_Deflation *def = (PC_Deflation*)pc->data; 563 PetscBool iascii; 564 PetscErrorCode ierr; 565 566 PetscFunctionBegin; 567 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 568 if (iascii) { 569 //if (cg->adaptiveconv) { 570 // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 571 //} 572 if (def->correct) { 573 ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 574 } 575 if (!def->nestedlvl) { 576 ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 577 ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 578 } 579 580 ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 581 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 582 ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 583 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 584 } 585 PetscFunctionReturn(0); 586 } 587 588 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 589 { 590 PC_Deflation *def = (PC_Deflation*)pc->data; 591 PetscErrorCode ierr; 592 PetscBool flg; 593 PCDeflationType deflt,type; 594 595 PetscFunctionBegin; 596 ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr); 597 ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 598 ierr = PetscOptionsEnum("-pc_deflation_type","Determine type of deflation","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 599 if (flg) { 600 ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr); 601 } 602 ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 603 ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 604 ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 605 //TODO add set function and fix manpages 606 ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 607 ierr = PetscOptionsBool("-pc_deflation_correct","Add Qr to descent direction","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 608 ierr = PetscOptionsBool("-pc_deflation_adaptive","Adaptive stopping criteria","PCDeflation",def->adaptiveconv,&def->adaptiveconv,NULL);CHKERRQ(ierr); 609 ierr = PetscOptionsReal("-pc_deflation_adaptive_const","Adaptive stopping criteria constant","PCDeflation",def->adaptiveconst,&def->adaptiveconst,NULL);CHKERRQ(ierr); 610 ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 611 ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 612 ierr = PetscOptionsTail();CHKERRQ(ierr); 613 PetscFunctionReturn(0); 614 } 615 616 /*MC 617 PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 618 or to a predefined value 619 620 Options Database Key: 621 + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 622 - -pc_jacobi_abs - use the absolute value of the diagonal entry 623 624 Level: beginner 625 626 Notes: 627 todo 628 629 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 630 PCDeflationSetType(), PCDeflationSetSpace() 631 M*/ 632 633 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 634 { 635 PC_Deflation *def; 636 PetscErrorCode ierr; 637 638 PetscFunctionBegin; 639 ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 640 pc->data = (void*)def; 641 642 def->init = PETSC_FALSE; 643 def->pre = PETSC_TRUE; 644 def->correct = PETSC_FALSE; 645 def->truenorm = PETSC_TRUE; 646 def->reductionfact = -1; 647 def->spacetype = PC_DEFLATION_SPACE_HAAR; 648 def->spacesize = 1; 649 def->extendsp = PETSC_FALSE; 650 def->nestedlvl = 0; 651 def->maxnestedlvl = 0; 652 def->adaptiveconv = PETSC_FALSE; 653 def->adaptiveconst = 1.0; 654 655 /* 656 Set the pointers for the functions that are provided above. 657 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 658 are called, they will automatically call these functions. Note we 659 choose not to provide a couple of these functions since they are 660 not needed. 661 */ 662 pc->ops->apply = PCApply_Deflation; 663 pc->ops->applytranspose = PCApply_Deflation; 664 pc->ops->presolve = PCPreSolve_Deflation; 665 pc->ops->postsolve = 0; 666 pc->ops->setup = PCSetUp_Deflation; 667 pc->ops->reset = PCReset_Deflation; 668 pc->ops->destroy = PCDestroy_Deflation; 669 pc->ops->setfromoptions = PCSetFromOptions_Deflation; 670 pc->ops->view = PCView_Deflation; 671 pc->ops->applyrichardson = 0; 672 pc->ops->applysymmetricleft = 0; 673 pc->ops->applysymmetricright = 0; 674 675 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr); 676 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr); 677 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 678 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 682