137eeb815SJakub Kruzik 237eeb815SJakub Kruzik /* -------------------------------------------------------------------- 337eeb815SJakub Kruzik 437eeb815SJakub Kruzik This file implements a Deflation preconditioner in PETSc as part of PC. 537eeb815SJakub Kruzik You can use this as a starting point for implementing your own 637eeb815SJakub Kruzik preconditioner that is not provided with PETSc. (You might also consider 737eeb815SJakub Kruzik just using PCSHELL) 837eeb815SJakub Kruzik 937eeb815SJakub Kruzik The following basic routines are required for each preconditioner. 1037eeb815SJakub Kruzik PCCreate_XXX() - Creates a preconditioner context 1137eeb815SJakub Kruzik PCSetFromOptions_XXX() - Sets runtime options 1237eeb815SJakub Kruzik PCApply_XXX() - Applies the preconditioner 1337eeb815SJakub Kruzik PCDestroy_XXX() - Destroys the preconditioner context 1437eeb815SJakub Kruzik where the suffix "_XXX" denotes a particular implementation, in 1537eeb815SJakub Kruzik this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation). 1637eeb815SJakub Kruzik These routines are actually called via the common user interface 1737eeb815SJakub Kruzik routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 1837eeb815SJakub Kruzik so the application code interface remains identical for all 1937eeb815SJakub Kruzik preconditioners. 2037eeb815SJakub Kruzik 2137eeb815SJakub Kruzik Another key routine is: 2237eeb815SJakub Kruzik PCSetUp_XXX() - Prepares for the use of a preconditioner 2337eeb815SJakub Kruzik by setting data structures and options. The interface routine PCSetUp() 2437eeb815SJakub Kruzik is not usually called directly by the user, but instead is called by 2537eeb815SJakub Kruzik PCApply() if necessary. 2637eeb815SJakub Kruzik 2737eeb815SJakub Kruzik Additional basic routines are: 2837eeb815SJakub Kruzik PCView_XXX() - Prints details of runtime options that 2937eeb815SJakub Kruzik have actually been used. 3037eeb815SJakub Kruzik These are called by application codes via the interface routines 3137eeb815SJakub Kruzik PCView(). 3237eeb815SJakub Kruzik 3337eeb815SJakub Kruzik The various types of solvers (preconditioners, Krylov subspace methods, 3437eeb815SJakub Kruzik nonlinear solvers, timesteppers) are all organized similarly, so the 3537eeb815SJakub Kruzik above description applies to these categories also. One exception is 3637eeb815SJakub Kruzik that the analogues of PCApply() for these components are KSPSolve(), 3737eeb815SJakub Kruzik SNESSolve(), and TSSolve(). 3837eeb815SJakub Kruzik 3937eeb815SJakub Kruzik Additional optional functionality unique to preconditioners is left and 4037eeb815SJakub Kruzik right symmetric preconditioner application via PCApplySymmetricLeft() 4137eeb815SJakub Kruzik and PCApplySymmetricRight(). The Deflation implementation is 4237eeb815SJakub Kruzik PCApplySymmetricLeftOrRight_Deflation(). 4337eeb815SJakub Kruzik 4437eeb815SJakub Kruzik -------------------------------------------------------------------- */ 4537eeb815SJakub Kruzik 4637eeb815SJakub Kruzik /* 4737eeb815SJakub Kruzik Include files needed for the Deflation preconditioner: 4837eeb815SJakub Kruzik pcimpl.h - private include file intended for use by all preconditioners 4937eeb815SJakub Kruzik */ 5037eeb815SJakub Kruzik 51292e2e67SJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 5237eeb815SJakub Kruzik 5337eeb815SJakub Kruzik const char *const PCDeflationTypes[] = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0}; 5437eeb815SJakub Kruzik 558513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 568513960bSJakub Kruzik "haar", 578513960bSJakub Kruzik "jacket-haar", 588513960bSJakub Kruzik "db2", 598513960bSJakub Kruzik "db4", 608513960bSJakub Kruzik "db8", 618513960bSJakub Kruzik "db16", 628513960bSJakub Kruzik "biorth22", 638513960bSJakub Kruzik "meyer", 648513960bSJakub Kruzik "aggregation", 658513960bSJakub Kruzik "slepc", 668513960bSJakub Kruzik "slepc-cheap", 678513960bSJakub Kruzik "user", 688513960bSJakub Kruzik "DdefSpaceType", 698513960bSJakub Kruzik "Ddef_SPACE_", 708513960bSJakub Kruzik 0 718513960bSJakub Kruzik }; 728513960bSJakub Kruzik 738513960bSJakub Kruzik 74e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 75e662bc50SJakub Kruzik { 76e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 77e662bc50SJakub Kruzik PetscErrorCode ierr; 78e662bc50SJakub Kruzik 79e662bc50SJakub Kruzik PetscFunctionBegin; 80e662bc50SJakub Kruzik if (transpose) { 81e662bc50SJakub Kruzik def->Wt = W; 82e662bc50SJakub Kruzik def->W = NULL; 83e662bc50SJakub Kruzik } else { 84e662bc50SJakub Kruzik def->W = W; 85e662bc50SJakub Kruzik } 86e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 87e662bc50SJakub Kruzik PetscFunctionReturn(0); 88e662bc50SJakub Kruzik } 89e662bc50SJakub Kruzik 90e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 91e662bc50SJakub Kruzik /*@ 92e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 93e662bc50SJakub Kruzik 94e662bc50SJakub Kruzik Logically Collective on PC 95e662bc50SJakub Kruzik 96e662bc50SJakub Kruzik Input Parameters: 97e662bc50SJakub Kruzik + pc - the preconditioner context 98e662bc50SJakub Kruzik . W - deflation matrix 99e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 100e662bc50SJakub Kruzik 101e662bc50SJakub Kruzik Level: intermediate 102e662bc50SJakub Kruzik 103e662bc50SJakub Kruzik .seealso: PCDEFLATION 104e662bc50SJakub Kruzik @*/ 105e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 106e662bc50SJakub Kruzik { 107e662bc50SJakub Kruzik PetscErrorCode ierr; 108e662bc50SJakub Kruzik 109e662bc50SJakub Kruzik PetscFunctionBegin; 110e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 111e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 112e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 113e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 114e662bc50SJakub Kruzik PetscFunctionReturn(0); 115e662bc50SJakub Kruzik } 116e662bc50SJakub Kruzik 1175c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 1185c0c31c5SJakub Kruzik { 1195c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1205c0c31c5SJakub Kruzik 1215c0c31c5SJakub Kruzik PetscFunctionBegin; 1225c0c31c5SJakub Kruzik def->nestedlvl = current; 1235c0c31c5SJakub Kruzik def->maxnestedlvl = max; 1245c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1255c0c31c5SJakub Kruzik } 1265c0c31c5SJakub Kruzik 1275c0c31c5SJakub Kruzik /*@ 1285c0c31c5SJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 1295c0c31c5SJakub Kruzik 1305c0c31c5SJakub Kruzik Logically Collective on PC 1315c0c31c5SJakub Kruzik 1325c0c31c5SJakub Kruzik Input Parameters: 1335c0c31c5SJakub Kruzik + pc - the preconditioner context 134*22b0793eSJakub Kruzik - max - maximum deflation level 1355c0c31c5SJakub Kruzik 1365c0c31c5SJakub Kruzik Level: intermediate 1375c0c31c5SJakub Kruzik 1385c0c31c5SJakub Kruzik .seealso: PCDEFLATION 1395c0c31c5SJakub Kruzik @*/ 1405c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 1415c0c31c5SJakub Kruzik { 1425c0c31c5SJakub Kruzik PetscErrorCode ierr; 1435c0c31c5SJakub Kruzik 1445c0c31c5SJakub Kruzik PetscFunctionBegin; 145*22b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1465c0c31c5SJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 147*22b0793eSJakub Kruzik /* TODO allow setting only for the top level */ 1485c0c31c5SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 1495c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1505c0c31c5SJakub Kruzik } 151e662bc50SJakub Kruzik 152*22b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 153*22b0793eSJakub Kruzik { 154*22b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 155*22b0793eSJakub Kruzik PetscErrorCode ierr; 156*22b0793eSJakub Kruzik 157*22b0793eSJakub Kruzik PetscFunctionBegin; 158*22b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 159*22b0793eSJakub Kruzik def->pc = apc; 160*22b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 161*22b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 162*22b0793eSJakub Kruzik PetscFunctionReturn(0); 163*22b0793eSJakub Kruzik } 164*22b0793eSJakub Kruzik 165*22b0793eSJakub Kruzik /*@ 166*22b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 167*22b0793eSJakub Kruzik 168*22b0793eSJakub Kruzik Logically Collective on PC 169*22b0793eSJakub Kruzik 170*22b0793eSJakub Kruzik Input Parameters: 171*22b0793eSJakub Kruzik + pc - the preconditioner context 172*22b0793eSJakub Kruzik - apc - additional preconditioner 173*22b0793eSJakub Kruzik 174*22b0793eSJakub Kruzik Level: advanced 175*22b0793eSJakub Kruzik 176*22b0793eSJakub Kruzik .seealso: PCDEFLATION 177*22b0793eSJakub Kruzik @*/ 178*22b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 179*22b0793eSJakub Kruzik { 180*22b0793eSJakub Kruzik PetscErrorCode ierr; 181*22b0793eSJakub Kruzik 182*22b0793eSJakub Kruzik PetscFunctionBegin; 183*22b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 184*22b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 185*22b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 186*22b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 187*22b0793eSJakub Kruzik PetscFunctionReturn(0); 188*22b0793eSJakub Kruzik } 189*22b0793eSJakub Kruzik 19037eeb815SJakub Kruzik /* 191b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 192b27c8092SJakub Kruzik */ 193b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 194b27c8092SJakub Kruzik { 195b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 196b27c8092SJakub Kruzik Mat A; 197b27c8092SJakub Kruzik Vec r,w1,w2; 198b27c8092SJakub Kruzik PetscBool nonzero; 199b27c8092SJakub Kruzik PetscErrorCode ierr; 20037eeb815SJakub Kruzik 201b27c8092SJakub Kruzik PetscFunctionBegin; 202b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 203b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 204b27c8092SJakub Kruzik r = def->work; 205b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 20637eeb815SJakub Kruzik 207b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 208b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 209b27c8092SJakub Kruzik if (nonzero) { 210b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 211b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 212b27c8092SJakub Kruzik } else { 213b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 214b27c8092SJakub Kruzik } 215b27c8092SJakub Kruzik 216b27c8092SJakub Kruzik ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 217b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 218b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 219b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 220b27c8092SJakub Kruzik PetscFunctionReturn(0); 221b27c8092SJakub Kruzik } 22237eeb815SJakub Kruzik 223f8bf229cSJakub Kruzik /* 224f8bf229cSJakub Kruzik if (def->correct) { 225*22b0793eSJakub Kruzik z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*M^{-1}*r 226f8bf229cSJakub Kruzik } else { 227*22b0793eSJakub Kruzik z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*M^{-1}*r 228f8bf229cSJakub Kruzik } 22937eeb815SJakub Kruzik */ 230f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 231f8bf229cSJakub Kruzik { 232f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 233f8bf229cSJakub Kruzik Mat A; 234f8bf229cSJakub Kruzik Vec u,w1,w2; 235f8bf229cSJakub Kruzik PetscErrorCode ierr; 236f8bf229cSJakub Kruzik 237f8bf229cSJakub Kruzik PetscFunctionBegin; 238f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 239f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 240f8bf229cSJakub Kruzik u = def->work; 241f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 242f8bf229cSJakub Kruzik 243*22b0793eSJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); 244*22b0793eSJakub Kruzik if (!def->init) { 245f8bf229cSJakub Kruzik if (!def->AW) { 246*22b0793eSJakub Kruzik ierr = MatMult(A,z,u);CHKERRQ(ierr); /* u <- A*z */ 247*22b0793eSJakub Kruzik /* TODO correct const */ 248*22b0793eSJakub Kruzik if (def->correct) ierr = VecAXPY(u,-1.0,z);CHKERRQ(ierr); /* u <- A*z -z */ 249f8bf229cSJakub Kruzik ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr); /* w1 <- W'*u */ 250f8bf229cSJakub Kruzik } else { 251*22b0793eSJakub Kruzik /* ONLY if A SYM */ 252*22b0793eSJakub Kruzik ierr = MatMultTranspose(def->AW,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*r */ 253f8bf229cSJakub Kruzik if (def->correct) { 254*22b0793eSJakub Kruzik ierr = MatMultTranspose(def->W,z,w2);CHKERRQ(ierr); /* w2 <- W'*u */ 255f8bf229cSJakub Kruzik ierr = VecAXPY(w1,-1.0,w2);CHKERRQ(ierr); /* w1 <- w1 - w2 */ 256f8bf229cSJakub Kruzik } 257f8bf229cSJakub Kruzik } 258f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 259f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 260*22b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 261*22b0793eSJakub Kruzik } 262f8bf229cSJakub Kruzik PetscFunctionReturn(0); 263f8bf229cSJakub Kruzik } 264f8bf229cSJakub Kruzik 265b27c8092SJakub Kruzik static PetscErrorCode PCDeflationSetType_Deflation(PC pc,PCDeflationType type) 266b27c8092SJakub Kruzik { 267b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 268b27c8092SJakub Kruzik 269b27c8092SJakub Kruzik PetscFunctionBegin; 270b27c8092SJakub Kruzik def->init = PETSC_FALSE; 271b27c8092SJakub Kruzik def->pre = PETSC_FALSE; 272b27c8092SJakub Kruzik if (type == PC_DEFLATION_POST) { 273b27c8092SJakub Kruzik //pc->ops->postsolve = PCPostSolve_Deflation; 274b27c8092SJakub Kruzik pc->ops->presolve = 0; 275b27c8092SJakub Kruzik } else { 276b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 277b27c8092SJakub Kruzik pc->ops->postsolve = 0; 278b27c8092SJakub Kruzik if (type == PC_DEFLATION_INIT) { 279b27c8092SJakub Kruzik def->init = PETSC_TRUE; 280b27c8092SJakub Kruzik pc->ops->apply = 0; 281b27c8092SJakub Kruzik } else { 282b27c8092SJakub Kruzik def->pre = PETSC_TRUE; 283b27c8092SJakub Kruzik } 284b27c8092SJakub Kruzik } 285b27c8092SJakub Kruzik PetscFunctionReturn(0); 286b27c8092SJakub Kruzik } 287b27c8092SJakub Kruzik 288b27c8092SJakub Kruzik /*@ 289b27c8092SJakub Kruzik PCDeflationSetType - Causes the deflation preconditioner to use only a special 290b27c8092SJakub Kruzik initial gues or pre/post solve solution update 291b27c8092SJakub Kruzik 292b27c8092SJakub Kruzik Logically Collective on PC 293b27c8092SJakub Kruzik 294b27c8092SJakub Kruzik Input Parameters: 295b27c8092SJakub Kruzik + pc - the preconditioner context 296b27c8092SJakub Kruzik - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 297b27c8092SJakub Kruzik 298b27c8092SJakub Kruzik Options Database Key: 299b27c8092SJakub Kruzik . -pc_deflation_type <pre,init,post> 300b27c8092SJakub Kruzik 301b27c8092SJakub Kruzik Level: intermediate 302b27c8092SJakub Kruzik 303b27c8092SJakub Kruzik Concepts: Deflation preconditioner 304b27c8092SJakub Kruzik 305b27c8092SJakub Kruzik .seealso: PCDeflationGetType() 306b27c8092SJakub Kruzik @*/ 307b27c8092SJakub Kruzik PetscErrorCode PCDeflationSetType(PC pc,PCDeflationType type) 308b27c8092SJakub Kruzik { 309b27c8092SJakub Kruzik PetscErrorCode ierr; 310b27c8092SJakub Kruzik 311b27c8092SJakub Kruzik PetscFunctionBegin; 312b27c8092SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 313b27c8092SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr); 314b27c8092SJakub Kruzik PetscFunctionReturn(0); 315b27c8092SJakub Kruzik } 316b27c8092SJakub Kruzik 317b27c8092SJakub Kruzik static PetscErrorCode PCDeflationGetType_Deflation(PC pc,PCDeflationType *type) 318b27c8092SJakub Kruzik { 319b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 320b27c8092SJakub Kruzik 321b27c8092SJakub Kruzik PetscFunctionBegin; 322b27c8092SJakub Kruzik if (def->init) { 323b27c8092SJakub Kruzik *type = PC_DEFLATION_INIT; 324b27c8092SJakub Kruzik } else if (def->pre) { 325b27c8092SJakub Kruzik *type = PC_DEFLATION_PRE; 326b27c8092SJakub Kruzik } else { 327b27c8092SJakub Kruzik *type = PC_DEFLATION_POST; 328b27c8092SJakub Kruzik } 329b27c8092SJakub Kruzik PetscFunctionReturn(0); 330b27c8092SJakub Kruzik } 331b27c8092SJakub Kruzik 332b27c8092SJakub Kruzik /*@ 333b27c8092SJakub Kruzik PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner 334b27c8092SJakub Kruzik 335b27c8092SJakub Kruzik Not Collective on PC 336b27c8092SJakub Kruzik 337b27c8092SJakub Kruzik Input Parameter: 338b27c8092SJakub Kruzik . pc - the preconditioner context 339b27c8092SJakub Kruzik 340b27c8092SJakub Kruzik Output Parameter: 341b27c8092SJakub Kruzik - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 342b27c8092SJakub Kruzik 343b27c8092SJakub Kruzik Level: intermediate 344b27c8092SJakub Kruzik 345b27c8092SJakub Kruzik Concepts: Deflation preconditioner 346b27c8092SJakub Kruzik 347b27c8092SJakub Kruzik .seealso: PCDeflationSetType() 348b27c8092SJakub Kruzik @*/ 349b27c8092SJakub Kruzik PetscErrorCode PCDeflationGetType(PC pc,PCDeflationType *type) 350b27c8092SJakub Kruzik { 351b27c8092SJakub Kruzik PetscErrorCode ierr; 352b27c8092SJakub Kruzik 353b27c8092SJakub Kruzik PetscFunctionBegin; 354b27c8092SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 355b27c8092SJakub Kruzik ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr); 356b27c8092SJakub Kruzik PetscFunctionReturn(0); 357b27c8092SJakub Kruzik } 358b27c8092SJakub Kruzik 35937eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 36037eeb815SJakub Kruzik { 361409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 362409a357bSJakub Kruzik KSP innerksp; 363409a357bSJakub Kruzik PC pcinner; 364409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 365409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 366409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 367409a357bSJakub Kruzik MatCompositeType ctype; 368409a357bSJakub Kruzik MPI_Comm comm; 369409a357bSJakub Kruzik const char *prefix; 37037eeb815SJakub Kruzik PetscErrorCode ierr; 37137eeb815SJakub Kruzik 37237eeb815SJakub Kruzik PetscFunctionBegin; 373409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 374409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 375f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 376f8bf229cSJakub Kruzik 377f8bf229cSJakub Kruzik /* compute a deflation space */ 378409a357bSJakub Kruzik if (def->W || def->Wt) { 379409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 380409a357bSJakub Kruzik } else { 381e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 38237eeb815SJakub Kruzik } 38337eeb815SJakub Kruzik 384409a357bSJakub Kruzik /* nested deflation */ 385409a357bSJakub Kruzik if (def->W) { 386409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 387409a357bSJakub Kruzik if (match) { 388409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 389409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 39037eeb815SJakub Kruzik } 391409a357bSJakub Kruzik } else { 392409a357bSJakub Kruzik ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr); 393409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 394409a357bSJakub Kruzik if (match) { 395409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 396409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 397409a357bSJakub Kruzik } 398409a357bSJakub Kruzik transp = PETSC_TRUE; 399409a357bSJakub Kruzik } 400409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 401409a357bSJakub Kruzik if (!transp) { 40228da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 40328da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 404409a357bSJakub Kruzik for (i=0; i<size; i++) { 405409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 406409a357bSJakub Kruzik } 407409a357bSJakub Kruzik size -= 1; 408409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 409409a357bSJakub Kruzik def->W = mats[size]; 410409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 411409a357bSJakub Kruzik if (size > 1) { 412409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 413409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 414409a357bSJakub Kruzik } else { 415409a357bSJakub Kruzik nextDef = mats[0]; 416409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 417409a357bSJakub Kruzik } 41828da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 419409a357bSJakub Kruzik } else { 420409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 421409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 422409a357bSJakub Kruzik } 42328da0170SJakub Kruzik } else { 42428da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 42528da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 42628da0170SJakub Kruzik for (i=0; i<size; i++) { 42728da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 42828da0170SJakub Kruzik } 42928da0170SJakub Kruzik size -= 1; 43028da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 43128da0170SJakub Kruzik def->Wt = mats[0]; 43228da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 43328da0170SJakub Kruzik if (size > 1) { 43428da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 43528da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 43628da0170SJakub Kruzik } else { 43728da0170SJakub Kruzik nextDef = mats[1]; 43828da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 439409a357bSJakub Kruzik } 440409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 44128da0170SJakub Kruzik } else { 44228da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 44328da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 44428da0170SJakub Kruzik } 44528da0170SJakub Kruzik } 44628da0170SJakub Kruzik } 44728da0170SJakub Kruzik 44828da0170SJakub Kruzik if (transp) { 44928da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 45028da0170SJakub Kruzik ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 451409a357bSJakub Kruzik } 452409a357bSJakub Kruzik 453*22b0793eSJakub Kruzik 454*22b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 455*22b0793eSJakub Kruzik 456409a357bSJakub Kruzik /* setup coarse problem */ 457409a357bSJakub Kruzik if (!def->WtAWinv) { 45828da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 459409a357bSJakub Kruzik if (!def->WtAW) { 460409a357bSJakub Kruzik /* TODO add implicit product version ? */ 461409a357bSJakub Kruzik if (!def->AW) { 462409a357bSJakub Kruzik ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 463409a357bSJakub Kruzik } else { 464409a357bSJakub Kruzik ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 465409a357bSJakub Kruzik } 466409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 467409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 468409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 469409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 470409a357bSJakub Kruzik /* Check WtAW is not sigular */ 471409a357bSJakub Kruzik PetscReal *norms; 472409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 473409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 474409a357bSJakub Kruzik for (i=0; i<m; i++) { 475409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 476409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 477409a357bSJakub Kruzik } 478409a357bSJakub Kruzik } 479409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 480409a357bSJakub Kruzik #endif 481409a357bSJakub Kruzik } else { 482409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 483409a357bSJakub Kruzik } 484409a357bSJakub Kruzik /* TODO use MATINV */ 485409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 486409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 487409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 488557a2f7dSJakub Kruzik /* Setup KSP and PC */ 489557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 490557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 491557a2f7dSJakub Kruzik /* set default KSPtype */ 492557a2f7dSJakub Kruzik if (!def->ksptype) { 493557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 494557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 495557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 496557a2f7dSJakub Kruzik } 497557a2f7dSJakub Kruzik } 498557a2f7dSJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */ 499557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 500557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 501557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 502557a2f7dSJakub Kruzik /* inherit options TODO if not set */ 503557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 504557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 505557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->adaptiveconv = def->adaptiveconv; 506557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->adaptiveconst = def->adaptiveconst; 507557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 508557a2f7dSJakub Kruzik } else { /* the last level */ 509557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 510409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 511ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 512ac66f006SJakub Kruzik if (prefix) { 513ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 514ac66f006SJakub Kruzik } 515*22b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 516ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 517409a357bSJakub Kruzik /* Reduction factor choice */ 518409a357bSJakub Kruzik red = def->reductionfact; 519409a357bSJakub Kruzik if (red < 0) { 520409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 521409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 522409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 523409a357bSJakub Kruzik if (match) red = commsize; 524409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 525409a357bSJakub Kruzik } 526409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 527ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 528409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 529ac66f006SJakub Kruzik if (innerksp) { 530409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 531409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 532ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 533409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 534409a357bSJakub Kruzik if (commsize == red) { 535ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 536409a357bSJakub Kruzik } else { 537ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 538409a357bSJakub Kruzik } 539409a357bSJakub Kruzik } 540557a2f7dSJakub Kruzik } 541557a2f7dSJakub Kruzik 542557a2f7dSJakub Kruzik if (innerksp) { 543409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 544409a357bSJakub Kruzik if (prefix) { 545409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 546409a357bSJakub Kruzik } 547*22b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 548557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 549557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 550ac66f006SJakub Kruzik } 551409a357bSJakub Kruzik /* TODO: check if WtAWinv is KSP and move following from this if */ 552409a357bSJakub Kruzik //if (def->adaptiveconv) { 553409a357bSJakub Kruzik // PetscReal *rnorm; 554409a357bSJakub Kruzik // PetscNew(&rnorm); 555409a357bSJakub Kruzik // ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr); 556409a357bSJakub Kruzik //} 557409a357bSJakub Kruzik } 558557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 559557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 560409a357bSJakub Kruzik 561*22b0793eSJakub Kruzik if (!def->pc) { 562*22b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 563*22b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 564*22b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 565*22b0793eSJakub Kruzik if (prefix) { 566*22b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 567*22b0793eSJakub Kruzik } 568*22b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 569*22b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 570*22b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 571*22b0793eSJakub Kruzik } 572*22b0793eSJakub Kruzik 573f8bf229cSJakub Kruzik /* create work vecs */ 574f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 575f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 57637eeb815SJakub Kruzik PetscFunctionReturn(0); 57737eeb815SJakub Kruzik } 578b30d4775SJakub Kruzik 57937eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 58037eeb815SJakub Kruzik { 581b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 58237eeb815SJakub Kruzik PetscErrorCode ierr; 58337eeb815SJakub Kruzik 58437eeb815SJakub Kruzik PetscFunctionBegin; 585b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 586b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 587b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 588b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 589b30d4775SJakub Kruzik ierr = MatDestroy(&def->AW);CHKERRQ(ierr); 590b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 591b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 592*22b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 59337eeb815SJakub Kruzik PetscFunctionReturn(0); 59437eeb815SJakub Kruzik } 59537eeb815SJakub Kruzik 59637eeb815SJakub Kruzik /* 59737eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 59837eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 59937eeb815SJakub Kruzik 60037eeb815SJakub Kruzik Input Parameter: 60137eeb815SJakub Kruzik . pc - the preconditioner context 60237eeb815SJakub Kruzik 60337eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 60437eeb815SJakub Kruzik */ 60537eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 60637eeb815SJakub Kruzik { 60737eeb815SJakub Kruzik PetscErrorCode ierr; 60837eeb815SJakub Kruzik 60937eeb815SJakub Kruzik PetscFunctionBegin; 61037eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 61137eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 61237eeb815SJakub Kruzik PetscFunctionReturn(0); 61337eeb815SJakub Kruzik } 61437eeb815SJakub Kruzik 6158513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 6168513960bSJakub Kruzik { 6178513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 6188513960bSJakub Kruzik PetscBool iascii; 6198513960bSJakub Kruzik PetscErrorCode ierr; 6208513960bSJakub Kruzik 6218513960bSJakub Kruzik PetscFunctionBegin; 6228513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6238513960bSJakub Kruzik if (iascii) { 6248513960bSJakub Kruzik //if (cg->adaptiveconv) { 6258513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 6268513960bSJakub Kruzik //} 6278513960bSJakub Kruzik if (def->correct) { 6288513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 6298513960bSJakub Kruzik } 6308513960bSJakub Kruzik if (!def->nestedlvl) { 6318513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 6328513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 6338513960bSJakub Kruzik } 6348513960bSJakub Kruzik 635*22b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 636*22b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 637*22b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 638*22b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 639*22b0793eSJakub Kruzik 6408513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 6418513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6428513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 6438513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 6448513960bSJakub Kruzik } 6458513960bSJakub Kruzik PetscFunctionReturn(0); 6468513960bSJakub Kruzik } 6478513960bSJakub Kruzik 64837eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 64937eeb815SJakub Kruzik { 650880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 65137eeb815SJakub Kruzik PetscErrorCode ierr; 65237eeb815SJakub Kruzik PetscBool flg; 65337eeb815SJakub Kruzik PCDeflationType deflt,type; 65437eeb815SJakub Kruzik 65537eeb815SJakub Kruzik PetscFunctionBegin; 65637eeb815SJakub Kruzik ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr); 65737eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 658880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_type","Determine type of deflation","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 65937eeb815SJakub Kruzik if (flg) { 66037eeb815SJakub Kruzik ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr); 66137eeb815SJakub Kruzik } 662880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 663880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 664880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 665880d05ceSJakub Kruzik //TODO add set function and fix manpages 666880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 667880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add Qr to descent direction","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 668880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_adaptive","Adaptive stopping criteria","PCDeflation",def->adaptiveconv,&def->adaptiveconv,NULL);CHKERRQ(ierr); 669880d05ceSJakub Kruzik ierr = PetscOptionsReal("-pc_deflation_adaptive_const","Adaptive stopping criteria constant","PCDeflation",def->adaptiveconst,&def->adaptiveconst,NULL);CHKERRQ(ierr); 670880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 671880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 67237eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 67337eeb815SJakub Kruzik PetscFunctionReturn(0); 67437eeb815SJakub Kruzik } 67537eeb815SJakub Kruzik 67637eeb815SJakub Kruzik /*MC 677e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 678e361f8b5SJakub Kruzik or to a predefined value 67937eeb815SJakub Kruzik 68037eeb815SJakub Kruzik Options Database Key: 681e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 68237eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 68337eeb815SJakub Kruzik 68437eeb815SJakub Kruzik Level: beginner 68537eeb815SJakub Kruzik 68637eeb815SJakub Kruzik Notes: 687e361f8b5SJakub Kruzik todo 68837eeb815SJakub Kruzik 68937eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 690e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 69137eeb815SJakub Kruzik M*/ 69237eeb815SJakub Kruzik 69337eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 69437eeb815SJakub Kruzik { 69537eeb815SJakub Kruzik PC_Deflation *def; 69637eeb815SJakub Kruzik PetscErrorCode ierr; 69737eeb815SJakub Kruzik 69837eeb815SJakub Kruzik PetscFunctionBegin; 69937eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 70037eeb815SJakub Kruzik pc->data = (void*)def; 70137eeb815SJakub Kruzik 702e662bc50SJakub Kruzik def->init = PETSC_FALSE; 703e662bc50SJakub Kruzik def->pre = PETSC_TRUE; 704e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 705e662bc50SJakub Kruzik def->truenorm = PETSC_TRUE; 706e662bc50SJakub Kruzik def->reductionfact = -1; 707e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 708e662bc50SJakub Kruzik def->spacesize = 1; 709e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 710e662bc50SJakub Kruzik def->nestedlvl = 0; 711e662bc50SJakub Kruzik def->maxnestedlvl = 0; 712e662bc50SJakub Kruzik def->adaptiveconv = PETSC_FALSE; 713e662bc50SJakub Kruzik def->adaptiveconst = 1.0; 71437eeb815SJakub Kruzik 71537eeb815SJakub Kruzik /* 71637eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 71737eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 71837eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 71937eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 72037eeb815SJakub Kruzik not needed. 72137eeb815SJakub Kruzik */ 72237eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 72337eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 724b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 725b27c8092SJakub Kruzik pc->ops->postsolve = 0; 72637eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 72737eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 72837eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 72937eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 7308513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 73137eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 73237eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 73337eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 73437eeb815SJakub Kruzik 73537eeb815SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr); 73637eeb815SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr); 737e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 7385c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 739*22b0793eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 74037eeb815SJakub Kruzik PetscFunctionReturn(0); 74137eeb815SJakub Kruzik } 74237eeb815SJakub Kruzik 743