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*a122ebaeSJakub Kruzik static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc,PetscBool flg) 72*a122ebaeSJakub Kruzik { 73*a122ebaeSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 74*a122ebaeSJakub Kruzik 75*a122ebaeSJakub Kruzik PetscFunctionBegin; 76*a122ebaeSJakub Kruzik def->init = flg; 77*a122ebaeSJakub Kruzik PetscFunctionReturn(0); 78*a122ebaeSJakub Kruzik } 79*a122ebaeSJakub Kruzik 80*a122ebaeSJakub Kruzik /*@ 81*a122ebaeSJakub Kruzik PCDeflationSetInitOnly - Do only initialization step. 82*a122ebaeSJakub Kruzik Sets initial guess to the solution on the deflation space but do not apply deflation preconditioner. 83*a122ebaeSJakub Kruzik The additional preconditioner is still applied. 84*a122ebaeSJakub Kruzik 85*a122ebaeSJakub Kruzik Logically Collective on PC 86*a122ebaeSJakub Kruzik 87*a122ebaeSJakub Kruzik Input Parameters: 88*a122ebaeSJakub Kruzik + pc - the preconditioner context 89*a122ebaeSJakub Kruzik - flg - default PETSC_FALSE 90*a122ebaeSJakub Kruzik 91*a122ebaeSJakub Kruzik Level: intermediate 92*a122ebaeSJakub Kruzik 93*a122ebaeSJakub Kruzik .seealso: PCDEFLATION 94*a122ebaeSJakub Kruzik @*/ 95*a122ebaeSJakub Kruzik PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg) 96*a122ebaeSJakub Kruzik { 97*a122ebaeSJakub Kruzik PetscErrorCode ierr; 98*a122ebaeSJakub Kruzik 99*a122ebaeSJakub Kruzik PetscFunctionBegin; 100*a122ebaeSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 101*a122ebaeSJakub Kruzik PetscValidLogicalCollectiveBool(pc,flg,2); 102*a122ebaeSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 103*a122ebaeSJakub Kruzik PetscFunctionReturn(0); 104*a122ebaeSJakub Kruzik } 105*a122ebaeSJakub Kruzik 106*a122ebaeSJakub Kruzik 10798d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 10898d6e3deSJakub Kruzik { 10998d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 11098d6e3deSJakub Kruzik 11198d6e3deSJakub Kruzik PetscFunctionBegin; 11298d6e3deSJakub Kruzik if (current) def->nestedlvl = current; 11398d6e3deSJakub Kruzik def->maxnestedlvl = max; 11498d6e3deSJakub Kruzik PetscFunctionReturn(0); 11598d6e3deSJakub Kruzik } 11698d6e3deSJakub Kruzik 11798d6e3deSJakub Kruzik /*@ 11898d6e3deSJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 11998d6e3deSJakub Kruzik 12098d6e3deSJakub Kruzik Logically Collective on PC 12198d6e3deSJakub Kruzik 12298d6e3deSJakub Kruzik Input Parameters: 12398d6e3deSJakub Kruzik + pc - the preconditioner context 12498d6e3deSJakub Kruzik - max - maximum deflation level 12598d6e3deSJakub Kruzik 12698d6e3deSJakub Kruzik Level: intermediate 12798d6e3deSJakub Kruzik 12898d6e3deSJakub Kruzik .seealso: PCDEFLATION 12998d6e3deSJakub Kruzik @*/ 13098d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 13198d6e3deSJakub Kruzik { 13298d6e3deSJakub Kruzik PetscErrorCode ierr; 13398d6e3deSJakub Kruzik 13498d6e3deSJakub Kruzik PetscFunctionBegin; 13598d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13698d6e3deSJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 13798d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 13898d6e3deSJakub Kruzik PetscFunctionReturn(0); 13998d6e3deSJakub Kruzik } 14098d6e3deSJakub Kruzik 141859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red) 142859a873cSJakub Kruzik { 143859a873cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 144859a873cSJakub Kruzik 145859a873cSJakub Kruzik PetscFunctionBegin; 146859a873cSJakub Kruzik def->reductionfact = red; 147859a873cSJakub Kruzik PetscFunctionReturn(0); 148859a873cSJakub Kruzik } 149859a873cSJakub Kruzik 150859a873cSJakub Kruzik /*@ 151859a873cSJakub Kruzik PCDeflationSetReductionFactor - Set reduction factor for PCTELESCOPE coarse problem solver. 152859a873cSJakub Kruzik 153859a873cSJakub Kruzik Logically Collective on PC 154859a873cSJakub Kruzik 155859a873cSJakub Kruzik Input Parameters: 156859a873cSJakub Kruzik + pc - the preconditioner context 157859a873cSJakub Kruzik - red - reduction factor (or PETSC_DETERMINE) 158859a873cSJakub Kruzik 159859a873cSJakub Kruzik Level: intermediate 160859a873cSJakub Kruzik 161859a873cSJakub Kruzik .seealso: PCDEFLATION 162859a873cSJakub Kruzik @*/ 163859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red) 164859a873cSJakub Kruzik { 165859a873cSJakub Kruzik PetscErrorCode ierr; 166859a873cSJakub Kruzik 167859a873cSJakub Kruzik PetscFunctionBegin; 168859a873cSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 169859a873cSJakub Kruzik PetscValidLogicalCollectiveInt(pc,red,2); 170859a873cSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr); 171859a873cSJakub Kruzik PetscFunctionReturn(0); 172859a873cSJakub Kruzik } 173859a873cSJakub Kruzik 1748a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact) 1758a71cb68SJakub Kruzik { 1768a71cb68SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1778a71cb68SJakub Kruzik 1788a71cb68SJakub Kruzik PetscFunctionBegin; 1798a71cb68SJakub Kruzik /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 1808a71cb68SJakub Kruzik def->correct = PETSC_TRUE; 1818a71cb68SJakub Kruzik def->correctfact = fact; 1828a71cb68SJakub Kruzik if (fact == PETSC_DEFAULT) { 1838a71cb68SJakub Kruzik def->correctfact = 1.0; 1848a71cb68SJakub Kruzik } else if (def->correct == 0.0) { 1858a71cb68SJakub Kruzik def->correct = PETSC_FALSE; 1868a71cb68SJakub Kruzik } 1878a71cb68SJakub Kruzik PetscFunctionReturn(0); 1888a71cb68SJakub Kruzik } 1898a71cb68SJakub Kruzik 1908a71cb68SJakub Kruzik /*@ 1918a71cb68SJakub Kruzik PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 1928a71cb68SJakub Kruzik The Preconditioner becomes P*M^{-1} + factor*Q. 1938a71cb68SJakub Kruzik 1948a71cb68SJakub Kruzik Logically Collective on PC 1958a71cb68SJakub Kruzik 1968a71cb68SJakub Kruzik Input Parameters: 1978a71cb68SJakub Kruzik + pc - the preconditioner context 1988a71cb68SJakub Kruzik - fact - correction factor (set 0.0 to disable, PETSC_DEFAULT = 1.0) 1998a71cb68SJakub Kruzik 2008a71cb68SJakub Kruzik Level: intermediate 2018a71cb68SJakub Kruzik 2028a71cb68SJakub Kruzik .seealso: PCDEFLATION 2038a71cb68SJakub Kruzik @*/ 2048a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact) 2058a71cb68SJakub Kruzik { 2068a71cb68SJakub Kruzik PetscErrorCode ierr; 2078a71cb68SJakub Kruzik 2088a71cb68SJakub Kruzik PetscFunctionBegin; 2098a71cb68SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2108a71cb68SJakub Kruzik PetscValidLogicalCollectiveScalar(pc,fact,2); 2118a71cb68SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr); 2128a71cb68SJakub Kruzik PetscFunctionReturn(0); 2138a71cb68SJakub Kruzik } 2148a71cb68SJakub Kruzik 21539417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 21639417d7eSJakub Kruzik { 21739417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 21839417d7eSJakub Kruzik 21939417d7eSJakub Kruzik PetscFunctionBegin; 22039417d7eSJakub Kruzik if (type) def->spacetype = type; 22139417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 22239417d7eSJakub Kruzik PetscFunctionReturn(0); 22339417d7eSJakub Kruzik } 22439417d7eSJakub Kruzik 22539417d7eSJakub Kruzik /*@ 22639417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 22739417d7eSJakub Kruzik 22839417d7eSJakub Kruzik Logically Collective on PC 22939417d7eSJakub Kruzik 23039417d7eSJakub Kruzik Input Parameters: 23139417d7eSJakub Kruzik + pc - the preconditioner context 23239417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 23339417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 23439417d7eSJakub Kruzik 23539417d7eSJakub Kruzik Notes: 23639417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 23739417d7eSJakub Kruzik The deflation space is computed in PCSetUP(). 23839417d7eSJakub Kruzik 23939417d7eSJakub Kruzik Level: intermediate 24039417d7eSJakub Kruzik 24139417d7eSJakub Kruzik .seealso: PCDEFLATION 24239417d7eSJakub Kruzik @*/ 24339417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 24439417d7eSJakub Kruzik { 24539417d7eSJakub Kruzik PetscErrorCode ierr; 24639417d7eSJakub Kruzik 24739417d7eSJakub Kruzik PetscFunctionBegin; 24839417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 24939417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 25039417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 25139417d7eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 25239417d7eSJakub Kruzik PetscFunctionReturn(0); 25339417d7eSJakub Kruzik } 2548513960bSJakub Kruzik 255e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 256e662bc50SJakub Kruzik { 257e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 258e662bc50SJakub Kruzik PetscErrorCode ierr; 259e662bc50SJakub Kruzik 260e662bc50SJakub Kruzik PetscFunctionBegin; 261e662bc50SJakub Kruzik if (transpose) { 262e662bc50SJakub Kruzik def->Wt = W; 263e662bc50SJakub Kruzik def->W = NULL; 264e662bc50SJakub Kruzik } else { 265e662bc50SJakub Kruzik def->W = W; 266e662bc50SJakub Kruzik } 267e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 268e662bc50SJakub Kruzik PetscFunctionReturn(0); 269e662bc50SJakub Kruzik } 270e662bc50SJakub Kruzik 271e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 272e662bc50SJakub Kruzik /*@ 273e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 274e662bc50SJakub Kruzik 275e662bc50SJakub Kruzik Logically Collective on PC 276e662bc50SJakub Kruzik 277e662bc50SJakub Kruzik Input Parameters: 278e662bc50SJakub Kruzik + pc - the preconditioner context 279e662bc50SJakub Kruzik . W - deflation matrix 280e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 281e662bc50SJakub Kruzik 282e662bc50SJakub Kruzik Level: intermediate 283e662bc50SJakub Kruzik 284e662bc50SJakub Kruzik .seealso: PCDEFLATION 285e662bc50SJakub Kruzik @*/ 286e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 287e662bc50SJakub Kruzik { 288e662bc50SJakub Kruzik PetscErrorCode ierr; 289e662bc50SJakub Kruzik 290e662bc50SJakub Kruzik PetscFunctionBegin; 291e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 292e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 293e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 294e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 295e662bc50SJakub Kruzik PetscFunctionReturn(0); 296e662bc50SJakub Kruzik } 297e662bc50SJakub Kruzik 2983cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 2993cf3a049SJakub Kruzik { 3003cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3013cf3a049SJakub Kruzik PetscErrorCode ierr; 3023cf3a049SJakub Kruzik 3033cf3a049SJakub Kruzik PetscFunctionBegin; 3043cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 3053cf3a049SJakub Kruzik def->WtA = mat; 3063cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 3073cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 3083cf3a049SJakub Kruzik PetscFunctionReturn(0); 3093cf3a049SJakub Kruzik } 3103cf3a049SJakub Kruzik 3113cf3a049SJakub Kruzik /*@ 3123cf3a049SJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 3133cf3a049SJakub Kruzik 3143cf3a049SJakub Kruzik Not Collective 3153cf3a049SJakub Kruzik 3163cf3a049SJakub Kruzik Input Parameters: 3173cf3a049SJakub Kruzik + pc - preconditioner context 3183cf3a049SJakub Kruzik - mat - projection null space matrix 3193cf3a049SJakub Kruzik 3203cf3a049SJakub Kruzik Level: developer 3213cf3a049SJakub Kruzik 3223cf3a049SJakub Kruzik .seealso: PCDEFLATION 3233cf3a049SJakub Kruzik @*/ 3243cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 3253cf3a049SJakub Kruzik { 3263cf3a049SJakub Kruzik PetscErrorCode ierr; 3273cf3a049SJakub Kruzik 3283cf3a049SJakub Kruzik PetscFunctionBegin; 3293cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3303cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 3313cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 3323cf3a049SJakub Kruzik PetscFunctionReturn(0); 3333cf3a049SJakub Kruzik } 3343cf3a049SJakub Kruzik 335e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 336e924b002SJakub Kruzik { 337e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 338e924b002SJakub Kruzik PetscErrorCode ierr; 339e924b002SJakub Kruzik 340e924b002SJakub Kruzik PetscFunctionBegin; 341e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 342e924b002SJakub Kruzik def->WtAW = mat; 343e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 344e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 345e924b002SJakub Kruzik PetscFunctionReturn(0); 346e924b002SJakub Kruzik } 347e924b002SJakub Kruzik 348e924b002SJakub Kruzik /*@ 349e924b002SJakub Kruzik PCDeflationSetCoarseMat - Set coarse problem Mat. 350e924b002SJakub Kruzik 351e924b002SJakub Kruzik Not Collective 352e924b002SJakub Kruzik 353e924b002SJakub Kruzik Input Parameters: 354e924b002SJakub Kruzik + pc - preconditioner context 355e924b002SJakub Kruzik - mat - coarse problem mat 356e924b002SJakub Kruzik 357e924b002SJakub Kruzik Level: developer 358e924b002SJakub Kruzik 359e924b002SJakub Kruzik .seealso: PCDEFLATION 360e924b002SJakub Kruzik @*/ 361e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 362e924b002SJakub Kruzik { 363e924b002SJakub Kruzik PetscErrorCode ierr; 364e924b002SJakub Kruzik 365e924b002SJakub Kruzik PetscFunctionBegin; 366e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 367e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 368e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 369e924b002SJakub Kruzik PetscFunctionReturn(0); 370e924b002SJakub Kruzik } 371e924b002SJakub Kruzik 37298d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 3735c0c31c5SJakub Kruzik { 3745c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3755c0c31c5SJakub Kruzik 3765c0c31c5SJakub Kruzik PetscFunctionBegin; 37798d6e3deSJakub Kruzik *ksp = def->WtAWinv; 3785c0c31c5SJakub Kruzik PetscFunctionReturn(0); 3795c0c31c5SJakub Kruzik } 3805c0c31c5SJakub Kruzik 3815c0c31c5SJakub Kruzik /*@ 38298d6e3deSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 3835c0c31c5SJakub Kruzik 38498d6e3deSJakub Kruzik Not Collective 3855c0c31c5SJakub Kruzik 3865c0c31c5SJakub Kruzik Input Parameters: 38798d6e3deSJakub Kruzik . pc - preconditioner context 3885c0c31c5SJakub Kruzik 38998d6e3deSJakub Kruzik Output Parameter: 39098d6e3deSJakub Kruzik . ksp - coarse problem KSP context 3915c0c31c5SJakub Kruzik 39298d6e3deSJakub Kruzik Level: developer 39398d6e3deSJakub Kruzik 39498d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 3955c0c31c5SJakub Kruzik @*/ 39698d6e3deSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3975c0c31c5SJakub Kruzik { 3985c0c31c5SJakub Kruzik PetscErrorCode ierr; 3995c0c31c5SJakub Kruzik 4005c0c31c5SJakub Kruzik PetscFunctionBegin; 40122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40298d6e3deSJakub Kruzik PetscValidPointer(ksp,2); 40398d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 40498d6e3deSJakub Kruzik PetscFunctionReturn(0); 40598d6e3deSJakub Kruzik } 40698d6e3deSJakub Kruzik 40798d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 40898d6e3deSJakub Kruzik { 40998d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 41098d6e3deSJakub Kruzik PetscErrorCode ierr; 41198d6e3deSJakub Kruzik 41298d6e3deSJakub Kruzik PetscFunctionBegin; 41398d6e3deSJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 41498d6e3deSJakub Kruzik def->WtAWinv = ksp; 41598d6e3deSJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 41698d6e3deSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 41798d6e3deSJakub Kruzik PetscFunctionReturn(0); 41898d6e3deSJakub Kruzik } 41998d6e3deSJakub Kruzik 42098d6e3deSJakub Kruzik /*@ 42198d6e3deSJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 42298d6e3deSJakub Kruzik 42398d6e3deSJakub Kruzik Collective on PC 42498d6e3deSJakub Kruzik 42598d6e3deSJakub Kruzik Input Parameters: 42698d6e3deSJakub Kruzik + pc - preconditioner context 42798d6e3deSJakub Kruzik - ksp - coarse problem KSP context 42898d6e3deSJakub Kruzik 42998d6e3deSJakub Kruzik Level: developer 43098d6e3deSJakub Kruzik 43198d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 43298d6e3deSJakub Kruzik @*/ 43398d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 43498d6e3deSJakub Kruzik { 43598d6e3deSJakub Kruzik PetscErrorCode ierr; 43698d6e3deSJakub Kruzik 43798d6e3deSJakub Kruzik PetscFunctionBegin; 43898d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 43998d6e3deSJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 44098d6e3deSJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 44198d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 4425c0c31c5SJakub Kruzik PetscFunctionReturn(0); 4435c0c31c5SJakub Kruzik } 444e662bc50SJakub Kruzik 445268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 446268c6673SJakub Kruzik { 447268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 448268c6673SJakub Kruzik 449268c6673SJakub Kruzik PetscFunctionBegin; 450268c6673SJakub Kruzik *apc = def->pc; 451268c6673SJakub Kruzik PetscFunctionReturn(0); 452268c6673SJakub Kruzik } 453268c6673SJakub Kruzik 454268c6673SJakub Kruzik /*@ 455268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 456268c6673SJakub Kruzik 457268c6673SJakub Kruzik Not Collective 458268c6673SJakub Kruzik 459268c6673SJakub Kruzik Input Parameters: 460268c6673SJakub Kruzik . pc - the preconditioner context 461268c6673SJakub Kruzik 462268c6673SJakub Kruzik Output Parameter: 463268c6673SJakub Kruzik . apc - additional preconditioner 464268c6673SJakub Kruzik 465268c6673SJakub Kruzik Level: advanced 466268c6673SJakub Kruzik 467268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 468268c6673SJakub Kruzik @*/ 469268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 470268c6673SJakub Kruzik { 471268c6673SJakub Kruzik PetscErrorCode ierr; 472268c6673SJakub Kruzik 473268c6673SJakub Kruzik PetscFunctionBegin; 474268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 475268c6673SJakub Kruzik PetscValidPointer(pc,2); 476268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 477268c6673SJakub Kruzik PetscFunctionReturn(0); 478268c6673SJakub Kruzik } 479268c6673SJakub Kruzik 48022b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 48122b0793eSJakub Kruzik { 48222b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 48322b0793eSJakub Kruzik PetscErrorCode ierr; 48422b0793eSJakub Kruzik 48522b0793eSJakub Kruzik PetscFunctionBegin; 48622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 48722b0793eSJakub Kruzik def->pc = apc; 48822b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 48922b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 49022b0793eSJakub Kruzik PetscFunctionReturn(0); 49122b0793eSJakub Kruzik } 49222b0793eSJakub Kruzik 49322b0793eSJakub Kruzik /*@ 49422b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 49522b0793eSJakub Kruzik 496268c6673SJakub Kruzik Collective on PC 49722b0793eSJakub Kruzik 49822b0793eSJakub Kruzik Input Parameters: 49922b0793eSJakub Kruzik + pc - the preconditioner context 50022b0793eSJakub Kruzik - apc - additional preconditioner 50122b0793eSJakub Kruzik 502268c6673SJakub Kruzik Level: developer 50322b0793eSJakub Kruzik 504268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 50522b0793eSJakub Kruzik @*/ 50622b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 50722b0793eSJakub Kruzik { 50822b0793eSJakub Kruzik PetscErrorCode ierr; 50922b0793eSJakub Kruzik 51022b0793eSJakub Kruzik PetscFunctionBegin; 51122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 51222b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 51322b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 51422b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 51522b0793eSJakub Kruzik PetscFunctionReturn(0); 51622b0793eSJakub Kruzik } 51722b0793eSJakub Kruzik 51837eeb815SJakub Kruzik /* 519b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 520b27c8092SJakub Kruzik */ 521b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 522b27c8092SJakub Kruzik { 523b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 524b27c8092SJakub Kruzik Mat A; 525b27c8092SJakub Kruzik Vec r,w1,w2; 526b27c8092SJakub Kruzik PetscBool nonzero; 527b27c8092SJakub Kruzik PetscErrorCode ierr; 52837eeb815SJakub Kruzik 529b27c8092SJakub Kruzik PetscFunctionBegin; 530b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 531b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 532b27c8092SJakub Kruzik r = def->work; 533b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 53437eeb815SJakub Kruzik 535b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 536b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 537b27c8092SJakub Kruzik if (nonzero) { 538b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 539b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 540b27c8092SJakub Kruzik } else { 541b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 542b27c8092SJakub Kruzik } 543b27c8092SJakub Kruzik 544a3177931SJakub Kruzik if (def->Wt) { 545a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 546a3177931SJakub Kruzik } else { 547a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 548a3177931SJakub Kruzik } 549b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 550b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 551b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 552b27c8092SJakub Kruzik PetscFunctionReturn(0); 553b27c8092SJakub Kruzik } 55437eeb815SJakub Kruzik 555f8bf229cSJakub Kruzik /* 556f8bf229cSJakub Kruzik if (def->correct) { 557ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 558f8bf229cSJakub Kruzik } else { 559ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 560f8bf229cSJakub Kruzik } 56137eeb815SJakub Kruzik */ 562f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 563f8bf229cSJakub Kruzik { 564f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 565f8bf229cSJakub Kruzik Mat A; 566f8bf229cSJakub Kruzik Vec u,w1,w2; 567f8bf229cSJakub Kruzik PetscErrorCode ierr; 568f8bf229cSJakub Kruzik 569f8bf229cSJakub Kruzik PetscFunctionBegin; 570f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 571f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 572f8bf229cSJakub Kruzik u = def->work; 573f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 574f8bf229cSJakub Kruzik 575ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 57622b0793eSJakub Kruzik if (!def->init) { 577ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 578f8bf229cSJakub Kruzik if (def->correct) { 579ae029463SJakub Kruzik if (def->Wt) { 580ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 581ae029463SJakub Kruzik } else { 582a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 583f8bf229cSJakub Kruzik } 584ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 585f8bf229cSJakub Kruzik } 586f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 587f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 58822b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 58922b0793eSJakub Kruzik } 590f8bf229cSJakub Kruzik PetscFunctionReturn(0); 591f8bf229cSJakub Kruzik } 592f8bf229cSJakub Kruzik 59337eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 59437eeb815SJakub Kruzik { 595409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 596409a357bSJakub Kruzik KSP innerksp; 597409a357bSJakub Kruzik PC pcinner; 598409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 599409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 600409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 601409a357bSJakub Kruzik MatCompositeType ctype; 602409a357bSJakub Kruzik MPI_Comm comm; 603409a357bSJakub Kruzik const char *prefix; 60437eeb815SJakub Kruzik PetscErrorCode ierr; 60537eeb815SJakub Kruzik 60637eeb815SJakub Kruzik PetscFunctionBegin; 607409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 608409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 609f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 610f8bf229cSJakub Kruzik 611f8bf229cSJakub Kruzik /* compute a deflation space */ 612409a357bSJakub Kruzik if (def->W || def->Wt) { 613409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 614409a357bSJakub Kruzik } else { 615e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 61637eeb815SJakub Kruzik } 61737eeb815SJakub Kruzik 618409a357bSJakub Kruzik /* nested deflation */ 619409a357bSJakub Kruzik if (def->W) { 620409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 621409a357bSJakub Kruzik if (match) { 622409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 623409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 62437eeb815SJakub Kruzik } 625409a357bSJakub Kruzik } else { 626a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 627409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 628409a357bSJakub Kruzik if (match) { 629409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 630409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 631409a357bSJakub Kruzik } 632409a357bSJakub Kruzik transp = PETSC_TRUE; 633409a357bSJakub Kruzik } 634409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 635409a357bSJakub Kruzik if (!transp) { 63628da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 63728da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 638409a357bSJakub Kruzik for (i=0; i<size; i++) { 639409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 640409a357bSJakub Kruzik } 641409a357bSJakub Kruzik size -= 1; 642409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 643409a357bSJakub Kruzik def->W = mats[size]; 644409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 645409a357bSJakub Kruzik if (size > 1) { 646409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 647409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 648409a357bSJakub Kruzik } else { 649409a357bSJakub Kruzik nextDef = mats[0]; 650409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 651409a357bSJakub Kruzik } 65228da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 653409a357bSJakub Kruzik } else { 654409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 655409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 656409a357bSJakub Kruzik } 65728da0170SJakub Kruzik } else { 65828da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 65928da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 66028da0170SJakub Kruzik for (i=0; i<size; i++) { 66128da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 66228da0170SJakub Kruzik } 66328da0170SJakub Kruzik size -= 1; 66428da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 66528da0170SJakub Kruzik def->Wt = mats[0]; 66628da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 66728da0170SJakub Kruzik if (size > 1) { 66828da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 66928da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 67028da0170SJakub Kruzik } else { 67128da0170SJakub Kruzik nextDef = mats[1]; 67228da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 673409a357bSJakub Kruzik } 674409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 67528da0170SJakub Kruzik } else { 67628da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 67728da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 67828da0170SJakub Kruzik } 67928da0170SJakub Kruzik } 68028da0170SJakub Kruzik } 68128da0170SJakub Kruzik 68228da0170SJakub Kruzik if (transp) { 68328da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 684a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 685409a357bSJakub Kruzik } 686409a357bSJakub Kruzik 68722b0793eSJakub Kruzik 68822b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 68922b0793eSJakub Kruzik 690ae029463SJakub Kruzik /* assemble WtA */ 691ae029463SJakub Kruzik if (!def->WtA) { 692ae029463SJakub Kruzik if (def->Wt) { 693ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 694ae029463SJakub Kruzik } else { 695a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 696a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 697a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 698a3177931SJakub Kruzik #else 699ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 700a3177931SJakub Kruzik #endif 701ae029463SJakub Kruzik } 702ae029463SJakub Kruzik } 703409a357bSJakub Kruzik /* setup coarse problem */ 704409a357bSJakub Kruzik if (!def->WtAWinv) { 70528da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 706409a357bSJakub Kruzik if (!def->WtAW) { 707ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 708409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 709409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 710409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 711409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 712ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 713409a357bSJakub Kruzik PetscReal *norms; 714409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 715409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 716409a357bSJakub Kruzik for (i=0; i<m; i++) { 717409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 718409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 719409a357bSJakub Kruzik } 720409a357bSJakub Kruzik } 721409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 722409a357bSJakub Kruzik #endif 723409a357bSJakub Kruzik } else { 724409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 725409a357bSJakub Kruzik } 726409a357bSJakub Kruzik /* TODO use MATINV */ 727409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 728409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 729409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 730557a2f7dSJakub Kruzik /* Setup KSP and PC */ 731557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 732557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 733557a2f7dSJakub Kruzik /* set default KSPtype */ 734557a2f7dSJakub Kruzik if (!def->ksptype) { 735557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 736557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 737557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 738557a2f7dSJakub Kruzik } 739557a2f7dSJakub Kruzik } 740ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 741557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 742557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 743557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 744ae029463SJakub Kruzik /* inherit options */ 745557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 746557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 747557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 748557a2f7dSJakub Kruzik } else { /* the last level */ 749557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 750409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 751ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 752ac66f006SJakub Kruzik if (prefix) { 753ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 754ac66f006SJakub Kruzik } 75522b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 756ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 757409a357bSJakub Kruzik /* Reduction factor choice */ 758409a357bSJakub Kruzik red = def->reductionfact; 759409a357bSJakub Kruzik if (red < 0) { 760409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 761409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 762409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 763409a357bSJakub Kruzik if (match) red = commsize; 764409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 765409a357bSJakub Kruzik } 766409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 767ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 768409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 769ac66f006SJakub Kruzik if (innerksp) { 770409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 771409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 772ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 773409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 774409a357bSJakub Kruzik if (commsize == red) { 775ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 776409a357bSJakub Kruzik } else { 777ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 778409a357bSJakub Kruzik } 779409a357bSJakub Kruzik } 780557a2f7dSJakub Kruzik } 781557a2f7dSJakub Kruzik 782557a2f7dSJakub Kruzik if (innerksp) { 783409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 784409a357bSJakub Kruzik if (prefix) { 785409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 786409a357bSJakub Kruzik } 78722b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 788557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 789557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 790ac66f006SJakub Kruzik } 791409a357bSJakub Kruzik } 792557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 793557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 794409a357bSJakub Kruzik 79522b0793eSJakub Kruzik if (!def->pc) { 79622b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 79722b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 79822b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 79922b0793eSJakub Kruzik if (prefix) { 80022b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 80122b0793eSJakub Kruzik } 80222b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 80322b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 80422b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 80522b0793eSJakub Kruzik } 80622b0793eSJakub Kruzik 807f8bf229cSJakub Kruzik /* create work vecs */ 808f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 809f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 81037eeb815SJakub Kruzik PetscFunctionReturn(0); 81137eeb815SJakub Kruzik } 812b30d4775SJakub Kruzik 81337eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 81437eeb815SJakub Kruzik { 815b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 81637eeb815SJakub Kruzik PetscErrorCode ierr; 81737eeb815SJakub Kruzik 81837eeb815SJakub Kruzik PetscFunctionBegin; 819b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 820b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 821b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 822b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 823ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 824b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 825b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 82622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 82737eeb815SJakub Kruzik PetscFunctionReturn(0); 82837eeb815SJakub Kruzik } 82937eeb815SJakub Kruzik 83037eeb815SJakub Kruzik /* 83137eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 83237eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 83337eeb815SJakub Kruzik 83437eeb815SJakub Kruzik Input Parameter: 83537eeb815SJakub Kruzik . pc - the preconditioner context 83637eeb815SJakub Kruzik 83737eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 83837eeb815SJakub Kruzik */ 83937eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 84037eeb815SJakub Kruzik { 84137eeb815SJakub Kruzik PetscErrorCode ierr; 84237eeb815SJakub Kruzik 84337eeb815SJakub Kruzik PetscFunctionBegin; 84437eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 84537eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 84637eeb815SJakub Kruzik PetscFunctionReturn(0); 84737eeb815SJakub Kruzik } 84837eeb815SJakub Kruzik 8498513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 8508513960bSJakub Kruzik { 8518513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 8528513960bSJakub Kruzik PetscBool iascii; 8538513960bSJakub Kruzik PetscErrorCode ierr; 8548513960bSJakub Kruzik 8558513960bSJakub Kruzik PetscFunctionBegin; 8568513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8578513960bSJakub Kruzik if (iascii) { 8588513960bSJakub Kruzik //if (cg->adaptiveconv) { 8598513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 8608513960bSJakub Kruzik //} 8618513960bSJakub Kruzik if (def->correct) { 8628513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 8638513960bSJakub Kruzik } 8648513960bSJakub Kruzik if (!def->nestedlvl) { 8658513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 8668513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 8678513960bSJakub Kruzik } 8688513960bSJakub Kruzik 86922b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 87022b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 87122b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 87222b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 87322b0793eSJakub Kruzik 8748513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 8758513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 8768513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 8778513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 8788513960bSJakub Kruzik } 8798513960bSJakub Kruzik PetscFunctionReturn(0); 8808513960bSJakub Kruzik } 8818513960bSJakub Kruzik 88237eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 88337eeb815SJakub Kruzik { 884880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 88537eeb815SJakub Kruzik PetscErrorCode ierr; 88637eeb815SJakub Kruzik 88737eeb815SJakub Kruzik PetscFunctionBegin; 88837eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 889*a122ebaeSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr); 890859a873cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 891859a873cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 8928a71cb68SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr); 8938a71cb68SJakub Kruzik ierr = PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr); 894880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 895880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 896880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 897880d05ceSJakub Kruzik //TODO add set function and fix manpages 89837eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 89937eeb815SJakub Kruzik PetscFunctionReturn(0); 90037eeb815SJakub Kruzik } 90137eeb815SJakub Kruzik 90237eeb815SJakub Kruzik /*MC 903e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 904e361f8b5SJakub Kruzik or to a predefined value 90537eeb815SJakub Kruzik 90637eeb815SJakub Kruzik Options Database Key: 907e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 90837eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 90937eeb815SJakub Kruzik 91037eeb815SJakub Kruzik Level: beginner 91137eeb815SJakub Kruzik 91237eeb815SJakub Kruzik Notes: 913e361f8b5SJakub Kruzik todo 91437eeb815SJakub Kruzik 91537eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 916e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 91737eeb815SJakub Kruzik M*/ 91837eeb815SJakub Kruzik 91937eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 92037eeb815SJakub Kruzik { 92137eeb815SJakub Kruzik PC_Deflation *def; 92237eeb815SJakub Kruzik PetscErrorCode ierr; 92337eeb815SJakub Kruzik 92437eeb815SJakub Kruzik PetscFunctionBegin; 92537eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 92637eeb815SJakub Kruzik pc->data = (void*)def; 92737eeb815SJakub Kruzik 928e662bc50SJakub Kruzik def->init = PETSC_FALSE; 929e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 930fcb31d99SJakub Kruzik def->correctfact = 1.0; 931e662bc50SJakub Kruzik def->reductionfact = -1; 932e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 933e662bc50SJakub Kruzik def->spacesize = 1; 934e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 935e662bc50SJakub Kruzik def->nestedlvl = 0; 936e662bc50SJakub Kruzik def->maxnestedlvl = 0; 93737eeb815SJakub Kruzik 93837eeb815SJakub Kruzik /* 93937eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 94037eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 94137eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 94237eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 94337eeb815SJakub Kruzik not needed. 94437eeb815SJakub Kruzik */ 94537eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 94637eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 947b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 948b27c8092SJakub Kruzik pc->ops->postsolve = 0; 94937eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 95037eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 95137eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 95237eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 9538513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 95437eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 95537eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 95637eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 95737eeb815SJakub Kruzik 958*a122ebaeSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr); 95998d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 960859a873cSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr); 9618a71cb68SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr); 96239417d7eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 963e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 9643cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 965e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 9664a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 967f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 96898d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 96998d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 97037eeb815SJakub Kruzik PetscFunctionReturn(0); 97137eeb815SJakub Kruzik } 97237eeb815SJakub Kruzik 973