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 538513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 548513960bSJakub Kruzik "haar", 558513960bSJakub Kruzik "jacket-haar", 568513960bSJakub Kruzik "db2", 578513960bSJakub Kruzik "db4", 588513960bSJakub Kruzik "db8", 598513960bSJakub Kruzik "db16", 608513960bSJakub Kruzik "biorth22", 618513960bSJakub Kruzik "meyer", 628513960bSJakub Kruzik "aggregation", 638513960bSJakub Kruzik "slepc", 648513960bSJakub Kruzik "slepc-cheap", 658513960bSJakub Kruzik "user", 668513960bSJakub Kruzik "DdefSpaceType", 678513960bSJakub Kruzik "Ddef_SPACE_", 688513960bSJakub Kruzik 0 698513960bSJakub Kruzik }; 708513960bSJakub Kruzik 718513960bSJakub Kruzik 72e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 73e662bc50SJakub Kruzik { 74e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 75e662bc50SJakub Kruzik PetscErrorCode ierr; 76e662bc50SJakub Kruzik 77e662bc50SJakub Kruzik PetscFunctionBegin; 78e662bc50SJakub Kruzik if (transpose) { 79e662bc50SJakub Kruzik def->Wt = W; 80e662bc50SJakub Kruzik def->W = NULL; 81e662bc50SJakub Kruzik } else { 82e662bc50SJakub Kruzik def->W = W; 83e662bc50SJakub Kruzik } 84e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 85e662bc50SJakub Kruzik PetscFunctionReturn(0); 86e662bc50SJakub Kruzik } 87e662bc50SJakub Kruzik 88e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 89e662bc50SJakub Kruzik /*@ 90e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 91e662bc50SJakub Kruzik 92e662bc50SJakub Kruzik Logically Collective on PC 93e662bc50SJakub Kruzik 94e662bc50SJakub Kruzik Input Parameters: 95e662bc50SJakub Kruzik + pc - the preconditioner context 96e662bc50SJakub Kruzik . W - deflation matrix 97e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 98e662bc50SJakub Kruzik 99e662bc50SJakub Kruzik Level: intermediate 100e662bc50SJakub Kruzik 101e662bc50SJakub Kruzik .seealso: PCDEFLATION 102e662bc50SJakub Kruzik @*/ 103e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 104e662bc50SJakub Kruzik { 105e662bc50SJakub Kruzik PetscErrorCode ierr; 106e662bc50SJakub Kruzik 107e662bc50SJakub Kruzik PetscFunctionBegin; 108e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 109e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 110e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 111e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 112e662bc50SJakub Kruzik PetscFunctionReturn(0); 113e662bc50SJakub Kruzik } 114e662bc50SJakub Kruzik 1155c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 1165c0c31c5SJakub Kruzik { 1175c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1185c0c31c5SJakub Kruzik 1195c0c31c5SJakub Kruzik PetscFunctionBegin; 12081e2347eSJakub Kruzik if (current) def->nestedlvl = current; 1215c0c31c5SJakub Kruzik def->maxnestedlvl = max; 1225c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1235c0c31c5SJakub Kruzik } 1245c0c31c5SJakub Kruzik 1255c0c31c5SJakub Kruzik /*@ 1265c0c31c5SJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 1275c0c31c5SJakub Kruzik 1285c0c31c5SJakub Kruzik Logically Collective on PC 1295c0c31c5SJakub Kruzik 1305c0c31c5SJakub Kruzik Input Parameters: 1315c0c31c5SJakub Kruzik + pc - the preconditioner context 13222b0793eSJakub Kruzik - max - maximum deflation level 1335c0c31c5SJakub Kruzik 1345c0c31c5SJakub Kruzik Level: intermediate 1355c0c31c5SJakub Kruzik 1365c0c31c5SJakub Kruzik .seealso: PCDEFLATION 1375c0c31c5SJakub Kruzik @*/ 1385c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 1395c0c31c5SJakub Kruzik { 1405c0c31c5SJakub Kruzik PetscErrorCode ierr; 1415c0c31c5SJakub Kruzik 1425c0c31c5SJakub Kruzik PetscFunctionBegin; 14322b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1445c0c31c5SJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 1455c0c31c5SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 1465c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1475c0c31c5SJakub Kruzik } 148e662bc50SJakub Kruzik 149*268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 150*268c6673SJakub Kruzik { 151*268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 152*268c6673SJakub Kruzik 153*268c6673SJakub Kruzik PetscFunctionBegin; 154*268c6673SJakub Kruzik *apc = def->pc; 155*268c6673SJakub Kruzik PetscFunctionReturn(0); 156*268c6673SJakub Kruzik } 157*268c6673SJakub Kruzik 158*268c6673SJakub Kruzik /*@ 159*268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 160*268c6673SJakub Kruzik 161*268c6673SJakub Kruzik Not Collective 162*268c6673SJakub Kruzik 163*268c6673SJakub Kruzik Input Parameters: 164*268c6673SJakub Kruzik . pc - the preconditioner context 165*268c6673SJakub Kruzik 166*268c6673SJakub Kruzik Output Parameter: 167*268c6673SJakub Kruzik . apc - additional preconditioner 168*268c6673SJakub Kruzik 169*268c6673SJakub Kruzik Level: advanced 170*268c6673SJakub Kruzik 171*268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 172*268c6673SJakub Kruzik @*/ 173*268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 174*268c6673SJakub Kruzik { 175*268c6673SJakub Kruzik PetscErrorCode ierr; 176*268c6673SJakub Kruzik 177*268c6673SJakub Kruzik PetscFunctionBegin; 178*268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 179*268c6673SJakub Kruzik PetscValidPointer(pc,2); 180*268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 181*268c6673SJakub Kruzik PetscFunctionReturn(0); 182*268c6673SJakub Kruzik } 183*268c6673SJakub Kruzik 18422b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 18522b0793eSJakub Kruzik { 18622b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 18722b0793eSJakub Kruzik PetscErrorCode ierr; 18822b0793eSJakub Kruzik 18922b0793eSJakub Kruzik PetscFunctionBegin; 19022b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 19122b0793eSJakub Kruzik def->pc = apc; 19222b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 19322b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 19422b0793eSJakub Kruzik PetscFunctionReturn(0); 19522b0793eSJakub Kruzik } 19622b0793eSJakub Kruzik 19722b0793eSJakub Kruzik /*@ 19822b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 19922b0793eSJakub Kruzik 200*268c6673SJakub Kruzik Collective on PC 20122b0793eSJakub Kruzik 20222b0793eSJakub Kruzik Input Parameters: 20322b0793eSJakub Kruzik + pc - the preconditioner context 20422b0793eSJakub Kruzik - apc - additional preconditioner 20522b0793eSJakub Kruzik 206*268c6673SJakub Kruzik Level: developer 20722b0793eSJakub Kruzik 208*268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 20922b0793eSJakub Kruzik @*/ 21022b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 21122b0793eSJakub Kruzik { 21222b0793eSJakub Kruzik PetscErrorCode ierr; 21322b0793eSJakub Kruzik 21422b0793eSJakub Kruzik PetscFunctionBegin; 21522b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 21622b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 21722b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 21822b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 21922b0793eSJakub Kruzik PetscFunctionReturn(0); 22022b0793eSJakub Kruzik } 22122b0793eSJakub Kruzik 2224a99276eSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 2234a99276eSJakub Kruzik { 2244a99276eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2254a99276eSJakub Kruzik 2264a99276eSJakub Kruzik PetscFunctionBegin; 2274a99276eSJakub Kruzik *ksp = def->WtAWinv; 2284a99276eSJakub Kruzik PetscFunctionReturn(0); 2294a99276eSJakub Kruzik } 2304a99276eSJakub Kruzik 2314a99276eSJakub Kruzik /*@ 2324a99276eSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 2334a99276eSJakub Kruzik 2344a99276eSJakub Kruzik Not Collective 2354a99276eSJakub Kruzik 2364a99276eSJakub Kruzik Input Parameters: 2374a99276eSJakub Kruzik . pc - preconditioner context 2384a99276eSJakub Kruzik 2394a99276eSJakub Kruzik Output Parameter: 2404a99276eSJakub Kruzik . ksp - coarse problem KSP context 2414a99276eSJakub Kruzik 2424a99276eSJakub Kruzik Level: developer 2434a99276eSJakub Kruzik 2444a99276eSJakub Kruzik .seealso: PCDEFLATION 2454a99276eSJakub Kruzik @*/ 2464a99276eSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 2474a99276eSJakub Kruzik { 2484a99276eSJakub Kruzik PetscErrorCode ierr; 2494a99276eSJakub Kruzik 2504a99276eSJakub Kruzik PetscFunctionBegin; 2514a99276eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2524a99276eSJakub Kruzik PetscValidPointer(ksp,2); 2534a99276eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 2544a99276eSJakub Kruzik PetscFunctionReturn(0); 2554a99276eSJakub Kruzik } 2564a99276eSJakub Kruzik 25737eeb815SJakub Kruzik /* 258b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 259b27c8092SJakub Kruzik */ 260b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 261b27c8092SJakub Kruzik { 262b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 263b27c8092SJakub Kruzik Mat A; 264b27c8092SJakub Kruzik Vec r,w1,w2; 265b27c8092SJakub Kruzik PetscBool nonzero; 266b27c8092SJakub Kruzik PetscErrorCode ierr; 26737eeb815SJakub Kruzik 268b27c8092SJakub Kruzik PetscFunctionBegin; 269b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 270b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 271b27c8092SJakub Kruzik r = def->work; 272b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 27337eeb815SJakub Kruzik 274b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 275b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 276b27c8092SJakub Kruzik if (nonzero) { 277b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 278b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 279b27c8092SJakub Kruzik } else { 280b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 281b27c8092SJakub Kruzik } 282b27c8092SJakub Kruzik 283b27c8092SJakub Kruzik ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 284b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 285b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 286b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 287b27c8092SJakub Kruzik PetscFunctionReturn(0); 288b27c8092SJakub Kruzik } 28937eeb815SJakub Kruzik 290f8bf229cSJakub Kruzik /* 291f8bf229cSJakub Kruzik if (def->correct) { 292ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 293f8bf229cSJakub Kruzik } else { 294ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 295f8bf229cSJakub Kruzik } 29637eeb815SJakub Kruzik */ 297f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 298f8bf229cSJakub Kruzik { 299f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 300f8bf229cSJakub Kruzik Mat A; 301f8bf229cSJakub Kruzik Vec u,w1,w2; 302f8bf229cSJakub Kruzik PetscErrorCode ierr; 303f8bf229cSJakub Kruzik 304f8bf229cSJakub Kruzik PetscFunctionBegin; 305f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 306f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 307f8bf229cSJakub Kruzik u = def->work; 308f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 309f8bf229cSJakub Kruzik 310ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 31122b0793eSJakub Kruzik if (!def->init) { 312ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 313f8bf229cSJakub Kruzik if (def->correct) { 314ae029463SJakub Kruzik if (def->Wt) { 315ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 316ae029463SJakub Kruzik } else { 317ae029463SJakub Kruzik ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 318f8bf229cSJakub Kruzik } 319ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 320f8bf229cSJakub Kruzik } 321f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 322f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 32322b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 32422b0793eSJakub Kruzik } 325f8bf229cSJakub Kruzik PetscFunctionReturn(0); 326f8bf229cSJakub Kruzik } 327f8bf229cSJakub Kruzik 32837eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 32937eeb815SJakub Kruzik { 330409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 331409a357bSJakub Kruzik KSP innerksp; 332409a357bSJakub Kruzik PC pcinner; 333409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 334409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 335409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 336409a357bSJakub Kruzik MatCompositeType ctype; 337409a357bSJakub Kruzik MPI_Comm comm; 338409a357bSJakub Kruzik const char *prefix; 33937eeb815SJakub Kruzik PetscErrorCode ierr; 34037eeb815SJakub Kruzik 34137eeb815SJakub Kruzik PetscFunctionBegin; 342409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 343409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 344f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 345f8bf229cSJakub Kruzik 346f8bf229cSJakub Kruzik /* compute a deflation space */ 347409a357bSJakub Kruzik if (def->W || def->Wt) { 348409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 349409a357bSJakub Kruzik } else { 350e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 35137eeb815SJakub Kruzik } 35237eeb815SJakub Kruzik 353409a357bSJakub Kruzik /* nested deflation */ 354409a357bSJakub Kruzik if (def->W) { 355409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 356409a357bSJakub Kruzik if (match) { 357409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 358409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 35937eeb815SJakub Kruzik } 360409a357bSJakub Kruzik } else { 361409a357bSJakub Kruzik ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr); 362409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 363409a357bSJakub Kruzik if (match) { 364409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 365409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 366409a357bSJakub Kruzik } 367409a357bSJakub Kruzik transp = PETSC_TRUE; 368409a357bSJakub Kruzik } 369409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 370409a357bSJakub Kruzik if (!transp) { 37128da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 37228da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 373409a357bSJakub Kruzik for (i=0; i<size; i++) { 374409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 375409a357bSJakub Kruzik } 376409a357bSJakub Kruzik size -= 1; 377409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 378409a357bSJakub Kruzik def->W = mats[size]; 379409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 380409a357bSJakub Kruzik if (size > 1) { 381409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 382409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 383409a357bSJakub Kruzik } else { 384409a357bSJakub Kruzik nextDef = mats[0]; 385409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 386409a357bSJakub Kruzik } 38728da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 388409a357bSJakub Kruzik } else { 389409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 390409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 391409a357bSJakub Kruzik } 39228da0170SJakub Kruzik } else { 39328da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 39428da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 39528da0170SJakub Kruzik for (i=0; i<size; i++) { 39628da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 39728da0170SJakub Kruzik } 39828da0170SJakub Kruzik size -= 1; 39928da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 40028da0170SJakub Kruzik def->Wt = mats[0]; 40128da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 40228da0170SJakub Kruzik if (size > 1) { 40328da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 40428da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 40528da0170SJakub Kruzik } else { 40628da0170SJakub Kruzik nextDef = mats[1]; 40728da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 408409a357bSJakub Kruzik } 409409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 41028da0170SJakub Kruzik } else { 41128da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 41228da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 41328da0170SJakub Kruzik } 41428da0170SJakub Kruzik } 41528da0170SJakub Kruzik } 41628da0170SJakub Kruzik 41728da0170SJakub Kruzik if (transp) { 41828da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 41928da0170SJakub Kruzik ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 420409a357bSJakub Kruzik } 421409a357bSJakub Kruzik 42222b0793eSJakub Kruzik 42322b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 42422b0793eSJakub Kruzik 425ae029463SJakub Kruzik /* assemble WtA */ 426ae029463SJakub Kruzik if (!def->WtA) { 427ae029463SJakub Kruzik if (def->Wt) { 428ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 429ae029463SJakub Kruzik } else { 430ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 431ae029463SJakub Kruzik } 432ae029463SJakub Kruzik } 433409a357bSJakub Kruzik /* setup coarse problem */ 434409a357bSJakub Kruzik if (!def->WtAWinv) { 43528da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 436409a357bSJakub Kruzik if (!def->WtAW) { 437ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 438409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 439409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 440409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 441409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 442ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 443409a357bSJakub Kruzik PetscReal *norms; 444409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 445409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 446409a357bSJakub Kruzik for (i=0; i<m; i++) { 447409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 448409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 449409a357bSJakub Kruzik } 450409a357bSJakub Kruzik } 451409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 452409a357bSJakub Kruzik #endif 453409a357bSJakub Kruzik } else { 454409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 455409a357bSJakub Kruzik } 456409a357bSJakub Kruzik /* TODO use MATINV */ 457409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 458409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 459409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 460557a2f7dSJakub Kruzik /* Setup KSP and PC */ 461557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 462557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 463557a2f7dSJakub Kruzik /* set default KSPtype */ 464557a2f7dSJakub Kruzik if (!def->ksptype) { 465557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 466557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 467557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 468557a2f7dSJakub Kruzik } 469557a2f7dSJakub Kruzik } 470ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 471557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 472557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 473557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 474ae029463SJakub Kruzik /* inherit options */ 475557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 476557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 477557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 478557a2f7dSJakub Kruzik } else { /* the last level */ 479557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 480409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 481ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 482ac66f006SJakub Kruzik if (prefix) { 483ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 484ac66f006SJakub Kruzik } 48522b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 486ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 487409a357bSJakub Kruzik /* Reduction factor choice */ 488409a357bSJakub Kruzik red = def->reductionfact; 489409a357bSJakub Kruzik if (red < 0) { 490409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 491409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 492409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 493409a357bSJakub Kruzik if (match) red = commsize; 494409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 495409a357bSJakub Kruzik } 496409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 497ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 498409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 499ac66f006SJakub Kruzik if (innerksp) { 500409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 501409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 502ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 503409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 504409a357bSJakub Kruzik if (commsize == red) { 505ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 506409a357bSJakub Kruzik } else { 507ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 508409a357bSJakub Kruzik } 509409a357bSJakub Kruzik } 510557a2f7dSJakub Kruzik } 511557a2f7dSJakub Kruzik 512557a2f7dSJakub Kruzik if (innerksp) { 513409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 514409a357bSJakub Kruzik if (prefix) { 515409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 516409a357bSJakub Kruzik } 51722b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 518557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 519557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 520ac66f006SJakub Kruzik } 521409a357bSJakub Kruzik } 522557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 523557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 524409a357bSJakub Kruzik 52522b0793eSJakub Kruzik if (!def->pc) { 52622b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 52722b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 52822b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 52922b0793eSJakub Kruzik if (prefix) { 53022b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 53122b0793eSJakub Kruzik } 53222b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 53322b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 53422b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 53522b0793eSJakub Kruzik } 53622b0793eSJakub Kruzik 537f8bf229cSJakub Kruzik /* create work vecs */ 538f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 539f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 54037eeb815SJakub Kruzik PetscFunctionReturn(0); 54137eeb815SJakub Kruzik } 542b30d4775SJakub Kruzik 54337eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 54437eeb815SJakub Kruzik { 545b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 54637eeb815SJakub Kruzik PetscErrorCode ierr; 54737eeb815SJakub Kruzik 54837eeb815SJakub Kruzik PetscFunctionBegin; 549b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 550b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 551b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 552b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 553ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 554b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 555b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 55622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 55737eeb815SJakub Kruzik PetscFunctionReturn(0); 55837eeb815SJakub Kruzik } 55937eeb815SJakub Kruzik 56037eeb815SJakub Kruzik /* 56137eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 56237eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 56337eeb815SJakub Kruzik 56437eeb815SJakub Kruzik Input Parameter: 56537eeb815SJakub Kruzik . pc - the preconditioner context 56637eeb815SJakub Kruzik 56737eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 56837eeb815SJakub Kruzik */ 56937eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 57037eeb815SJakub Kruzik { 57137eeb815SJakub Kruzik PetscErrorCode ierr; 57237eeb815SJakub Kruzik 57337eeb815SJakub Kruzik PetscFunctionBegin; 57437eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 57537eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 57637eeb815SJakub Kruzik PetscFunctionReturn(0); 57737eeb815SJakub Kruzik } 57837eeb815SJakub Kruzik 5798513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 5808513960bSJakub Kruzik { 5818513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 5828513960bSJakub Kruzik PetscBool iascii; 5838513960bSJakub Kruzik PetscErrorCode ierr; 5848513960bSJakub Kruzik 5858513960bSJakub Kruzik PetscFunctionBegin; 5868513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5878513960bSJakub Kruzik if (iascii) { 5888513960bSJakub Kruzik //if (cg->adaptiveconv) { 5898513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 5908513960bSJakub Kruzik //} 5918513960bSJakub Kruzik if (def->correct) { 5928513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 5938513960bSJakub Kruzik } 5948513960bSJakub Kruzik if (!def->nestedlvl) { 5958513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 5968513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 5978513960bSJakub Kruzik } 5988513960bSJakub Kruzik 59922b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 60022b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 60122b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 60222b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 60322b0793eSJakub Kruzik 6048513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 6058513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6068513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 6078513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 6088513960bSJakub Kruzik } 6098513960bSJakub Kruzik PetscFunctionReturn(0); 6108513960bSJakub Kruzik } 6118513960bSJakub Kruzik 61237eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 61337eeb815SJakub Kruzik { 614880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 61537eeb815SJakub Kruzik PetscErrorCode ierr; 61637eeb815SJakub Kruzik 61737eeb815SJakub Kruzik PetscFunctionBegin; 61837eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 619880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 620880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 621880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 622880d05ceSJakub Kruzik //TODO add set function and fix manpages 623880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 624fcb31d99SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 625fcb31d99SJakub Kruzik ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr); 626880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 627880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 62837eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 62937eeb815SJakub Kruzik PetscFunctionReturn(0); 63037eeb815SJakub Kruzik } 63137eeb815SJakub Kruzik 63237eeb815SJakub Kruzik /*MC 633e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 634e361f8b5SJakub Kruzik or to a predefined value 63537eeb815SJakub Kruzik 63637eeb815SJakub Kruzik Options Database Key: 637e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 63837eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 63937eeb815SJakub Kruzik 64037eeb815SJakub Kruzik Level: beginner 64137eeb815SJakub Kruzik 64237eeb815SJakub Kruzik Notes: 643e361f8b5SJakub Kruzik todo 64437eeb815SJakub Kruzik 64537eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 646e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 64737eeb815SJakub Kruzik M*/ 64837eeb815SJakub Kruzik 64937eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 65037eeb815SJakub Kruzik { 65137eeb815SJakub Kruzik PC_Deflation *def; 65237eeb815SJakub Kruzik PetscErrorCode ierr; 65337eeb815SJakub Kruzik 65437eeb815SJakub Kruzik PetscFunctionBegin; 65537eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 65637eeb815SJakub Kruzik pc->data = (void*)def; 65737eeb815SJakub Kruzik 658e662bc50SJakub Kruzik def->init = PETSC_FALSE; 659e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 660fcb31d99SJakub Kruzik def->correctfact = 1.0; 661e662bc50SJakub Kruzik def->reductionfact = -1; 662e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 663e662bc50SJakub Kruzik def->spacesize = 1; 664e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 665e662bc50SJakub Kruzik def->nestedlvl = 0; 666e662bc50SJakub Kruzik def->maxnestedlvl = 0; 66737eeb815SJakub Kruzik 66837eeb815SJakub Kruzik /* 66937eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 67037eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 67137eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 67237eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 67337eeb815SJakub Kruzik not needed. 67437eeb815SJakub Kruzik */ 67537eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 67637eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 677b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 678b27c8092SJakub Kruzik pc->ops->postsolve = 0; 67937eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 68037eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 68137eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 68237eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 6838513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 68437eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 68537eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 68637eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 68737eeb815SJakub Kruzik 688e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 6895c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 690*268c6673SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 69122b0793eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 6924a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 69337eeb815SJakub Kruzik PetscFunctionReturn(0); 69437eeb815SJakub Kruzik } 69537eeb815SJakub Kruzik 696