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