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 115*3cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 116*3cf3a049SJakub Kruzik { 117*3cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 118*3cf3a049SJakub Kruzik PetscErrorCode ierr; 119*3cf3a049SJakub Kruzik 120*3cf3a049SJakub Kruzik PetscFunctionBegin; 121*3cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 122*3cf3a049SJakub Kruzik def->WtA = mat; 123*3cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 124*3cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 125*3cf3a049SJakub Kruzik PetscFunctionReturn(0); 126*3cf3a049SJakub Kruzik } 127*3cf3a049SJakub Kruzik 128*3cf3a049SJakub Kruzik /*@ 129*3cf3a049SJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 130*3cf3a049SJakub Kruzik 131*3cf3a049SJakub Kruzik Not Collective 132*3cf3a049SJakub Kruzik 133*3cf3a049SJakub Kruzik Input Parameters: 134*3cf3a049SJakub Kruzik + pc - preconditioner context 135*3cf3a049SJakub Kruzik - mat - projection null space matrix 136*3cf3a049SJakub Kruzik 137*3cf3a049SJakub Kruzik Level: developer 138*3cf3a049SJakub Kruzik 139*3cf3a049SJakub Kruzik .seealso: PCDEFLATION 140*3cf3a049SJakub Kruzik @*/ 141*3cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 142*3cf3a049SJakub Kruzik { 143*3cf3a049SJakub Kruzik PetscErrorCode ierr; 144*3cf3a049SJakub Kruzik 145*3cf3a049SJakub Kruzik PetscFunctionBegin; 146*3cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 147*3cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 148*3cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 149*3cf3a049SJakub Kruzik PetscFunctionReturn(0); 150*3cf3a049SJakub Kruzik } 151*3cf3a049SJakub Kruzik 152e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 153e924b002SJakub Kruzik { 154e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 155e924b002SJakub Kruzik PetscErrorCode ierr; 156e924b002SJakub Kruzik 157e924b002SJakub Kruzik PetscFunctionBegin; 158e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 159e924b002SJakub Kruzik def->WtAW = mat; 160e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 161e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 162e924b002SJakub Kruzik PetscFunctionReturn(0); 163e924b002SJakub Kruzik } 164e924b002SJakub Kruzik 165e924b002SJakub Kruzik /*@ 166e924b002SJakub Kruzik PCDeflationSetCoarseMat - Set coarse problem Mat. 167e924b002SJakub Kruzik 168e924b002SJakub Kruzik Not Collective 169e924b002SJakub Kruzik 170e924b002SJakub Kruzik Input Parameters: 171e924b002SJakub Kruzik + pc - preconditioner context 172e924b002SJakub Kruzik - mat - coarse problem mat 173e924b002SJakub Kruzik 174e924b002SJakub Kruzik Level: developer 175e924b002SJakub Kruzik 176e924b002SJakub Kruzik .seealso: PCDEFLATION 177e924b002SJakub Kruzik @*/ 178e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 179e924b002SJakub Kruzik { 180e924b002SJakub Kruzik PetscErrorCode ierr; 181e924b002SJakub Kruzik 182e924b002SJakub Kruzik PetscFunctionBegin; 183e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 184e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 185e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 186e924b002SJakub Kruzik PetscFunctionReturn(0); 187e924b002SJakub Kruzik } 188e924b002SJakub Kruzik 1895c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 1905c0c31c5SJakub Kruzik { 1915c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1925c0c31c5SJakub Kruzik 1935c0c31c5SJakub Kruzik PetscFunctionBegin; 19481e2347eSJakub Kruzik if (current) def->nestedlvl = current; 1955c0c31c5SJakub Kruzik def->maxnestedlvl = max; 1965c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1975c0c31c5SJakub Kruzik } 1985c0c31c5SJakub Kruzik 1995c0c31c5SJakub Kruzik /*@ 2005c0c31c5SJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 2015c0c31c5SJakub Kruzik 2025c0c31c5SJakub Kruzik Logically Collective on PC 2035c0c31c5SJakub Kruzik 2045c0c31c5SJakub Kruzik Input Parameters: 2055c0c31c5SJakub Kruzik + pc - the preconditioner context 20622b0793eSJakub Kruzik - max - maximum deflation level 2075c0c31c5SJakub Kruzik 2085c0c31c5SJakub Kruzik Level: intermediate 2095c0c31c5SJakub Kruzik 2105c0c31c5SJakub Kruzik .seealso: PCDEFLATION 2115c0c31c5SJakub Kruzik @*/ 2125c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 2135c0c31c5SJakub Kruzik { 2145c0c31c5SJakub Kruzik PetscErrorCode ierr; 2155c0c31c5SJakub Kruzik 2165c0c31c5SJakub Kruzik PetscFunctionBegin; 21722b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2185c0c31c5SJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 2195c0c31c5SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 2205c0c31c5SJakub Kruzik PetscFunctionReturn(0); 2215c0c31c5SJakub Kruzik } 222e662bc50SJakub Kruzik 223268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 224268c6673SJakub Kruzik { 225268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 226268c6673SJakub Kruzik 227268c6673SJakub Kruzik PetscFunctionBegin; 228268c6673SJakub Kruzik *apc = def->pc; 229268c6673SJakub Kruzik PetscFunctionReturn(0); 230268c6673SJakub Kruzik } 231268c6673SJakub Kruzik 232268c6673SJakub Kruzik /*@ 233268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 234268c6673SJakub Kruzik 235268c6673SJakub Kruzik Not Collective 236268c6673SJakub Kruzik 237268c6673SJakub Kruzik Input Parameters: 238268c6673SJakub Kruzik . pc - the preconditioner context 239268c6673SJakub Kruzik 240268c6673SJakub Kruzik Output Parameter: 241268c6673SJakub Kruzik . apc - additional preconditioner 242268c6673SJakub Kruzik 243268c6673SJakub Kruzik Level: advanced 244268c6673SJakub Kruzik 245268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 246268c6673SJakub Kruzik @*/ 247268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 248268c6673SJakub Kruzik { 249268c6673SJakub Kruzik PetscErrorCode ierr; 250268c6673SJakub Kruzik 251268c6673SJakub Kruzik PetscFunctionBegin; 252268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 253268c6673SJakub Kruzik PetscValidPointer(pc,2); 254268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 255268c6673SJakub Kruzik PetscFunctionReturn(0); 256268c6673SJakub Kruzik } 257268c6673SJakub Kruzik 25822b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 25922b0793eSJakub Kruzik { 26022b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 26122b0793eSJakub Kruzik PetscErrorCode ierr; 26222b0793eSJakub Kruzik 26322b0793eSJakub Kruzik PetscFunctionBegin; 26422b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 26522b0793eSJakub Kruzik def->pc = apc; 26622b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 26722b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 26822b0793eSJakub Kruzik PetscFunctionReturn(0); 26922b0793eSJakub Kruzik } 27022b0793eSJakub Kruzik 27122b0793eSJakub Kruzik /*@ 27222b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 27322b0793eSJakub Kruzik 274268c6673SJakub Kruzik Collective on PC 27522b0793eSJakub Kruzik 27622b0793eSJakub Kruzik Input Parameters: 27722b0793eSJakub Kruzik + pc - the preconditioner context 27822b0793eSJakub Kruzik - apc - additional preconditioner 27922b0793eSJakub Kruzik 280268c6673SJakub Kruzik Level: developer 28122b0793eSJakub Kruzik 282268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 28322b0793eSJakub Kruzik @*/ 28422b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 28522b0793eSJakub Kruzik { 28622b0793eSJakub Kruzik PetscErrorCode ierr; 28722b0793eSJakub Kruzik 28822b0793eSJakub Kruzik PetscFunctionBegin; 28922b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 29022b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 29122b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 29222b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 29322b0793eSJakub Kruzik PetscFunctionReturn(0); 29422b0793eSJakub Kruzik } 29522b0793eSJakub Kruzik 2964a99276eSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 2974a99276eSJakub Kruzik { 2984a99276eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2994a99276eSJakub Kruzik 3004a99276eSJakub Kruzik PetscFunctionBegin; 3014a99276eSJakub Kruzik *ksp = def->WtAWinv; 3024a99276eSJakub Kruzik PetscFunctionReturn(0); 3034a99276eSJakub Kruzik } 3044a99276eSJakub Kruzik 3054a99276eSJakub Kruzik /*@ 3064a99276eSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 3074a99276eSJakub Kruzik 3084a99276eSJakub Kruzik Not Collective 3094a99276eSJakub Kruzik 3104a99276eSJakub Kruzik Input Parameters: 3114a99276eSJakub Kruzik . pc - preconditioner context 3124a99276eSJakub Kruzik 3134a99276eSJakub Kruzik Output Parameter: 3144a99276eSJakub Kruzik . ksp - coarse problem KSP context 3154a99276eSJakub Kruzik 3164a99276eSJakub Kruzik Level: developer 3174a99276eSJakub Kruzik 318f3eaa4f0SJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 3194a99276eSJakub Kruzik @*/ 3204a99276eSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3214a99276eSJakub Kruzik { 3224a99276eSJakub Kruzik PetscErrorCode ierr; 3234a99276eSJakub Kruzik 3244a99276eSJakub Kruzik PetscFunctionBegin; 3254a99276eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3264a99276eSJakub Kruzik PetscValidPointer(ksp,2); 3274a99276eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 3284a99276eSJakub Kruzik PetscFunctionReturn(0); 3294a99276eSJakub Kruzik } 3304a99276eSJakub Kruzik 331f3eaa4f0SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 332f3eaa4f0SJakub Kruzik { 333f3eaa4f0SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 334f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 335f3eaa4f0SJakub Kruzik 336f3eaa4f0SJakub Kruzik PetscFunctionBegin; 337f3eaa4f0SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 338f3eaa4f0SJakub Kruzik def->WtAWinv = ksp; 339f3eaa4f0SJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 340f3eaa4f0SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 341f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 342f3eaa4f0SJakub Kruzik } 343f3eaa4f0SJakub Kruzik 344f3eaa4f0SJakub Kruzik /*@ 345f3eaa4f0SJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 346f3eaa4f0SJakub Kruzik 347f3eaa4f0SJakub Kruzik Collective on PC 348f3eaa4f0SJakub Kruzik 349f3eaa4f0SJakub Kruzik Input Parameters: 350f3eaa4f0SJakub Kruzik + pc - preconditioner context 351f3eaa4f0SJakub Kruzik - ksp - coarse problem KSP context 352f3eaa4f0SJakub Kruzik 353f3eaa4f0SJakub Kruzik Level: developer 354f3eaa4f0SJakub Kruzik 355f3eaa4f0SJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 356f3eaa4f0SJakub Kruzik @*/ 357f3eaa4f0SJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 358f3eaa4f0SJakub Kruzik { 359f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 360f3eaa4f0SJakub Kruzik 361f3eaa4f0SJakub Kruzik PetscFunctionBegin; 362f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 363f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 364f3eaa4f0SJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 365f3eaa4f0SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 366f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 367f3eaa4f0SJakub Kruzik } 368f3eaa4f0SJakub Kruzik 36937eeb815SJakub Kruzik /* 370b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 371b27c8092SJakub Kruzik */ 372b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 373b27c8092SJakub Kruzik { 374b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 375b27c8092SJakub Kruzik Mat A; 376b27c8092SJakub Kruzik Vec r,w1,w2; 377b27c8092SJakub Kruzik PetscBool nonzero; 378b27c8092SJakub Kruzik PetscErrorCode ierr; 37937eeb815SJakub Kruzik 380b27c8092SJakub Kruzik PetscFunctionBegin; 381b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 382b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 383b27c8092SJakub Kruzik r = def->work; 384b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 38537eeb815SJakub Kruzik 386b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 387b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 388b27c8092SJakub Kruzik if (nonzero) { 389b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 390b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 391b27c8092SJakub Kruzik } else { 392b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 393b27c8092SJakub Kruzik } 394b27c8092SJakub Kruzik 395a3177931SJakub Kruzik if (def->Wt) { 396a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 397a3177931SJakub Kruzik } else { 398a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 399a3177931SJakub Kruzik } 400b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 401b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 402b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 403b27c8092SJakub Kruzik PetscFunctionReturn(0); 404b27c8092SJakub Kruzik } 40537eeb815SJakub Kruzik 406f8bf229cSJakub Kruzik /* 407f8bf229cSJakub Kruzik if (def->correct) { 408ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 409f8bf229cSJakub Kruzik } else { 410ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 411f8bf229cSJakub Kruzik } 41237eeb815SJakub Kruzik */ 413f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 414f8bf229cSJakub Kruzik { 415f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 416f8bf229cSJakub Kruzik Mat A; 417f8bf229cSJakub Kruzik Vec u,w1,w2; 418f8bf229cSJakub Kruzik PetscErrorCode ierr; 419f8bf229cSJakub Kruzik 420f8bf229cSJakub Kruzik PetscFunctionBegin; 421f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 422f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 423f8bf229cSJakub Kruzik u = def->work; 424f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 425f8bf229cSJakub Kruzik 426ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 42722b0793eSJakub Kruzik if (!def->init) { 428ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 429f8bf229cSJakub Kruzik if (def->correct) { 430ae029463SJakub Kruzik if (def->Wt) { 431ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 432ae029463SJakub Kruzik } else { 433a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 434f8bf229cSJakub Kruzik } 435ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 436f8bf229cSJakub Kruzik } 437f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 438f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 43922b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 44022b0793eSJakub Kruzik } 441f8bf229cSJakub Kruzik PetscFunctionReturn(0); 442f8bf229cSJakub Kruzik } 443f8bf229cSJakub Kruzik 44437eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 44537eeb815SJakub Kruzik { 446409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 447409a357bSJakub Kruzik KSP innerksp; 448409a357bSJakub Kruzik PC pcinner; 449409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 450409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 451409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 452409a357bSJakub Kruzik MatCompositeType ctype; 453409a357bSJakub Kruzik MPI_Comm comm; 454409a357bSJakub Kruzik const char *prefix; 45537eeb815SJakub Kruzik PetscErrorCode ierr; 45637eeb815SJakub Kruzik 45737eeb815SJakub Kruzik PetscFunctionBegin; 458409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 459409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 460f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 461f8bf229cSJakub Kruzik 462f8bf229cSJakub Kruzik /* compute a deflation space */ 463409a357bSJakub Kruzik if (def->W || def->Wt) { 464409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 465409a357bSJakub Kruzik } else { 466e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 46737eeb815SJakub Kruzik } 46837eeb815SJakub Kruzik 469409a357bSJakub Kruzik /* nested deflation */ 470409a357bSJakub Kruzik if (def->W) { 471409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 472409a357bSJakub Kruzik if (match) { 473409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 474409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 47537eeb815SJakub Kruzik } 476409a357bSJakub Kruzik } else { 477a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 478409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 479409a357bSJakub Kruzik if (match) { 480409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 481409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 482409a357bSJakub Kruzik } 483409a357bSJakub Kruzik transp = PETSC_TRUE; 484409a357bSJakub Kruzik } 485409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 486409a357bSJakub Kruzik if (!transp) { 48728da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 48828da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 489409a357bSJakub Kruzik for (i=0; i<size; i++) { 490409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 491409a357bSJakub Kruzik } 492409a357bSJakub Kruzik size -= 1; 493409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 494409a357bSJakub Kruzik def->W = mats[size]; 495409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 496409a357bSJakub Kruzik if (size > 1) { 497409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 498409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 499409a357bSJakub Kruzik } else { 500409a357bSJakub Kruzik nextDef = mats[0]; 501409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 502409a357bSJakub Kruzik } 50328da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 504409a357bSJakub Kruzik } else { 505409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 506409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 507409a357bSJakub Kruzik } 50828da0170SJakub Kruzik } else { 50928da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 51028da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 51128da0170SJakub Kruzik for (i=0; i<size; i++) { 51228da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 51328da0170SJakub Kruzik } 51428da0170SJakub Kruzik size -= 1; 51528da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 51628da0170SJakub Kruzik def->Wt = mats[0]; 51728da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 51828da0170SJakub Kruzik if (size > 1) { 51928da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 52028da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 52128da0170SJakub Kruzik } else { 52228da0170SJakub Kruzik nextDef = mats[1]; 52328da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 524409a357bSJakub Kruzik } 525409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 52628da0170SJakub Kruzik } else { 52728da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 52828da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 52928da0170SJakub Kruzik } 53028da0170SJakub Kruzik } 53128da0170SJakub Kruzik } 53228da0170SJakub Kruzik 53328da0170SJakub Kruzik if (transp) { 53428da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 535a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 536409a357bSJakub Kruzik } 537409a357bSJakub Kruzik 53822b0793eSJakub Kruzik 53922b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 54022b0793eSJakub Kruzik 541ae029463SJakub Kruzik /* assemble WtA */ 542ae029463SJakub Kruzik if (!def->WtA) { 543ae029463SJakub Kruzik if (def->Wt) { 544ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 545ae029463SJakub Kruzik } else { 546a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 547a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 548a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 549a3177931SJakub Kruzik #else 550ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 551a3177931SJakub Kruzik #endif 552ae029463SJakub Kruzik } 553ae029463SJakub Kruzik } 554409a357bSJakub Kruzik /* setup coarse problem */ 555409a357bSJakub Kruzik if (!def->WtAWinv) { 55628da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 557409a357bSJakub Kruzik if (!def->WtAW) { 558ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 559409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 560409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 561409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 562409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 563ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 564409a357bSJakub Kruzik PetscReal *norms; 565409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 566409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 567409a357bSJakub Kruzik for (i=0; i<m; i++) { 568409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 569409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 570409a357bSJakub Kruzik } 571409a357bSJakub Kruzik } 572409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 573409a357bSJakub Kruzik #endif 574409a357bSJakub Kruzik } else { 575409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 576409a357bSJakub Kruzik } 577409a357bSJakub Kruzik /* TODO use MATINV */ 578409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 579409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 580409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 581557a2f7dSJakub Kruzik /* Setup KSP and PC */ 582557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 583557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 584557a2f7dSJakub Kruzik /* set default KSPtype */ 585557a2f7dSJakub Kruzik if (!def->ksptype) { 586557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 587557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 588557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 589557a2f7dSJakub Kruzik } 590557a2f7dSJakub Kruzik } 591ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 592557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 593557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 594557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 595ae029463SJakub Kruzik /* inherit options */ 596557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 597557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 598557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 599557a2f7dSJakub Kruzik } else { /* the last level */ 600557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 601409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 602ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 603ac66f006SJakub Kruzik if (prefix) { 604ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 605ac66f006SJakub Kruzik } 60622b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 607ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 608409a357bSJakub Kruzik /* Reduction factor choice */ 609409a357bSJakub Kruzik red = def->reductionfact; 610409a357bSJakub Kruzik if (red < 0) { 611409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 612409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 613409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 614409a357bSJakub Kruzik if (match) red = commsize; 615409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 616409a357bSJakub Kruzik } 617409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 618ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 619409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 620ac66f006SJakub Kruzik if (innerksp) { 621409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 622409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 623ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 624409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 625409a357bSJakub Kruzik if (commsize == red) { 626ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 627409a357bSJakub Kruzik } else { 628ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 629409a357bSJakub Kruzik } 630409a357bSJakub Kruzik } 631557a2f7dSJakub Kruzik } 632557a2f7dSJakub Kruzik 633557a2f7dSJakub Kruzik if (innerksp) { 634409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 635409a357bSJakub Kruzik if (prefix) { 636409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 637409a357bSJakub Kruzik } 63822b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 639557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 640557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 641ac66f006SJakub Kruzik } 642409a357bSJakub Kruzik } 643557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 644557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 645409a357bSJakub Kruzik 64622b0793eSJakub Kruzik if (!def->pc) { 64722b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 64822b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 64922b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 65022b0793eSJakub Kruzik if (prefix) { 65122b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 65222b0793eSJakub Kruzik } 65322b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 65422b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 65522b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 65622b0793eSJakub Kruzik } 65722b0793eSJakub Kruzik 658f8bf229cSJakub Kruzik /* create work vecs */ 659f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 660f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 66137eeb815SJakub Kruzik PetscFunctionReturn(0); 66237eeb815SJakub Kruzik } 663b30d4775SJakub Kruzik 66437eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 66537eeb815SJakub Kruzik { 666b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 66737eeb815SJakub Kruzik PetscErrorCode ierr; 66837eeb815SJakub Kruzik 66937eeb815SJakub Kruzik PetscFunctionBegin; 670b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 671b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 672b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 673b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 674ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 675b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 676b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 67722b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 67837eeb815SJakub Kruzik PetscFunctionReturn(0); 67937eeb815SJakub Kruzik } 68037eeb815SJakub Kruzik 68137eeb815SJakub Kruzik /* 68237eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 68337eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 68437eeb815SJakub Kruzik 68537eeb815SJakub Kruzik Input Parameter: 68637eeb815SJakub Kruzik . pc - the preconditioner context 68737eeb815SJakub Kruzik 68837eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 68937eeb815SJakub Kruzik */ 69037eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 69137eeb815SJakub Kruzik { 69237eeb815SJakub Kruzik PetscErrorCode ierr; 69337eeb815SJakub Kruzik 69437eeb815SJakub Kruzik PetscFunctionBegin; 69537eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 69637eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 69737eeb815SJakub Kruzik PetscFunctionReturn(0); 69837eeb815SJakub Kruzik } 69937eeb815SJakub Kruzik 7008513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 7018513960bSJakub Kruzik { 7028513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 7038513960bSJakub Kruzik PetscBool iascii; 7048513960bSJakub Kruzik PetscErrorCode ierr; 7058513960bSJakub Kruzik 7068513960bSJakub Kruzik PetscFunctionBegin; 7078513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7088513960bSJakub Kruzik if (iascii) { 7098513960bSJakub Kruzik //if (cg->adaptiveconv) { 7108513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 7118513960bSJakub Kruzik //} 7128513960bSJakub Kruzik if (def->correct) { 7138513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 7148513960bSJakub Kruzik } 7158513960bSJakub Kruzik if (!def->nestedlvl) { 7168513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 7178513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 7188513960bSJakub Kruzik } 7198513960bSJakub Kruzik 72022b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 72122b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 72222b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 72322b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 72422b0793eSJakub Kruzik 7258513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 7268513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7278513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 7288513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7298513960bSJakub Kruzik } 7308513960bSJakub Kruzik PetscFunctionReturn(0); 7318513960bSJakub Kruzik } 7328513960bSJakub Kruzik 73337eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 73437eeb815SJakub Kruzik { 735880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 73637eeb815SJakub Kruzik PetscErrorCode ierr; 73737eeb815SJakub Kruzik 73837eeb815SJakub Kruzik PetscFunctionBegin; 73937eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 740880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 741880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 742880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 743880d05ceSJakub Kruzik //TODO add set function and fix manpages 744880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 745fcb31d99SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 746fcb31d99SJakub 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); 747880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 748880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 74937eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 75037eeb815SJakub Kruzik PetscFunctionReturn(0); 75137eeb815SJakub Kruzik } 75237eeb815SJakub Kruzik 75337eeb815SJakub Kruzik /*MC 754e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 755e361f8b5SJakub Kruzik or to a predefined value 75637eeb815SJakub Kruzik 75737eeb815SJakub Kruzik Options Database Key: 758e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 75937eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 76037eeb815SJakub Kruzik 76137eeb815SJakub Kruzik Level: beginner 76237eeb815SJakub Kruzik 76337eeb815SJakub Kruzik Notes: 764e361f8b5SJakub Kruzik todo 76537eeb815SJakub Kruzik 76637eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 767e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 76837eeb815SJakub Kruzik M*/ 76937eeb815SJakub Kruzik 77037eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 77137eeb815SJakub Kruzik { 77237eeb815SJakub Kruzik PC_Deflation *def; 77337eeb815SJakub Kruzik PetscErrorCode ierr; 77437eeb815SJakub Kruzik 77537eeb815SJakub Kruzik PetscFunctionBegin; 77637eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 77737eeb815SJakub Kruzik pc->data = (void*)def; 77837eeb815SJakub Kruzik 779e662bc50SJakub Kruzik def->init = PETSC_FALSE; 780e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 781fcb31d99SJakub Kruzik def->correctfact = 1.0; 782e662bc50SJakub Kruzik def->reductionfact = -1; 783e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 784e662bc50SJakub Kruzik def->spacesize = 1; 785e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 786e662bc50SJakub Kruzik def->nestedlvl = 0; 787e662bc50SJakub Kruzik def->maxnestedlvl = 0; 78837eeb815SJakub Kruzik 78937eeb815SJakub Kruzik /* 79037eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 79137eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 79237eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 79337eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 79437eeb815SJakub Kruzik not needed. 79537eeb815SJakub Kruzik */ 79637eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 79737eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 798b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 799b27c8092SJakub Kruzik pc->ops->postsolve = 0; 80037eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 80137eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 80237eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 80337eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 8048513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 80537eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 80637eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 80737eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 80837eeb815SJakub Kruzik 809e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 810*3cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 8115c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 812268c6673SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 81322b0793eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 814e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 8154a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 816f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 81737eeb815SJakub Kruzik PetscFunctionReturn(0); 81837eeb815SJakub Kruzik } 81937eeb815SJakub Kruzik 820