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 7198d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 7298d6e3deSJakub Kruzik { 7398d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 7498d6e3deSJakub Kruzik 7598d6e3deSJakub Kruzik PetscFunctionBegin; 7698d6e3deSJakub Kruzik if (current) def->nestedlvl = current; 7798d6e3deSJakub Kruzik def->maxnestedlvl = max; 7898d6e3deSJakub Kruzik PetscFunctionReturn(0); 7998d6e3deSJakub Kruzik } 8098d6e3deSJakub Kruzik 8198d6e3deSJakub Kruzik /*@ 8298d6e3deSJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 8398d6e3deSJakub Kruzik 8498d6e3deSJakub Kruzik Logically Collective on PC 8598d6e3deSJakub Kruzik 8698d6e3deSJakub Kruzik Input Parameters: 8798d6e3deSJakub Kruzik + pc - the preconditioner context 8898d6e3deSJakub Kruzik - max - maximum deflation level 8998d6e3deSJakub Kruzik 9098d6e3deSJakub Kruzik Level: intermediate 9198d6e3deSJakub Kruzik 9298d6e3deSJakub Kruzik .seealso: PCDEFLATION 9398d6e3deSJakub Kruzik @*/ 9498d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 9598d6e3deSJakub Kruzik { 9698d6e3deSJakub Kruzik PetscErrorCode ierr; 9798d6e3deSJakub Kruzik 9898d6e3deSJakub Kruzik PetscFunctionBegin; 9998d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10098d6e3deSJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 10198d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 10298d6e3deSJakub Kruzik PetscFunctionReturn(0); 10398d6e3deSJakub Kruzik } 10498d6e3deSJakub Kruzik 105859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red) 106859a873cSJakub Kruzik { 107859a873cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 108859a873cSJakub Kruzik 109859a873cSJakub Kruzik PetscFunctionBegin; 110859a873cSJakub Kruzik def->reductionfact = red; 111859a873cSJakub Kruzik PetscFunctionReturn(0); 112859a873cSJakub Kruzik } 113859a873cSJakub Kruzik 114859a873cSJakub Kruzik /*@ 115859a873cSJakub Kruzik PCDeflationSetReductionFactor - Set reduction factor for PCTELESCOPE coarse problem solver. 116859a873cSJakub Kruzik 117859a873cSJakub Kruzik Logically Collective on PC 118859a873cSJakub Kruzik 119859a873cSJakub Kruzik Input Parameters: 120859a873cSJakub Kruzik + pc - the preconditioner context 121859a873cSJakub Kruzik - red - reduction factor (or PETSC_DETERMINE) 122859a873cSJakub Kruzik 123859a873cSJakub Kruzik Level: intermediate 124859a873cSJakub Kruzik 125859a873cSJakub Kruzik .seealso: PCDEFLATION 126859a873cSJakub Kruzik @*/ 127859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red) 128859a873cSJakub Kruzik { 129859a873cSJakub Kruzik PetscErrorCode ierr; 130859a873cSJakub Kruzik 131859a873cSJakub Kruzik PetscFunctionBegin; 132859a873cSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 133859a873cSJakub Kruzik PetscValidLogicalCollectiveInt(pc,red,2); 134859a873cSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr); 135859a873cSJakub Kruzik PetscFunctionReturn(0); 136859a873cSJakub Kruzik } 137859a873cSJakub Kruzik 138*8a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact) 139*8a71cb68SJakub Kruzik { 140*8a71cb68SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 141*8a71cb68SJakub Kruzik 142*8a71cb68SJakub Kruzik PetscFunctionBegin; 143*8a71cb68SJakub Kruzik /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 144*8a71cb68SJakub Kruzik def->correct = PETSC_TRUE; 145*8a71cb68SJakub Kruzik def->correctfact = fact; 146*8a71cb68SJakub Kruzik if (fact == PETSC_DEFAULT) { 147*8a71cb68SJakub Kruzik def->correctfact = 1.0; 148*8a71cb68SJakub Kruzik } else if (def->correct == 0.0) { 149*8a71cb68SJakub Kruzik def->correct = PETSC_FALSE; 150*8a71cb68SJakub Kruzik } 151*8a71cb68SJakub Kruzik PetscFunctionReturn(0); 152*8a71cb68SJakub Kruzik } 153*8a71cb68SJakub Kruzik 154*8a71cb68SJakub Kruzik /*@ 155*8a71cb68SJakub Kruzik PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 156*8a71cb68SJakub Kruzik The Preconditioner becomes P*M^{-1} + factor*Q. 157*8a71cb68SJakub Kruzik 158*8a71cb68SJakub Kruzik Logically Collective on PC 159*8a71cb68SJakub Kruzik 160*8a71cb68SJakub Kruzik Input Parameters: 161*8a71cb68SJakub Kruzik + pc - the preconditioner context 162*8a71cb68SJakub Kruzik - fact - correction factor (set 0.0 to disable, PETSC_DEFAULT = 1.0) 163*8a71cb68SJakub Kruzik 164*8a71cb68SJakub Kruzik Level: intermediate 165*8a71cb68SJakub Kruzik 166*8a71cb68SJakub Kruzik .seealso: PCDEFLATION 167*8a71cb68SJakub Kruzik @*/ 168*8a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact) 169*8a71cb68SJakub Kruzik { 170*8a71cb68SJakub Kruzik PetscErrorCode ierr; 171*8a71cb68SJakub Kruzik 172*8a71cb68SJakub Kruzik PetscFunctionBegin; 173*8a71cb68SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 174*8a71cb68SJakub Kruzik PetscValidLogicalCollectiveScalar(pc,fact,2); 175*8a71cb68SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr); 176*8a71cb68SJakub Kruzik PetscFunctionReturn(0); 177*8a71cb68SJakub Kruzik } 178*8a71cb68SJakub Kruzik 17939417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 18039417d7eSJakub Kruzik { 18139417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 18239417d7eSJakub Kruzik 18339417d7eSJakub Kruzik PetscFunctionBegin; 18439417d7eSJakub Kruzik if (type) def->spacetype = type; 18539417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 18639417d7eSJakub Kruzik PetscFunctionReturn(0); 18739417d7eSJakub Kruzik } 18839417d7eSJakub Kruzik 18939417d7eSJakub Kruzik /*@ 19039417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 19139417d7eSJakub Kruzik 19239417d7eSJakub Kruzik Logically Collective on PC 19339417d7eSJakub Kruzik 19439417d7eSJakub Kruzik Input Parameters: 19539417d7eSJakub Kruzik + pc - the preconditioner context 19639417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 19739417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 19839417d7eSJakub Kruzik 19939417d7eSJakub Kruzik Notes: 20039417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 20139417d7eSJakub Kruzik The deflation space is computed in PCSetUP(). 20239417d7eSJakub Kruzik 20339417d7eSJakub Kruzik Level: intermediate 20439417d7eSJakub Kruzik 20539417d7eSJakub Kruzik .seealso: PCDEFLATION 20639417d7eSJakub Kruzik @*/ 20739417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 20839417d7eSJakub Kruzik { 20939417d7eSJakub Kruzik PetscErrorCode ierr; 21039417d7eSJakub Kruzik 21139417d7eSJakub Kruzik PetscFunctionBegin; 21239417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 21339417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 21439417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 21539417d7eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 21639417d7eSJakub Kruzik PetscFunctionReturn(0); 21739417d7eSJakub Kruzik } 2188513960bSJakub Kruzik 219e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 220e662bc50SJakub Kruzik { 221e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 222e662bc50SJakub Kruzik PetscErrorCode ierr; 223e662bc50SJakub Kruzik 224e662bc50SJakub Kruzik PetscFunctionBegin; 225e662bc50SJakub Kruzik if (transpose) { 226e662bc50SJakub Kruzik def->Wt = W; 227e662bc50SJakub Kruzik def->W = NULL; 228e662bc50SJakub Kruzik } else { 229e662bc50SJakub Kruzik def->W = W; 230e662bc50SJakub Kruzik } 231e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 232e662bc50SJakub Kruzik PetscFunctionReturn(0); 233e662bc50SJakub Kruzik } 234e662bc50SJakub Kruzik 235e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 236e662bc50SJakub Kruzik /*@ 237e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 238e662bc50SJakub Kruzik 239e662bc50SJakub Kruzik Logically Collective on PC 240e662bc50SJakub Kruzik 241e662bc50SJakub Kruzik Input Parameters: 242e662bc50SJakub Kruzik + pc - the preconditioner context 243e662bc50SJakub Kruzik . W - deflation matrix 244e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 245e662bc50SJakub Kruzik 246e662bc50SJakub Kruzik Level: intermediate 247e662bc50SJakub Kruzik 248e662bc50SJakub Kruzik .seealso: PCDEFLATION 249e662bc50SJakub Kruzik @*/ 250e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 251e662bc50SJakub Kruzik { 252e662bc50SJakub Kruzik PetscErrorCode ierr; 253e662bc50SJakub Kruzik 254e662bc50SJakub Kruzik PetscFunctionBegin; 255e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 256e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 257e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 258e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 259e662bc50SJakub Kruzik PetscFunctionReturn(0); 260e662bc50SJakub Kruzik } 261e662bc50SJakub Kruzik 2623cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 2633cf3a049SJakub Kruzik { 2643cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2653cf3a049SJakub Kruzik PetscErrorCode ierr; 2663cf3a049SJakub Kruzik 2673cf3a049SJakub Kruzik PetscFunctionBegin; 2683cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 2693cf3a049SJakub Kruzik def->WtA = mat; 2703cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2713cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 2723cf3a049SJakub Kruzik PetscFunctionReturn(0); 2733cf3a049SJakub Kruzik } 2743cf3a049SJakub Kruzik 2753cf3a049SJakub Kruzik /*@ 2763cf3a049SJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 2773cf3a049SJakub Kruzik 2783cf3a049SJakub Kruzik Not Collective 2793cf3a049SJakub Kruzik 2803cf3a049SJakub Kruzik Input Parameters: 2813cf3a049SJakub Kruzik + pc - preconditioner context 2823cf3a049SJakub Kruzik - mat - projection null space matrix 2833cf3a049SJakub Kruzik 2843cf3a049SJakub Kruzik Level: developer 2853cf3a049SJakub Kruzik 2863cf3a049SJakub Kruzik .seealso: PCDEFLATION 2873cf3a049SJakub Kruzik @*/ 2883cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 2893cf3a049SJakub Kruzik { 2903cf3a049SJakub Kruzik PetscErrorCode ierr; 2913cf3a049SJakub Kruzik 2923cf3a049SJakub Kruzik PetscFunctionBegin; 2933cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2943cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 2953cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 2963cf3a049SJakub Kruzik PetscFunctionReturn(0); 2973cf3a049SJakub Kruzik } 2983cf3a049SJakub Kruzik 299e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 300e924b002SJakub Kruzik { 301e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 302e924b002SJakub Kruzik PetscErrorCode ierr; 303e924b002SJakub Kruzik 304e924b002SJakub Kruzik PetscFunctionBegin; 305e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 306e924b002SJakub Kruzik def->WtAW = mat; 307e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 308e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 309e924b002SJakub Kruzik PetscFunctionReturn(0); 310e924b002SJakub Kruzik } 311e924b002SJakub Kruzik 312e924b002SJakub Kruzik /*@ 313e924b002SJakub Kruzik PCDeflationSetCoarseMat - Set coarse problem Mat. 314e924b002SJakub Kruzik 315e924b002SJakub Kruzik Not Collective 316e924b002SJakub Kruzik 317e924b002SJakub Kruzik Input Parameters: 318e924b002SJakub Kruzik + pc - preconditioner context 319e924b002SJakub Kruzik - mat - coarse problem mat 320e924b002SJakub Kruzik 321e924b002SJakub Kruzik Level: developer 322e924b002SJakub Kruzik 323e924b002SJakub Kruzik .seealso: PCDEFLATION 324e924b002SJakub Kruzik @*/ 325e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 326e924b002SJakub Kruzik { 327e924b002SJakub Kruzik PetscErrorCode ierr; 328e924b002SJakub Kruzik 329e924b002SJakub Kruzik PetscFunctionBegin; 330e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 331e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 332e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 333e924b002SJakub Kruzik PetscFunctionReturn(0); 334e924b002SJakub Kruzik } 335e924b002SJakub Kruzik 33698d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 3375c0c31c5SJakub Kruzik { 3385c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3395c0c31c5SJakub Kruzik 3405c0c31c5SJakub Kruzik PetscFunctionBegin; 34198d6e3deSJakub Kruzik *ksp = def->WtAWinv; 3425c0c31c5SJakub Kruzik PetscFunctionReturn(0); 3435c0c31c5SJakub Kruzik } 3445c0c31c5SJakub Kruzik 3455c0c31c5SJakub Kruzik /*@ 34698d6e3deSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 3475c0c31c5SJakub Kruzik 34898d6e3deSJakub Kruzik Not Collective 3495c0c31c5SJakub Kruzik 3505c0c31c5SJakub Kruzik Input Parameters: 35198d6e3deSJakub Kruzik . pc - preconditioner context 3525c0c31c5SJakub Kruzik 35398d6e3deSJakub Kruzik Output Parameter: 35498d6e3deSJakub Kruzik . ksp - coarse problem KSP context 3555c0c31c5SJakub Kruzik 35698d6e3deSJakub Kruzik Level: developer 35798d6e3deSJakub Kruzik 35898d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 3595c0c31c5SJakub Kruzik @*/ 36098d6e3deSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3615c0c31c5SJakub Kruzik { 3625c0c31c5SJakub Kruzik PetscErrorCode ierr; 3635c0c31c5SJakub Kruzik 3645c0c31c5SJakub Kruzik PetscFunctionBegin; 36522b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 36698d6e3deSJakub Kruzik PetscValidPointer(ksp,2); 36798d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 36898d6e3deSJakub Kruzik PetscFunctionReturn(0); 36998d6e3deSJakub Kruzik } 37098d6e3deSJakub Kruzik 37198d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 37298d6e3deSJakub Kruzik { 37398d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 37498d6e3deSJakub Kruzik PetscErrorCode ierr; 37598d6e3deSJakub Kruzik 37698d6e3deSJakub Kruzik PetscFunctionBegin; 37798d6e3deSJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 37898d6e3deSJakub Kruzik def->WtAWinv = ksp; 37998d6e3deSJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 38098d6e3deSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 38198d6e3deSJakub Kruzik PetscFunctionReturn(0); 38298d6e3deSJakub Kruzik } 38398d6e3deSJakub Kruzik 38498d6e3deSJakub Kruzik /*@ 38598d6e3deSJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 38698d6e3deSJakub Kruzik 38798d6e3deSJakub Kruzik Collective on PC 38898d6e3deSJakub Kruzik 38998d6e3deSJakub Kruzik Input Parameters: 39098d6e3deSJakub Kruzik + pc - preconditioner context 39198d6e3deSJakub Kruzik - ksp - coarse problem KSP context 39298d6e3deSJakub Kruzik 39398d6e3deSJakub Kruzik Level: developer 39498d6e3deSJakub Kruzik 39598d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 39698d6e3deSJakub Kruzik @*/ 39798d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 39898d6e3deSJakub Kruzik { 39998d6e3deSJakub Kruzik PetscErrorCode ierr; 40098d6e3deSJakub Kruzik 40198d6e3deSJakub Kruzik PetscFunctionBegin; 40298d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40398d6e3deSJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 40498d6e3deSJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 40598d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 4065c0c31c5SJakub Kruzik PetscFunctionReturn(0); 4075c0c31c5SJakub Kruzik } 408e662bc50SJakub Kruzik 409268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 410268c6673SJakub Kruzik { 411268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 412268c6673SJakub Kruzik 413268c6673SJakub Kruzik PetscFunctionBegin; 414268c6673SJakub Kruzik *apc = def->pc; 415268c6673SJakub Kruzik PetscFunctionReturn(0); 416268c6673SJakub Kruzik } 417268c6673SJakub Kruzik 418268c6673SJakub Kruzik /*@ 419268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 420268c6673SJakub Kruzik 421268c6673SJakub Kruzik Not Collective 422268c6673SJakub Kruzik 423268c6673SJakub Kruzik Input Parameters: 424268c6673SJakub Kruzik . pc - the preconditioner context 425268c6673SJakub Kruzik 426268c6673SJakub Kruzik Output Parameter: 427268c6673SJakub Kruzik . apc - additional preconditioner 428268c6673SJakub Kruzik 429268c6673SJakub Kruzik Level: advanced 430268c6673SJakub Kruzik 431268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 432268c6673SJakub Kruzik @*/ 433268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 434268c6673SJakub Kruzik { 435268c6673SJakub Kruzik PetscErrorCode ierr; 436268c6673SJakub Kruzik 437268c6673SJakub Kruzik PetscFunctionBegin; 438268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 439268c6673SJakub Kruzik PetscValidPointer(pc,2); 440268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 441268c6673SJakub Kruzik PetscFunctionReturn(0); 442268c6673SJakub Kruzik } 443268c6673SJakub Kruzik 44422b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 44522b0793eSJakub Kruzik { 44622b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 44722b0793eSJakub Kruzik PetscErrorCode ierr; 44822b0793eSJakub Kruzik 44922b0793eSJakub Kruzik PetscFunctionBegin; 45022b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 45122b0793eSJakub Kruzik def->pc = apc; 45222b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 45322b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 45422b0793eSJakub Kruzik PetscFunctionReturn(0); 45522b0793eSJakub Kruzik } 45622b0793eSJakub Kruzik 45722b0793eSJakub Kruzik /*@ 45822b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 45922b0793eSJakub Kruzik 460268c6673SJakub Kruzik Collective on PC 46122b0793eSJakub Kruzik 46222b0793eSJakub Kruzik Input Parameters: 46322b0793eSJakub Kruzik + pc - the preconditioner context 46422b0793eSJakub Kruzik - apc - additional preconditioner 46522b0793eSJakub Kruzik 466268c6673SJakub Kruzik Level: developer 46722b0793eSJakub Kruzik 468268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 46922b0793eSJakub Kruzik @*/ 47022b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 47122b0793eSJakub Kruzik { 47222b0793eSJakub Kruzik PetscErrorCode ierr; 47322b0793eSJakub Kruzik 47422b0793eSJakub Kruzik PetscFunctionBegin; 47522b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 47622b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 47722b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 47822b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 47922b0793eSJakub Kruzik PetscFunctionReturn(0); 48022b0793eSJakub Kruzik } 48122b0793eSJakub Kruzik 48237eeb815SJakub Kruzik /* 483b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 484b27c8092SJakub Kruzik */ 485b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 486b27c8092SJakub Kruzik { 487b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 488b27c8092SJakub Kruzik Mat A; 489b27c8092SJakub Kruzik Vec r,w1,w2; 490b27c8092SJakub Kruzik PetscBool nonzero; 491b27c8092SJakub Kruzik PetscErrorCode ierr; 49237eeb815SJakub Kruzik 493b27c8092SJakub Kruzik PetscFunctionBegin; 494b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 495b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 496b27c8092SJakub Kruzik r = def->work; 497b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 49837eeb815SJakub Kruzik 499b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 500b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 501b27c8092SJakub Kruzik if (nonzero) { 502b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 503b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 504b27c8092SJakub Kruzik } else { 505b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 506b27c8092SJakub Kruzik } 507b27c8092SJakub Kruzik 508a3177931SJakub Kruzik if (def->Wt) { 509a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 510a3177931SJakub Kruzik } else { 511a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 512a3177931SJakub Kruzik } 513b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 514b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 515b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 516b27c8092SJakub Kruzik PetscFunctionReturn(0); 517b27c8092SJakub Kruzik } 51837eeb815SJakub Kruzik 519f8bf229cSJakub Kruzik /* 520f8bf229cSJakub Kruzik if (def->correct) { 521ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 522f8bf229cSJakub Kruzik } else { 523ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 524f8bf229cSJakub Kruzik } 52537eeb815SJakub Kruzik */ 526f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 527f8bf229cSJakub Kruzik { 528f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 529f8bf229cSJakub Kruzik Mat A; 530f8bf229cSJakub Kruzik Vec u,w1,w2; 531f8bf229cSJakub Kruzik PetscErrorCode ierr; 532f8bf229cSJakub Kruzik 533f8bf229cSJakub Kruzik PetscFunctionBegin; 534f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 535f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 536f8bf229cSJakub Kruzik u = def->work; 537f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 538f8bf229cSJakub Kruzik 539ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 54022b0793eSJakub Kruzik if (!def->init) { 541ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 542f8bf229cSJakub Kruzik if (def->correct) { 543ae029463SJakub Kruzik if (def->Wt) { 544ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 545ae029463SJakub Kruzik } else { 546a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 547f8bf229cSJakub Kruzik } 548ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 549f8bf229cSJakub Kruzik } 550f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 551f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 55222b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 55322b0793eSJakub Kruzik } 554f8bf229cSJakub Kruzik PetscFunctionReturn(0); 555f8bf229cSJakub Kruzik } 556f8bf229cSJakub Kruzik 55737eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 55837eeb815SJakub Kruzik { 559409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 560409a357bSJakub Kruzik KSP innerksp; 561409a357bSJakub Kruzik PC pcinner; 562409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 563409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 564409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 565409a357bSJakub Kruzik MatCompositeType ctype; 566409a357bSJakub Kruzik MPI_Comm comm; 567409a357bSJakub Kruzik const char *prefix; 56837eeb815SJakub Kruzik PetscErrorCode ierr; 56937eeb815SJakub Kruzik 57037eeb815SJakub Kruzik PetscFunctionBegin; 571409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 572409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 573f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 574f8bf229cSJakub Kruzik 575f8bf229cSJakub Kruzik /* compute a deflation space */ 576409a357bSJakub Kruzik if (def->W || def->Wt) { 577409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 578409a357bSJakub Kruzik } else { 579e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 58037eeb815SJakub Kruzik } 58137eeb815SJakub Kruzik 582409a357bSJakub Kruzik /* nested deflation */ 583409a357bSJakub Kruzik if (def->W) { 584409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 585409a357bSJakub Kruzik if (match) { 586409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 587409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 58837eeb815SJakub Kruzik } 589409a357bSJakub Kruzik } else { 590a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 591409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 592409a357bSJakub Kruzik if (match) { 593409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 594409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 595409a357bSJakub Kruzik } 596409a357bSJakub Kruzik transp = PETSC_TRUE; 597409a357bSJakub Kruzik } 598409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 599409a357bSJakub Kruzik if (!transp) { 60028da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 60128da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 602409a357bSJakub Kruzik for (i=0; i<size; i++) { 603409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 604409a357bSJakub Kruzik } 605409a357bSJakub Kruzik size -= 1; 606409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 607409a357bSJakub Kruzik def->W = mats[size]; 608409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 609409a357bSJakub Kruzik if (size > 1) { 610409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 611409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 612409a357bSJakub Kruzik } else { 613409a357bSJakub Kruzik nextDef = mats[0]; 614409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 615409a357bSJakub Kruzik } 61628da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 617409a357bSJakub Kruzik } else { 618409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 619409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 620409a357bSJakub Kruzik } 62128da0170SJakub Kruzik } else { 62228da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 62328da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 62428da0170SJakub Kruzik for (i=0; i<size; i++) { 62528da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 62628da0170SJakub Kruzik } 62728da0170SJakub Kruzik size -= 1; 62828da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 62928da0170SJakub Kruzik def->Wt = mats[0]; 63028da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 63128da0170SJakub Kruzik if (size > 1) { 63228da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 63328da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 63428da0170SJakub Kruzik } else { 63528da0170SJakub Kruzik nextDef = mats[1]; 63628da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 637409a357bSJakub Kruzik } 638409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 63928da0170SJakub Kruzik } else { 64028da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 64128da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 64228da0170SJakub Kruzik } 64328da0170SJakub Kruzik } 64428da0170SJakub Kruzik } 64528da0170SJakub Kruzik 64628da0170SJakub Kruzik if (transp) { 64728da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 648a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 649409a357bSJakub Kruzik } 650409a357bSJakub Kruzik 65122b0793eSJakub Kruzik 65222b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 65322b0793eSJakub Kruzik 654ae029463SJakub Kruzik /* assemble WtA */ 655ae029463SJakub Kruzik if (!def->WtA) { 656ae029463SJakub Kruzik if (def->Wt) { 657ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 658ae029463SJakub Kruzik } else { 659a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 660a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 661a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 662a3177931SJakub Kruzik #else 663ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 664a3177931SJakub Kruzik #endif 665ae029463SJakub Kruzik } 666ae029463SJakub Kruzik } 667409a357bSJakub Kruzik /* setup coarse problem */ 668409a357bSJakub Kruzik if (!def->WtAWinv) { 66928da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 670409a357bSJakub Kruzik if (!def->WtAW) { 671ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 672409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 673409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 674409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 675409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 676ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 677409a357bSJakub Kruzik PetscReal *norms; 678409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 679409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 680409a357bSJakub Kruzik for (i=0; i<m; i++) { 681409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 682409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 683409a357bSJakub Kruzik } 684409a357bSJakub Kruzik } 685409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 686409a357bSJakub Kruzik #endif 687409a357bSJakub Kruzik } else { 688409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 689409a357bSJakub Kruzik } 690409a357bSJakub Kruzik /* TODO use MATINV */ 691409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 692409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 693409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 694557a2f7dSJakub Kruzik /* Setup KSP and PC */ 695557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 696557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 697557a2f7dSJakub Kruzik /* set default KSPtype */ 698557a2f7dSJakub Kruzik if (!def->ksptype) { 699557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 700557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 701557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 702557a2f7dSJakub Kruzik } 703557a2f7dSJakub Kruzik } 704ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 705557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 706557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 707557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 708ae029463SJakub Kruzik /* inherit options */ 709557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 710557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 711557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 712557a2f7dSJakub Kruzik } else { /* the last level */ 713557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 714409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 715ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 716ac66f006SJakub Kruzik if (prefix) { 717ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 718ac66f006SJakub Kruzik } 71922b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 720ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 721409a357bSJakub Kruzik /* Reduction factor choice */ 722409a357bSJakub Kruzik red = def->reductionfact; 723409a357bSJakub Kruzik if (red < 0) { 724409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 725409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 726409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 727409a357bSJakub Kruzik if (match) red = commsize; 728409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 729409a357bSJakub Kruzik } 730409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 731ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 732409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 733ac66f006SJakub Kruzik if (innerksp) { 734409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 735409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 736ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 737409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 738409a357bSJakub Kruzik if (commsize == red) { 739ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 740409a357bSJakub Kruzik } else { 741ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 742409a357bSJakub Kruzik } 743409a357bSJakub Kruzik } 744557a2f7dSJakub Kruzik } 745557a2f7dSJakub Kruzik 746557a2f7dSJakub Kruzik if (innerksp) { 747409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 748409a357bSJakub Kruzik if (prefix) { 749409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 750409a357bSJakub Kruzik } 75122b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 752557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 753557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 754ac66f006SJakub Kruzik } 755409a357bSJakub Kruzik } 756557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 757557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 758409a357bSJakub Kruzik 75922b0793eSJakub Kruzik if (!def->pc) { 76022b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 76122b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 76222b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 76322b0793eSJakub Kruzik if (prefix) { 76422b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 76522b0793eSJakub Kruzik } 76622b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 76722b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 76822b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 76922b0793eSJakub Kruzik } 77022b0793eSJakub Kruzik 771f8bf229cSJakub Kruzik /* create work vecs */ 772f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 773f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 77437eeb815SJakub Kruzik PetscFunctionReturn(0); 77537eeb815SJakub Kruzik } 776b30d4775SJakub Kruzik 77737eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 77837eeb815SJakub Kruzik { 779b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 78037eeb815SJakub Kruzik PetscErrorCode ierr; 78137eeb815SJakub Kruzik 78237eeb815SJakub Kruzik PetscFunctionBegin; 783b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 784b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 785b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 786b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 787ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 788b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 789b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 79022b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 79137eeb815SJakub Kruzik PetscFunctionReturn(0); 79237eeb815SJakub Kruzik } 79337eeb815SJakub Kruzik 79437eeb815SJakub Kruzik /* 79537eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 79637eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 79737eeb815SJakub Kruzik 79837eeb815SJakub Kruzik Input Parameter: 79937eeb815SJakub Kruzik . pc - the preconditioner context 80037eeb815SJakub Kruzik 80137eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 80237eeb815SJakub Kruzik */ 80337eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 80437eeb815SJakub Kruzik { 80537eeb815SJakub Kruzik PetscErrorCode ierr; 80637eeb815SJakub Kruzik 80737eeb815SJakub Kruzik PetscFunctionBegin; 80837eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 80937eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 81037eeb815SJakub Kruzik PetscFunctionReturn(0); 81137eeb815SJakub Kruzik } 81237eeb815SJakub Kruzik 8138513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 8148513960bSJakub Kruzik { 8158513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 8168513960bSJakub Kruzik PetscBool iascii; 8178513960bSJakub Kruzik PetscErrorCode ierr; 8188513960bSJakub Kruzik 8198513960bSJakub Kruzik PetscFunctionBegin; 8208513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8218513960bSJakub Kruzik if (iascii) { 8228513960bSJakub Kruzik //if (cg->adaptiveconv) { 8238513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 8248513960bSJakub Kruzik //} 8258513960bSJakub Kruzik if (def->correct) { 8268513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 8278513960bSJakub Kruzik } 8288513960bSJakub Kruzik if (!def->nestedlvl) { 8298513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 8308513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 8318513960bSJakub Kruzik } 8328513960bSJakub Kruzik 83322b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 83422b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 83522b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 83622b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 83722b0793eSJakub Kruzik 8388513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 8398513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 8408513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 8418513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 8428513960bSJakub Kruzik } 8438513960bSJakub Kruzik PetscFunctionReturn(0); 8448513960bSJakub Kruzik } 8458513960bSJakub Kruzik 84637eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 84737eeb815SJakub Kruzik { 848880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 84937eeb815SJakub Kruzik PetscErrorCode ierr; 85037eeb815SJakub Kruzik 85137eeb815SJakub Kruzik PetscFunctionBegin; 85237eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 853859a873cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 854859a873cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 855*8a71cb68SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr); 856*8a71cb68SJakub 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); 857880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 858880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 859880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 860880d05ceSJakub Kruzik //TODO add set function and fix manpages 861880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 86237eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 86337eeb815SJakub Kruzik PetscFunctionReturn(0); 86437eeb815SJakub Kruzik } 86537eeb815SJakub Kruzik 86637eeb815SJakub Kruzik /*MC 867e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 868e361f8b5SJakub Kruzik or to a predefined value 86937eeb815SJakub Kruzik 87037eeb815SJakub Kruzik Options Database Key: 871e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 87237eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 87337eeb815SJakub Kruzik 87437eeb815SJakub Kruzik Level: beginner 87537eeb815SJakub Kruzik 87637eeb815SJakub Kruzik Notes: 877e361f8b5SJakub Kruzik todo 87837eeb815SJakub Kruzik 87937eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 880e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 88137eeb815SJakub Kruzik M*/ 88237eeb815SJakub Kruzik 88337eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 88437eeb815SJakub Kruzik { 88537eeb815SJakub Kruzik PC_Deflation *def; 88637eeb815SJakub Kruzik PetscErrorCode ierr; 88737eeb815SJakub Kruzik 88837eeb815SJakub Kruzik PetscFunctionBegin; 88937eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 89037eeb815SJakub Kruzik pc->data = (void*)def; 89137eeb815SJakub Kruzik 892e662bc50SJakub Kruzik def->init = PETSC_FALSE; 893e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 894fcb31d99SJakub Kruzik def->correctfact = 1.0; 895e662bc50SJakub Kruzik def->reductionfact = -1; 896e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 897e662bc50SJakub Kruzik def->spacesize = 1; 898e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 899e662bc50SJakub Kruzik def->nestedlvl = 0; 900e662bc50SJakub Kruzik def->maxnestedlvl = 0; 90137eeb815SJakub Kruzik 90237eeb815SJakub Kruzik /* 90337eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 90437eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 90537eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 90637eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 90737eeb815SJakub Kruzik not needed. 90837eeb815SJakub Kruzik */ 90937eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 91037eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 911b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 912b27c8092SJakub Kruzik pc->ops->postsolve = 0; 91337eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 91437eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 91537eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 91637eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 9178513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 91837eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 91937eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 92037eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 92137eeb815SJakub Kruzik 92298d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 923859a873cSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr); 924*8a71cb68SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr); 92539417d7eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 926e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 9273cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 928e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 9294a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 930f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 93198d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 93298d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 93337eeb815SJakub Kruzik PetscFunctionReturn(0); 93437eeb815SJakub Kruzik } 93537eeb815SJakub Kruzik 936