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 71*98d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 72*98d6e3deSJakub Kruzik { 73*98d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 74*98d6e3deSJakub Kruzik 75*98d6e3deSJakub Kruzik PetscFunctionBegin; 76*98d6e3deSJakub Kruzik if (current) def->nestedlvl = current; 77*98d6e3deSJakub Kruzik def->maxnestedlvl = max; 78*98d6e3deSJakub Kruzik PetscFunctionReturn(0); 79*98d6e3deSJakub Kruzik } 80*98d6e3deSJakub Kruzik 81*98d6e3deSJakub Kruzik /*@ 82*98d6e3deSJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 83*98d6e3deSJakub Kruzik 84*98d6e3deSJakub Kruzik Logically Collective on PC 85*98d6e3deSJakub Kruzik 86*98d6e3deSJakub Kruzik Input Parameters: 87*98d6e3deSJakub Kruzik + pc - the preconditioner context 88*98d6e3deSJakub Kruzik - max - maximum deflation level 89*98d6e3deSJakub Kruzik 90*98d6e3deSJakub Kruzik Level: intermediate 91*98d6e3deSJakub Kruzik 92*98d6e3deSJakub Kruzik .seealso: PCDEFLATION 93*98d6e3deSJakub Kruzik @*/ 94*98d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 95*98d6e3deSJakub Kruzik { 96*98d6e3deSJakub Kruzik PetscErrorCode ierr; 97*98d6e3deSJakub Kruzik 98*98d6e3deSJakub Kruzik PetscFunctionBegin; 99*98d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 100*98d6e3deSJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 101*98d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 102*98d6e3deSJakub Kruzik PetscFunctionReturn(0); 103*98d6e3deSJakub Kruzik } 104*98d6e3deSJakub Kruzik 10539417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 10639417d7eSJakub Kruzik { 10739417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 10839417d7eSJakub Kruzik 10939417d7eSJakub Kruzik PetscFunctionBegin; 11039417d7eSJakub Kruzik if (type) def->spacetype = type; 11139417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 11239417d7eSJakub Kruzik PetscFunctionReturn(0); 11339417d7eSJakub Kruzik } 11439417d7eSJakub Kruzik 11539417d7eSJakub Kruzik /*@ 11639417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 11739417d7eSJakub Kruzik 11839417d7eSJakub Kruzik Logically Collective on PC 11939417d7eSJakub Kruzik 12039417d7eSJakub Kruzik Input Parameters: 12139417d7eSJakub Kruzik + pc - the preconditioner context 12239417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 12339417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 12439417d7eSJakub Kruzik 12539417d7eSJakub Kruzik Notes: 12639417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 12739417d7eSJakub Kruzik The deflation space is computed in PCSetUP(). 12839417d7eSJakub Kruzik 12939417d7eSJakub Kruzik Level: intermediate 13039417d7eSJakub Kruzik 13139417d7eSJakub Kruzik .seealso: PCDEFLATION 13239417d7eSJakub Kruzik @*/ 13339417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 13439417d7eSJakub Kruzik { 13539417d7eSJakub Kruzik PetscErrorCode ierr; 13639417d7eSJakub Kruzik 13739417d7eSJakub Kruzik PetscFunctionBegin; 13839417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13939417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 14039417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 14139417d7eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 14239417d7eSJakub Kruzik PetscFunctionReturn(0); 14339417d7eSJakub Kruzik } 1448513960bSJakub Kruzik 145e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 146e662bc50SJakub Kruzik { 147e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 148e662bc50SJakub Kruzik PetscErrorCode ierr; 149e662bc50SJakub Kruzik 150e662bc50SJakub Kruzik PetscFunctionBegin; 151e662bc50SJakub Kruzik if (transpose) { 152e662bc50SJakub Kruzik def->Wt = W; 153e662bc50SJakub Kruzik def->W = NULL; 154e662bc50SJakub Kruzik } else { 155e662bc50SJakub Kruzik def->W = W; 156e662bc50SJakub Kruzik } 157e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 158e662bc50SJakub Kruzik PetscFunctionReturn(0); 159e662bc50SJakub Kruzik } 160e662bc50SJakub Kruzik 161e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 162e662bc50SJakub Kruzik /*@ 163e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 164e662bc50SJakub Kruzik 165e662bc50SJakub Kruzik Logically Collective on PC 166e662bc50SJakub Kruzik 167e662bc50SJakub Kruzik Input Parameters: 168e662bc50SJakub Kruzik + pc - the preconditioner context 169e662bc50SJakub Kruzik . W - deflation matrix 170e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 171e662bc50SJakub Kruzik 172e662bc50SJakub Kruzik Level: intermediate 173e662bc50SJakub Kruzik 174e662bc50SJakub Kruzik .seealso: PCDEFLATION 175e662bc50SJakub Kruzik @*/ 176e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 177e662bc50SJakub Kruzik { 178e662bc50SJakub Kruzik PetscErrorCode ierr; 179e662bc50SJakub Kruzik 180e662bc50SJakub Kruzik PetscFunctionBegin; 181e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 183e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 184e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 185e662bc50SJakub Kruzik PetscFunctionReturn(0); 186e662bc50SJakub Kruzik } 187e662bc50SJakub Kruzik 1883cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 1893cf3a049SJakub Kruzik { 1903cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1913cf3a049SJakub Kruzik PetscErrorCode ierr; 1923cf3a049SJakub Kruzik 1933cf3a049SJakub Kruzik PetscFunctionBegin; 1943cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 1953cf3a049SJakub Kruzik def->WtA = mat; 1963cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1973cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 1983cf3a049SJakub Kruzik PetscFunctionReturn(0); 1993cf3a049SJakub Kruzik } 2003cf3a049SJakub Kruzik 2013cf3a049SJakub Kruzik /*@ 2023cf3a049SJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 2033cf3a049SJakub Kruzik 2043cf3a049SJakub Kruzik Not Collective 2053cf3a049SJakub Kruzik 2063cf3a049SJakub Kruzik Input Parameters: 2073cf3a049SJakub Kruzik + pc - preconditioner context 2083cf3a049SJakub Kruzik - mat - projection null space matrix 2093cf3a049SJakub Kruzik 2103cf3a049SJakub Kruzik Level: developer 2113cf3a049SJakub Kruzik 2123cf3a049SJakub Kruzik .seealso: PCDEFLATION 2133cf3a049SJakub Kruzik @*/ 2143cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 2153cf3a049SJakub Kruzik { 2163cf3a049SJakub Kruzik PetscErrorCode ierr; 2173cf3a049SJakub Kruzik 2183cf3a049SJakub Kruzik PetscFunctionBegin; 2193cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2203cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 2213cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 2223cf3a049SJakub Kruzik PetscFunctionReturn(0); 2233cf3a049SJakub Kruzik } 2243cf3a049SJakub Kruzik 225e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 226e924b002SJakub Kruzik { 227e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 228e924b002SJakub Kruzik PetscErrorCode ierr; 229e924b002SJakub Kruzik 230e924b002SJakub Kruzik PetscFunctionBegin; 231e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 232e924b002SJakub Kruzik def->WtAW = mat; 233e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 234e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 235e924b002SJakub Kruzik PetscFunctionReturn(0); 236e924b002SJakub Kruzik } 237e924b002SJakub Kruzik 238e924b002SJakub Kruzik /*@ 239e924b002SJakub Kruzik PCDeflationSetCoarseMat - Set coarse problem Mat. 240e924b002SJakub Kruzik 241e924b002SJakub Kruzik Not Collective 242e924b002SJakub Kruzik 243e924b002SJakub Kruzik Input Parameters: 244e924b002SJakub Kruzik + pc - preconditioner context 245e924b002SJakub Kruzik - mat - coarse problem mat 246e924b002SJakub Kruzik 247e924b002SJakub Kruzik Level: developer 248e924b002SJakub Kruzik 249e924b002SJakub Kruzik .seealso: PCDEFLATION 250e924b002SJakub Kruzik @*/ 251e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 252e924b002SJakub Kruzik { 253e924b002SJakub Kruzik PetscErrorCode ierr; 254e924b002SJakub Kruzik 255e924b002SJakub Kruzik PetscFunctionBegin; 256e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 257e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 258e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 259e924b002SJakub Kruzik PetscFunctionReturn(0); 260e924b002SJakub Kruzik } 261e924b002SJakub Kruzik 262*98d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 2635c0c31c5SJakub Kruzik { 2645c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2655c0c31c5SJakub Kruzik 2665c0c31c5SJakub Kruzik PetscFunctionBegin; 267*98d6e3deSJakub Kruzik *ksp = def->WtAWinv; 2685c0c31c5SJakub Kruzik PetscFunctionReturn(0); 2695c0c31c5SJakub Kruzik } 2705c0c31c5SJakub Kruzik 2715c0c31c5SJakub Kruzik /*@ 272*98d6e3deSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 2735c0c31c5SJakub Kruzik 274*98d6e3deSJakub Kruzik Not Collective 2755c0c31c5SJakub Kruzik 2765c0c31c5SJakub Kruzik Input Parameters: 277*98d6e3deSJakub Kruzik . pc - preconditioner context 2785c0c31c5SJakub Kruzik 279*98d6e3deSJakub Kruzik Output Parameter: 280*98d6e3deSJakub Kruzik . ksp - coarse problem KSP context 2815c0c31c5SJakub Kruzik 282*98d6e3deSJakub Kruzik Level: developer 283*98d6e3deSJakub Kruzik 284*98d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 2855c0c31c5SJakub Kruzik @*/ 286*98d6e3deSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 2875c0c31c5SJakub Kruzik { 2885c0c31c5SJakub Kruzik PetscErrorCode ierr; 2895c0c31c5SJakub Kruzik 2905c0c31c5SJakub Kruzik PetscFunctionBegin; 29122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 292*98d6e3deSJakub Kruzik PetscValidPointer(ksp,2); 293*98d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 294*98d6e3deSJakub Kruzik PetscFunctionReturn(0); 295*98d6e3deSJakub Kruzik } 296*98d6e3deSJakub Kruzik 297*98d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 298*98d6e3deSJakub Kruzik { 299*98d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 300*98d6e3deSJakub Kruzik PetscErrorCode ierr; 301*98d6e3deSJakub Kruzik 302*98d6e3deSJakub Kruzik PetscFunctionBegin; 303*98d6e3deSJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 304*98d6e3deSJakub Kruzik def->WtAWinv = ksp; 305*98d6e3deSJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 306*98d6e3deSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 307*98d6e3deSJakub Kruzik PetscFunctionReturn(0); 308*98d6e3deSJakub Kruzik } 309*98d6e3deSJakub Kruzik 310*98d6e3deSJakub Kruzik /*@ 311*98d6e3deSJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 312*98d6e3deSJakub Kruzik 313*98d6e3deSJakub Kruzik Collective on PC 314*98d6e3deSJakub Kruzik 315*98d6e3deSJakub Kruzik Input Parameters: 316*98d6e3deSJakub Kruzik + pc - preconditioner context 317*98d6e3deSJakub Kruzik - ksp - coarse problem KSP context 318*98d6e3deSJakub Kruzik 319*98d6e3deSJakub Kruzik Level: developer 320*98d6e3deSJakub Kruzik 321*98d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 322*98d6e3deSJakub Kruzik @*/ 323*98d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 324*98d6e3deSJakub Kruzik { 325*98d6e3deSJakub Kruzik PetscErrorCode ierr; 326*98d6e3deSJakub Kruzik 327*98d6e3deSJakub Kruzik PetscFunctionBegin; 328*98d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 329*98d6e3deSJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 330*98d6e3deSJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 331*98d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 3325c0c31c5SJakub Kruzik PetscFunctionReturn(0); 3335c0c31c5SJakub Kruzik } 334e662bc50SJakub Kruzik 335268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 336268c6673SJakub Kruzik { 337268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 338268c6673SJakub Kruzik 339268c6673SJakub Kruzik PetscFunctionBegin; 340268c6673SJakub Kruzik *apc = def->pc; 341268c6673SJakub Kruzik PetscFunctionReturn(0); 342268c6673SJakub Kruzik } 343268c6673SJakub Kruzik 344268c6673SJakub Kruzik /*@ 345268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 346268c6673SJakub Kruzik 347268c6673SJakub Kruzik Not Collective 348268c6673SJakub Kruzik 349268c6673SJakub Kruzik Input Parameters: 350268c6673SJakub Kruzik . pc - the preconditioner context 351268c6673SJakub Kruzik 352268c6673SJakub Kruzik Output Parameter: 353268c6673SJakub Kruzik . apc - additional preconditioner 354268c6673SJakub Kruzik 355268c6673SJakub Kruzik Level: advanced 356268c6673SJakub Kruzik 357268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 358268c6673SJakub Kruzik @*/ 359268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 360268c6673SJakub Kruzik { 361268c6673SJakub Kruzik PetscErrorCode ierr; 362268c6673SJakub Kruzik 363268c6673SJakub Kruzik PetscFunctionBegin; 364268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 365268c6673SJakub Kruzik PetscValidPointer(pc,2); 366268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 367268c6673SJakub Kruzik PetscFunctionReturn(0); 368268c6673SJakub Kruzik } 369268c6673SJakub Kruzik 37022b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 37122b0793eSJakub Kruzik { 37222b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 37322b0793eSJakub Kruzik PetscErrorCode ierr; 37422b0793eSJakub Kruzik 37522b0793eSJakub Kruzik PetscFunctionBegin; 37622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 37722b0793eSJakub Kruzik def->pc = apc; 37822b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 37922b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 38022b0793eSJakub Kruzik PetscFunctionReturn(0); 38122b0793eSJakub Kruzik } 38222b0793eSJakub Kruzik 38322b0793eSJakub Kruzik /*@ 38422b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 38522b0793eSJakub Kruzik 386268c6673SJakub Kruzik Collective on PC 38722b0793eSJakub Kruzik 38822b0793eSJakub Kruzik Input Parameters: 38922b0793eSJakub Kruzik + pc - the preconditioner context 39022b0793eSJakub Kruzik - apc - additional preconditioner 39122b0793eSJakub Kruzik 392268c6673SJakub Kruzik Level: developer 39322b0793eSJakub Kruzik 394268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 39522b0793eSJakub Kruzik @*/ 39622b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 39722b0793eSJakub Kruzik { 39822b0793eSJakub Kruzik PetscErrorCode ierr; 39922b0793eSJakub Kruzik 40022b0793eSJakub Kruzik PetscFunctionBegin; 40122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40222b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 40322b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 40422b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 40522b0793eSJakub Kruzik PetscFunctionReturn(0); 40622b0793eSJakub Kruzik } 40722b0793eSJakub Kruzik 40837eeb815SJakub Kruzik /* 409b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 410b27c8092SJakub Kruzik */ 411b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 412b27c8092SJakub Kruzik { 413b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 414b27c8092SJakub Kruzik Mat A; 415b27c8092SJakub Kruzik Vec r,w1,w2; 416b27c8092SJakub Kruzik PetscBool nonzero; 417b27c8092SJakub Kruzik PetscErrorCode ierr; 41837eeb815SJakub Kruzik 419b27c8092SJakub Kruzik PetscFunctionBegin; 420b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 421b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 422b27c8092SJakub Kruzik r = def->work; 423b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 42437eeb815SJakub Kruzik 425b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 426b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 427b27c8092SJakub Kruzik if (nonzero) { 428b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 429b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 430b27c8092SJakub Kruzik } else { 431b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 432b27c8092SJakub Kruzik } 433b27c8092SJakub Kruzik 434a3177931SJakub Kruzik if (def->Wt) { 435a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 436a3177931SJakub Kruzik } else { 437a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 438a3177931SJakub Kruzik } 439b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 440b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 441b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 442b27c8092SJakub Kruzik PetscFunctionReturn(0); 443b27c8092SJakub Kruzik } 44437eeb815SJakub Kruzik 445f8bf229cSJakub Kruzik /* 446f8bf229cSJakub Kruzik if (def->correct) { 447ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 448f8bf229cSJakub Kruzik } else { 449ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 450f8bf229cSJakub Kruzik } 45137eeb815SJakub Kruzik */ 452f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 453f8bf229cSJakub Kruzik { 454f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 455f8bf229cSJakub Kruzik Mat A; 456f8bf229cSJakub Kruzik Vec u,w1,w2; 457f8bf229cSJakub Kruzik PetscErrorCode ierr; 458f8bf229cSJakub Kruzik 459f8bf229cSJakub Kruzik PetscFunctionBegin; 460f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 461f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 462f8bf229cSJakub Kruzik u = def->work; 463f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 464f8bf229cSJakub Kruzik 465ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 46622b0793eSJakub Kruzik if (!def->init) { 467ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 468f8bf229cSJakub Kruzik if (def->correct) { 469ae029463SJakub Kruzik if (def->Wt) { 470ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 471ae029463SJakub Kruzik } else { 472a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 473f8bf229cSJakub Kruzik } 474ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 475f8bf229cSJakub Kruzik } 476f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 477f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 47822b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 47922b0793eSJakub Kruzik } 480f8bf229cSJakub Kruzik PetscFunctionReturn(0); 481f8bf229cSJakub Kruzik } 482f8bf229cSJakub Kruzik 48337eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 48437eeb815SJakub Kruzik { 485409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 486409a357bSJakub Kruzik KSP innerksp; 487409a357bSJakub Kruzik PC pcinner; 488409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 489409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 490409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 491409a357bSJakub Kruzik MatCompositeType ctype; 492409a357bSJakub Kruzik MPI_Comm comm; 493409a357bSJakub Kruzik const char *prefix; 49437eeb815SJakub Kruzik PetscErrorCode ierr; 49537eeb815SJakub Kruzik 49637eeb815SJakub Kruzik PetscFunctionBegin; 497409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 498409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 499f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 500f8bf229cSJakub Kruzik 501f8bf229cSJakub Kruzik /* compute a deflation space */ 502409a357bSJakub Kruzik if (def->W || def->Wt) { 503409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 504409a357bSJakub Kruzik } else { 505e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 50637eeb815SJakub Kruzik } 50737eeb815SJakub Kruzik 508409a357bSJakub Kruzik /* nested deflation */ 509409a357bSJakub Kruzik if (def->W) { 510409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 511409a357bSJakub Kruzik if (match) { 512409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 513409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 51437eeb815SJakub Kruzik } 515409a357bSJakub Kruzik } else { 516a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 517409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 518409a357bSJakub Kruzik if (match) { 519409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 520409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 521409a357bSJakub Kruzik } 522409a357bSJakub Kruzik transp = PETSC_TRUE; 523409a357bSJakub Kruzik } 524409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 525409a357bSJakub Kruzik if (!transp) { 52628da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 52728da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 528409a357bSJakub Kruzik for (i=0; i<size; i++) { 529409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 530409a357bSJakub Kruzik } 531409a357bSJakub Kruzik size -= 1; 532409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 533409a357bSJakub Kruzik def->W = mats[size]; 534409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 535409a357bSJakub Kruzik if (size > 1) { 536409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 537409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 538409a357bSJakub Kruzik } else { 539409a357bSJakub Kruzik nextDef = mats[0]; 540409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 541409a357bSJakub Kruzik } 54228da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 543409a357bSJakub Kruzik } else { 544409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 545409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 546409a357bSJakub Kruzik } 54728da0170SJakub Kruzik } else { 54828da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 54928da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 55028da0170SJakub Kruzik for (i=0; i<size; i++) { 55128da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 55228da0170SJakub Kruzik } 55328da0170SJakub Kruzik size -= 1; 55428da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 55528da0170SJakub Kruzik def->Wt = mats[0]; 55628da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 55728da0170SJakub Kruzik if (size > 1) { 55828da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 55928da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 56028da0170SJakub Kruzik } else { 56128da0170SJakub Kruzik nextDef = mats[1]; 56228da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 563409a357bSJakub Kruzik } 564409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 56528da0170SJakub Kruzik } else { 56628da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 56728da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 56828da0170SJakub Kruzik } 56928da0170SJakub Kruzik } 57028da0170SJakub Kruzik } 57128da0170SJakub Kruzik 57228da0170SJakub Kruzik if (transp) { 57328da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 574a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 575409a357bSJakub Kruzik } 576409a357bSJakub Kruzik 57722b0793eSJakub Kruzik 57822b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 57922b0793eSJakub Kruzik 580ae029463SJakub Kruzik /* assemble WtA */ 581ae029463SJakub Kruzik if (!def->WtA) { 582ae029463SJakub Kruzik if (def->Wt) { 583ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 584ae029463SJakub Kruzik } else { 585a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 586a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 587a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 588a3177931SJakub Kruzik #else 589ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 590a3177931SJakub Kruzik #endif 591ae029463SJakub Kruzik } 592ae029463SJakub Kruzik } 593409a357bSJakub Kruzik /* setup coarse problem */ 594409a357bSJakub Kruzik if (!def->WtAWinv) { 59528da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 596409a357bSJakub Kruzik if (!def->WtAW) { 597ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 598409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 599409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 600409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 601409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 602ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 603409a357bSJakub Kruzik PetscReal *norms; 604409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 605409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 606409a357bSJakub Kruzik for (i=0; i<m; i++) { 607409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 608409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 609409a357bSJakub Kruzik } 610409a357bSJakub Kruzik } 611409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 612409a357bSJakub Kruzik #endif 613409a357bSJakub Kruzik } else { 614409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 615409a357bSJakub Kruzik } 616409a357bSJakub Kruzik /* TODO use MATINV */ 617409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 618409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 619409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 620557a2f7dSJakub Kruzik /* Setup KSP and PC */ 621557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 622557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 623557a2f7dSJakub Kruzik /* set default KSPtype */ 624557a2f7dSJakub Kruzik if (!def->ksptype) { 625557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 626557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 627557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 628557a2f7dSJakub Kruzik } 629557a2f7dSJakub Kruzik } 630ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 631557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 632557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 633557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 634ae029463SJakub Kruzik /* inherit options */ 635557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 636557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 637557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 638557a2f7dSJakub Kruzik } else { /* the last level */ 639557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 640409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 641ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 642ac66f006SJakub Kruzik if (prefix) { 643ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 644ac66f006SJakub Kruzik } 64522b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 646ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 647409a357bSJakub Kruzik /* Reduction factor choice */ 648409a357bSJakub Kruzik red = def->reductionfact; 649409a357bSJakub Kruzik if (red < 0) { 650409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 651409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 652409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 653409a357bSJakub Kruzik if (match) red = commsize; 654409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 655409a357bSJakub Kruzik } 656409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 657ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 658409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 659ac66f006SJakub Kruzik if (innerksp) { 660409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 661409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 662ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 663409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 664409a357bSJakub Kruzik if (commsize == red) { 665ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 666409a357bSJakub Kruzik } else { 667ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 668409a357bSJakub Kruzik } 669409a357bSJakub Kruzik } 670557a2f7dSJakub Kruzik } 671557a2f7dSJakub Kruzik 672557a2f7dSJakub Kruzik if (innerksp) { 673409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 674409a357bSJakub Kruzik if (prefix) { 675409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 676409a357bSJakub Kruzik } 67722b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 678557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 679557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 680ac66f006SJakub Kruzik } 681409a357bSJakub Kruzik } 682557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 683557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 684409a357bSJakub Kruzik 68522b0793eSJakub Kruzik if (!def->pc) { 68622b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 68722b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 68822b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 68922b0793eSJakub Kruzik if (prefix) { 69022b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 69122b0793eSJakub Kruzik } 69222b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 69322b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 69422b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 69522b0793eSJakub Kruzik } 69622b0793eSJakub Kruzik 697f8bf229cSJakub Kruzik /* create work vecs */ 698f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 699f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 70037eeb815SJakub Kruzik PetscFunctionReturn(0); 70137eeb815SJakub Kruzik } 702b30d4775SJakub Kruzik 70337eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 70437eeb815SJakub Kruzik { 705b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 70637eeb815SJakub Kruzik PetscErrorCode ierr; 70737eeb815SJakub Kruzik 70837eeb815SJakub Kruzik PetscFunctionBegin; 709b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 710b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 711b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 712b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 713ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 714b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 715b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 71622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 71737eeb815SJakub Kruzik PetscFunctionReturn(0); 71837eeb815SJakub Kruzik } 71937eeb815SJakub Kruzik 72037eeb815SJakub Kruzik /* 72137eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 72237eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 72337eeb815SJakub Kruzik 72437eeb815SJakub Kruzik Input Parameter: 72537eeb815SJakub Kruzik . pc - the preconditioner context 72637eeb815SJakub Kruzik 72737eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 72837eeb815SJakub Kruzik */ 72937eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 73037eeb815SJakub Kruzik { 73137eeb815SJakub Kruzik PetscErrorCode ierr; 73237eeb815SJakub Kruzik 73337eeb815SJakub Kruzik PetscFunctionBegin; 73437eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 73537eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 73637eeb815SJakub Kruzik PetscFunctionReturn(0); 73737eeb815SJakub Kruzik } 73837eeb815SJakub Kruzik 7398513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 7408513960bSJakub Kruzik { 7418513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 7428513960bSJakub Kruzik PetscBool iascii; 7438513960bSJakub Kruzik PetscErrorCode ierr; 7448513960bSJakub Kruzik 7458513960bSJakub Kruzik PetscFunctionBegin; 7468513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7478513960bSJakub Kruzik if (iascii) { 7488513960bSJakub Kruzik //if (cg->adaptiveconv) { 7498513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 7508513960bSJakub Kruzik //} 7518513960bSJakub Kruzik if (def->correct) { 7528513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 7538513960bSJakub Kruzik } 7548513960bSJakub Kruzik if (!def->nestedlvl) { 7558513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 7568513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 7578513960bSJakub Kruzik } 7588513960bSJakub Kruzik 75922b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 76022b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 76122b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 76222b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 76322b0793eSJakub Kruzik 7648513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 7658513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7668513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 7678513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7688513960bSJakub Kruzik } 7698513960bSJakub Kruzik PetscFunctionReturn(0); 7708513960bSJakub Kruzik } 7718513960bSJakub Kruzik 77237eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 77337eeb815SJakub Kruzik { 774880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 77537eeb815SJakub Kruzik PetscErrorCode ierr; 77637eeb815SJakub Kruzik 77737eeb815SJakub Kruzik PetscFunctionBegin; 77837eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 779880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 780880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 781880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 782880d05ceSJakub Kruzik //TODO add set function and fix manpages 783880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 784fcb31d99SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 785fcb31d99SJakub 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); 786880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 787880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 78837eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 78937eeb815SJakub Kruzik PetscFunctionReturn(0); 79037eeb815SJakub Kruzik } 79137eeb815SJakub Kruzik 79237eeb815SJakub Kruzik /*MC 793e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 794e361f8b5SJakub Kruzik or to a predefined value 79537eeb815SJakub Kruzik 79637eeb815SJakub Kruzik Options Database Key: 797e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 79837eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 79937eeb815SJakub Kruzik 80037eeb815SJakub Kruzik Level: beginner 80137eeb815SJakub Kruzik 80237eeb815SJakub Kruzik Notes: 803e361f8b5SJakub Kruzik todo 80437eeb815SJakub Kruzik 80537eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 806e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 80737eeb815SJakub Kruzik M*/ 80837eeb815SJakub Kruzik 80937eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 81037eeb815SJakub Kruzik { 81137eeb815SJakub Kruzik PC_Deflation *def; 81237eeb815SJakub Kruzik PetscErrorCode ierr; 81337eeb815SJakub Kruzik 81437eeb815SJakub Kruzik PetscFunctionBegin; 81537eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 81637eeb815SJakub Kruzik pc->data = (void*)def; 81737eeb815SJakub Kruzik 818e662bc50SJakub Kruzik def->init = PETSC_FALSE; 819e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 820fcb31d99SJakub Kruzik def->correctfact = 1.0; 821e662bc50SJakub Kruzik def->reductionfact = -1; 822e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 823e662bc50SJakub Kruzik def->spacesize = 1; 824e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 825e662bc50SJakub Kruzik def->nestedlvl = 0; 826e662bc50SJakub Kruzik def->maxnestedlvl = 0; 82737eeb815SJakub Kruzik 82837eeb815SJakub Kruzik /* 82937eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 83037eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 83137eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 83237eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 83337eeb815SJakub Kruzik not needed. 83437eeb815SJakub Kruzik */ 83537eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 83637eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 837b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 838b27c8092SJakub Kruzik pc->ops->postsolve = 0; 83937eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 84037eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 84137eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 84237eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 8438513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 84437eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 84537eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 84637eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 84737eeb815SJakub Kruzik 848*98d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 84939417d7eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 850e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 8513cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 852e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 8534a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 854f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 855*98d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 856*98d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 85737eeb815SJakub Kruzik PetscFunctionReturn(0); 85837eeb815SJakub Kruzik } 85937eeb815SJakub Kruzik 860