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 149268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 150268c6673SJakub Kruzik { 151268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 152268c6673SJakub Kruzik 153268c6673SJakub Kruzik PetscFunctionBegin; 154268c6673SJakub Kruzik *apc = def->pc; 155268c6673SJakub Kruzik PetscFunctionReturn(0); 156268c6673SJakub Kruzik } 157268c6673SJakub Kruzik 158268c6673SJakub Kruzik /*@ 159268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 160268c6673SJakub Kruzik 161268c6673SJakub Kruzik Not Collective 162268c6673SJakub Kruzik 163268c6673SJakub Kruzik Input Parameters: 164268c6673SJakub Kruzik . pc - the preconditioner context 165268c6673SJakub Kruzik 166268c6673SJakub Kruzik Output Parameter: 167268c6673SJakub Kruzik . apc - additional preconditioner 168268c6673SJakub Kruzik 169268c6673SJakub Kruzik Level: advanced 170268c6673SJakub Kruzik 171268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 172268c6673SJakub Kruzik @*/ 173268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 174268c6673SJakub Kruzik { 175268c6673SJakub Kruzik PetscErrorCode ierr; 176268c6673SJakub Kruzik 177268c6673SJakub Kruzik PetscFunctionBegin; 178268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 179268c6673SJakub Kruzik PetscValidPointer(pc,2); 180268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 181268c6673SJakub Kruzik PetscFunctionReturn(0); 182268c6673SJakub Kruzik } 183268c6673SJakub 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 200268c6673SJakub Kruzik Collective on PC 20122b0793eSJakub Kruzik 20222b0793eSJakub Kruzik Input Parameters: 20322b0793eSJakub Kruzik + pc - the preconditioner context 20422b0793eSJakub Kruzik - apc - additional preconditioner 20522b0793eSJakub Kruzik 206268c6673SJakub Kruzik Level: developer 20722b0793eSJakub Kruzik 208268c6673SJakub 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 244*f3eaa4f0SJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), 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 257*f3eaa4f0SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 258*f3eaa4f0SJakub Kruzik { 259*f3eaa4f0SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 260*f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 261*f3eaa4f0SJakub Kruzik 262*f3eaa4f0SJakub Kruzik PetscFunctionBegin; 263*f3eaa4f0SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 264*f3eaa4f0SJakub Kruzik def->WtAWinv = ksp; 265*f3eaa4f0SJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 266*f3eaa4f0SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 267*f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 268*f3eaa4f0SJakub Kruzik } 269*f3eaa4f0SJakub Kruzik 270*f3eaa4f0SJakub Kruzik /*@ 271*f3eaa4f0SJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 272*f3eaa4f0SJakub Kruzik 273*f3eaa4f0SJakub Kruzik Collective on PC 274*f3eaa4f0SJakub Kruzik 275*f3eaa4f0SJakub Kruzik Input Parameters: 276*f3eaa4f0SJakub Kruzik + pc - preconditioner context 277*f3eaa4f0SJakub Kruzik - ksp - coarse problem KSP context 278*f3eaa4f0SJakub Kruzik 279*f3eaa4f0SJakub Kruzik Level: developer 280*f3eaa4f0SJakub Kruzik 281*f3eaa4f0SJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 282*f3eaa4f0SJakub Kruzik @*/ 283*f3eaa4f0SJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 284*f3eaa4f0SJakub Kruzik { 285*f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 286*f3eaa4f0SJakub Kruzik 287*f3eaa4f0SJakub Kruzik PetscFunctionBegin; 288*f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 289*f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 290*f3eaa4f0SJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 291*f3eaa4f0SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 292*f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 293*f3eaa4f0SJakub Kruzik } 294*f3eaa4f0SJakub Kruzik 29537eeb815SJakub Kruzik /* 296b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 297b27c8092SJakub Kruzik */ 298b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 299b27c8092SJakub Kruzik { 300b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 301b27c8092SJakub Kruzik Mat A; 302b27c8092SJakub Kruzik Vec r,w1,w2; 303b27c8092SJakub Kruzik PetscBool nonzero; 304b27c8092SJakub Kruzik PetscErrorCode ierr; 30537eeb815SJakub Kruzik 306b27c8092SJakub Kruzik PetscFunctionBegin; 307b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 308b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 309b27c8092SJakub Kruzik r = def->work; 310b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 31137eeb815SJakub Kruzik 312b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 313b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 314b27c8092SJakub Kruzik if (nonzero) { 315b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 316b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 317b27c8092SJakub Kruzik } else { 318b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 319b27c8092SJakub Kruzik } 320b27c8092SJakub Kruzik 321b27c8092SJakub Kruzik ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 322b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 323b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 324b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 325b27c8092SJakub Kruzik PetscFunctionReturn(0); 326b27c8092SJakub Kruzik } 32737eeb815SJakub Kruzik 328f8bf229cSJakub Kruzik /* 329f8bf229cSJakub Kruzik if (def->correct) { 330ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 331f8bf229cSJakub Kruzik } else { 332ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 333f8bf229cSJakub Kruzik } 33437eeb815SJakub Kruzik */ 335f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 336f8bf229cSJakub Kruzik { 337f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 338f8bf229cSJakub Kruzik Mat A; 339f8bf229cSJakub Kruzik Vec u,w1,w2; 340f8bf229cSJakub Kruzik PetscErrorCode ierr; 341f8bf229cSJakub Kruzik 342f8bf229cSJakub Kruzik PetscFunctionBegin; 343f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 344f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 345f8bf229cSJakub Kruzik u = def->work; 346f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 347f8bf229cSJakub Kruzik 348ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 34922b0793eSJakub Kruzik if (!def->init) { 350ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 351f8bf229cSJakub Kruzik if (def->correct) { 352ae029463SJakub Kruzik if (def->Wt) { 353ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 354ae029463SJakub Kruzik } else { 355ae029463SJakub Kruzik ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 356f8bf229cSJakub Kruzik } 357ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 358f8bf229cSJakub Kruzik } 359f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 360f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 36122b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 36222b0793eSJakub Kruzik } 363f8bf229cSJakub Kruzik PetscFunctionReturn(0); 364f8bf229cSJakub Kruzik } 365f8bf229cSJakub Kruzik 36637eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 36737eeb815SJakub Kruzik { 368409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 369409a357bSJakub Kruzik KSP innerksp; 370409a357bSJakub Kruzik PC pcinner; 371409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 372409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 373409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 374409a357bSJakub Kruzik MatCompositeType ctype; 375409a357bSJakub Kruzik MPI_Comm comm; 376409a357bSJakub Kruzik const char *prefix; 37737eeb815SJakub Kruzik PetscErrorCode ierr; 37837eeb815SJakub Kruzik 37937eeb815SJakub Kruzik PetscFunctionBegin; 380409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 381409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 382f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 383f8bf229cSJakub Kruzik 384f8bf229cSJakub Kruzik /* compute a deflation space */ 385409a357bSJakub Kruzik if (def->W || def->Wt) { 386409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 387409a357bSJakub Kruzik } else { 388e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 38937eeb815SJakub Kruzik } 39037eeb815SJakub Kruzik 391409a357bSJakub Kruzik /* nested deflation */ 392409a357bSJakub Kruzik if (def->W) { 393409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 394409a357bSJakub Kruzik if (match) { 395409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 396409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 39737eeb815SJakub Kruzik } 398409a357bSJakub Kruzik } else { 399409a357bSJakub Kruzik ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr); 400409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 401409a357bSJakub Kruzik if (match) { 402409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 403409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 404409a357bSJakub Kruzik } 405409a357bSJakub Kruzik transp = PETSC_TRUE; 406409a357bSJakub Kruzik } 407409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 408409a357bSJakub Kruzik if (!transp) { 40928da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 41028da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 411409a357bSJakub Kruzik for (i=0; i<size; i++) { 412409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 413409a357bSJakub Kruzik } 414409a357bSJakub Kruzik size -= 1; 415409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 416409a357bSJakub Kruzik def->W = mats[size]; 417409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 418409a357bSJakub Kruzik if (size > 1) { 419409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 420409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 421409a357bSJakub Kruzik } else { 422409a357bSJakub Kruzik nextDef = mats[0]; 423409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 424409a357bSJakub Kruzik } 42528da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 426409a357bSJakub Kruzik } else { 427409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 428409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 429409a357bSJakub Kruzik } 43028da0170SJakub Kruzik } else { 43128da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 43228da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 43328da0170SJakub Kruzik for (i=0; i<size; i++) { 43428da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 43528da0170SJakub Kruzik } 43628da0170SJakub Kruzik size -= 1; 43728da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 43828da0170SJakub Kruzik def->Wt = mats[0]; 43928da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 44028da0170SJakub Kruzik if (size > 1) { 44128da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 44228da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 44328da0170SJakub Kruzik } else { 44428da0170SJakub Kruzik nextDef = mats[1]; 44528da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 446409a357bSJakub Kruzik } 447409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 44828da0170SJakub Kruzik } else { 44928da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 45028da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 45128da0170SJakub Kruzik } 45228da0170SJakub Kruzik } 45328da0170SJakub Kruzik } 45428da0170SJakub Kruzik 45528da0170SJakub Kruzik if (transp) { 45628da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 45728da0170SJakub Kruzik ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 458409a357bSJakub Kruzik } 459409a357bSJakub Kruzik 46022b0793eSJakub Kruzik 46122b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 46222b0793eSJakub Kruzik 463ae029463SJakub Kruzik /* assemble WtA */ 464ae029463SJakub Kruzik if (!def->WtA) { 465ae029463SJakub Kruzik if (def->Wt) { 466ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 467ae029463SJakub Kruzik } else { 468ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 469ae029463SJakub Kruzik } 470ae029463SJakub Kruzik } 471409a357bSJakub Kruzik /* setup coarse problem */ 472409a357bSJakub Kruzik if (!def->WtAWinv) { 47328da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 474409a357bSJakub Kruzik if (!def->WtAW) { 475ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 476409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 477409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 478409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 479409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 480ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 481409a357bSJakub Kruzik PetscReal *norms; 482409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 483409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 484409a357bSJakub Kruzik for (i=0; i<m; i++) { 485409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 486409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 487409a357bSJakub Kruzik } 488409a357bSJakub Kruzik } 489409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 490409a357bSJakub Kruzik #endif 491409a357bSJakub Kruzik } else { 492409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 493409a357bSJakub Kruzik } 494409a357bSJakub Kruzik /* TODO use MATINV */ 495409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 496409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 497409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 498557a2f7dSJakub Kruzik /* Setup KSP and PC */ 499557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 500557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 501557a2f7dSJakub Kruzik /* set default KSPtype */ 502557a2f7dSJakub Kruzik if (!def->ksptype) { 503557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 504557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 505557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 506557a2f7dSJakub Kruzik } 507557a2f7dSJakub Kruzik } 508ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 509557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 510557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 511557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 512ae029463SJakub Kruzik /* inherit options */ 513557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 514557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 515557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 516557a2f7dSJakub Kruzik } else { /* the last level */ 517557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 518409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 519ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 520ac66f006SJakub Kruzik if (prefix) { 521ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 522ac66f006SJakub Kruzik } 52322b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 524ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 525409a357bSJakub Kruzik /* Reduction factor choice */ 526409a357bSJakub Kruzik red = def->reductionfact; 527409a357bSJakub Kruzik if (red < 0) { 528409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 529409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 530409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 531409a357bSJakub Kruzik if (match) red = commsize; 532409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 533409a357bSJakub Kruzik } 534409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 535ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 536409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 537ac66f006SJakub Kruzik if (innerksp) { 538409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 539409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 540ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 541409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 542409a357bSJakub Kruzik if (commsize == red) { 543ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 544409a357bSJakub Kruzik } else { 545ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 546409a357bSJakub Kruzik } 547409a357bSJakub Kruzik } 548557a2f7dSJakub Kruzik } 549557a2f7dSJakub Kruzik 550557a2f7dSJakub Kruzik if (innerksp) { 551409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 552409a357bSJakub Kruzik if (prefix) { 553409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 554409a357bSJakub Kruzik } 55522b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 556557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 557557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 558ac66f006SJakub Kruzik } 559409a357bSJakub Kruzik } 560557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 561557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 562409a357bSJakub Kruzik 56322b0793eSJakub Kruzik if (!def->pc) { 56422b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 56522b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 56622b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 56722b0793eSJakub Kruzik if (prefix) { 56822b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 56922b0793eSJakub Kruzik } 57022b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 57122b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 57222b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 57322b0793eSJakub Kruzik } 57422b0793eSJakub Kruzik 575f8bf229cSJakub Kruzik /* create work vecs */ 576f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 577f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 57837eeb815SJakub Kruzik PetscFunctionReturn(0); 57937eeb815SJakub Kruzik } 580b30d4775SJakub Kruzik 58137eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 58237eeb815SJakub Kruzik { 583b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 58437eeb815SJakub Kruzik PetscErrorCode ierr; 58537eeb815SJakub Kruzik 58637eeb815SJakub Kruzik PetscFunctionBegin; 587b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 588b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 589b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 590b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 591ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 592b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 593b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 59422b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 59537eeb815SJakub Kruzik PetscFunctionReturn(0); 59637eeb815SJakub Kruzik } 59737eeb815SJakub Kruzik 59837eeb815SJakub Kruzik /* 59937eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 60037eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 60137eeb815SJakub Kruzik 60237eeb815SJakub Kruzik Input Parameter: 60337eeb815SJakub Kruzik . pc - the preconditioner context 60437eeb815SJakub Kruzik 60537eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 60637eeb815SJakub Kruzik */ 60737eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 60837eeb815SJakub Kruzik { 60937eeb815SJakub Kruzik PetscErrorCode ierr; 61037eeb815SJakub Kruzik 61137eeb815SJakub Kruzik PetscFunctionBegin; 61237eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 61337eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 61437eeb815SJakub Kruzik PetscFunctionReturn(0); 61537eeb815SJakub Kruzik } 61637eeb815SJakub Kruzik 6178513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 6188513960bSJakub Kruzik { 6198513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 6208513960bSJakub Kruzik PetscBool iascii; 6218513960bSJakub Kruzik PetscErrorCode ierr; 6228513960bSJakub Kruzik 6238513960bSJakub Kruzik PetscFunctionBegin; 6248513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6258513960bSJakub Kruzik if (iascii) { 6268513960bSJakub Kruzik //if (cg->adaptiveconv) { 6278513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 6288513960bSJakub Kruzik //} 6298513960bSJakub Kruzik if (def->correct) { 6308513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 6318513960bSJakub Kruzik } 6328513960bSJakub Kruzik if (!def->nestedlvl) { 6338513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 6348513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 6358513960bSJakub Kruzik } 6368513960bSJakub Kruzik 63722b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 63822b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 63922b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 64022b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 64122b0793eSJakub Kruzik 6428513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 6438513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6448513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 6458513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 6468513960bSJakub Kruzik } 6478513960bSJakub Kruzik PetscFunctionReturn(0); 6488513960bSJakub Kruzik } 6498513960bSJakub Kruzik 65037eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 65137eeb815SJakub Kruzik { 652880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 65337eeb815SJakub Kruzik PetscErrorCode ierr; 65437eeb815SJakub Kruzik 65537eeb815SJakub Kruzik PetscFunctionBegin; 65637eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 657880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 658880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 659880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 660880d05ceSJakub Kruzik //TODO add set function and fix manpages 661880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 662fcb31d99SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 663fcb31d99SJakub 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); 664880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 665880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 66637eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 66737eeb815SJakub Kruzik PetscFunctionReturn(0); 66837eeb815SJakub Kruzik } 66937eeb815SJakub Kruzik 67037eeb815SJakub Kruzik /*MC 671e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 672e361f8b5SJakub Kruzik or to a predefined value 67337eeb815SJakub Kruzik 67437eeb815SJakub Kruzik Options Database Key: 675e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 67637eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 67737eeb815SJakub Kruzik 67837eeb815SJakub Kruzik Level: beginner 67937eeb815SJakub Kruzik 68037eeb815SJakub Kruzik Notes: 681e361f8b5SJakub Kruzik todo 68237eeb815SJakub Kruzik 68337eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 684e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 68537eeb815SJakub Kruzik M*/ 68637eeb815SJakub Kruzik 68737eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 68837eeb815SJakub Kruzik { 68937eeb815SJakub Kruzik PC_Deflation *def; 69037eeb815SJakub Kruzik PetscErrorCode ierr; 69137eeb815SJakub Kruzik 69237eeb815SJakub Kruzik PetscFunctionBegin; 69337eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 69437eeb815SJakub Kruzik pc->data = (void*)def; 69537eeb815SJakub Kruzik 696e662bc50SJakub Kruzik def->init = PETSC_FALSE; 697e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 698fcb31d99SJakub Kruzik def->correctfact = 1.0; 699e662bc50SJakub Kruzik def->reductionfact = -1; 700e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 701e662bc50SJakub Kruzik def->spacesize = 1; 702e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 703e662bc50SJakub Kruzik def->nestedlvl = 0; 704e662bc50SJakub Kruzik def->maxnestedlvl = 0; 70537eeb815SJakub Kruzik 70637eeb815SJakub Kruzik /* 70737eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 70837eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 70937eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 71037eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 71137eeb815SJakub Kruzik not needed. 71237eeb815SJakub Kruzik */ 71337eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 71437eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 715b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 716b27c8092SJakub Kruzik pc->ops->postsolve = 0; 71737eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 71837eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 71937eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 72037eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 7218513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 72237eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 72337eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 72437eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 72537eeb815SJakub Kruzik 726e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 7275c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 728268c6673SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 72922b0793eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 7304a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 731*f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 73237eeb815SJakub Kruzik PetscFunctionReturn(0); 73337eeb815SJakub Kruzik } 73437eeb815SJakub Kruzik 735