15170378fSJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/ /* includes for fortran wrappers */ 237eeb815SJakub Kruzik 38513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 48513960bSJakub Kruzik "haar", 58513960bSJakub Kruzik "db2", 68513960bSJakub Kruzik "db4", 78513960bSJakub Kruzik "db8", 88513960bSJakub Kruzik "db16", 98513960bSJakub Kruzik "biorth22", 108513960bSJakub Kruzik "meyer", 118513960bSJakub Kruzik "aggregation", 128513960bSJakub Kruzik "user", 130a78af22SJakub Kruzik "PCDeflationSpaceType", 140a78af22SJakub Kruzik "PC_DEFLATION_SPACE_", 150a545947SLisandro Dalcin NULL 168513960bSJakub Kruzik }; 178513960bSJakub Kruzik 18a122ebaeSJakub Kruzik static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc,PetscBool flg) 19a122ebaeSJakub Kruzik { 20a122ebaeSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 21a122ebaeSJakub Kruzik 22a122ebaeSJakub Kruzik PetscFunctionBegin; 23a122ebaeSJakub Kruzik def->init = flg; 24a122ebaeSJakub Kruzik PetscFunctionReturn(0); 25a122ebaeSJakub Kruzik } 26a122ebaeSJakub Kruzik 27a122ebaeSJakub Kruzik /*@ 28a122ebaeSJakub Kruzik PCDeflationSetInitOnly - Do only initialization step. 29e26ad46dSJakub Kruzik Sets initial guess to the solution on the deflation space but does not apply 30e26ad46dSJakub Kruzik the deflation preconditioner. The additional preconditioner is still applied. 31a122ebaeSJakub Kruzik 32e26ad46dSJakub Kruzik Logically Collective 33a122ebaeSJakub Kruzik 34a122ebaeSJakub Kruzik Input Parameters: 35a122ebaeSJakub Kruzik + pc - the preconditioner context 36a122ebaeSJakub Kruzik - flg - default PETSC_FALSE 37a122ebaeSJakub Kruzik 38e26ad46dSJakub Kruzik Options Database Keys: 39e26ad46dSJakub Kruzik . -pc_deflation_init_only <false> - if true computes only the special guess 40e26ad46dSJakub Kruzik 41a122ebaeSJakub Kruzik Level: intermediate 42a122ebaeSJakub Kruzik 43a122ebaeSJakub Kruzik .seealso: PCDEFLATION 44a122ebaeSJakub Kruzik @*/ 45a122ebaeSJakub Kruzik PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg) 46a122ebaeSJakub Kruzik { 47a122ebaeSJakub Kruzik PetscFunctionBegin; 48a122ebaeSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 49a122ebaeSJakub Kruzik PetscValidLogicalCollectiveBool(pc,flg,2); 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg))); 51a122ebaeSJakub Kruzik PetscFunctionReturn(0); 52a122ebaeSJakub Kruzik } 53a122ebaeSJakub Kruzik 5493b79a5aSJakub Kruzik static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc,PetscInt current,PetscInt max) 5598d6e3deSJakub Kruzik { 5698d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 5798d6e3deSJakub Kruzik 5898d6e3deSJakub Kruzik PetscFunctionBegin; 596c93e71cSJakub Kruzik if (current) def->lvl = current; 606c93e71cSJakub Kruzik def->maxlvl = max; 6198d6e3deSJakub Kruzik PetscFunctionReturn(0); 6298d6e3deSJakub Kruzik } 6398d6e3deSJakub Kruzik 6498d6e3deSJakub Kruzik /*@ 6593b79a5aSJakub Kruzik PCDeflationSetLevels - Set the maximum level of deflation nesting. 6698d6e3deSJakub Kruzik 67e26ad46dSJakub Kruzik Logically Collective 6898d6e3deSJakub Kruzik 6998d6e3deSJakub Kruzik Input Parameters: 7098d6e3deSJakub Kruzik + pc - the preconditioner context 7198d6e3deSJakub Kruzik - max - maximum deflation level 7298d6e3deSJakub Kruzik 73e26ad46dSJakub Kruzik Options Database Keys: 74e26ad46dSJakub Kruzik . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 75e26ad46dSJakub Kruzik 7698d6e3deSJakub Kruzik Level: intermediate 7798d6e3deSJakub Kruzik 78e26ad46dSJakub Kruzik .seealso: PCDeflationSetSpaceToCompute(), PCDeflationSetSpace(), PCDEFLATION 7998d6e3deSJakub Kruzik @*/ 8093b79a5aSJakub Kruzik PetscErrorCode PCDeflationSetLevels(PC pc,PetscInt max) 8198d6e3deSJakub Kruzik { 8298d6e3deSJakub Kruzik PetscFunctionBegin; 8398d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8498d6e3deSJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetLevels_C",(PC,PetscInt,PetscInt),(pc,0,max))); 8698d6e3deSJakub Kruzik PetscFunctionReturn(0); 8798d6e3deSJakub Kruzik } 8898d6e3deSJakub Kruzik 89859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red) 90859a873cSJakub Kruzik { 91859a873cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 92859a873cSJakub Kruzik 93859a873cSJakub Kruzik PetscFunctionBegin; 94859a873cSJakub Kruzik def->reductionfact = red; 95859a873cSJakub Kruzik PetscFunctionReturn(0); 96859a873cSJakub Kruzik } 97859a873cSJakub Kruzik 98859a873cSJakub Kruzik /*@ 99e26ad46dSJakub Kruzik PCDeflationSetReductionFactor - Set reduction factor for the bottom PCTELESCOPE coarse problem solver. 100859a873cSJakub Kruzik 101e26ad46dSJakub Kruzik Logically Collective 102859a873cSJakub Kruzik 103859a873cSJakub Kruzik Input Parameters: 104859a873cSJakub Kruzik + pc - the preconditioner context 105859a873cSJakub Kruzik - red - reduction factor (or PETSC_DETERMINE) 106859a873cSJakub Kruzik 107e26ad46dSJakub Kruzik Options Database Keys: 10871f2caa7Sprj- . -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for PCTELESCOPE 109e26ad46dSJakub Kruzik 110e26ad46dSJakub Kruzik Notes: 111e26ad46dSJakub Kruzik Default is computed based on the size of the coarse problem. 112e26ad46dSJakub Kruzik 113859a873cSJakub Kruzik Level: intermediate 114859a873cSJakub Kruzik 115e26ad46dSJakub Kruzik .seealso: PCTELESCOPE, PCDEFLATION 116859a873cSJakub Kruzik @*/ 117859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red) 118859a873cSJakub Kruzik { 119859a873cSJakub Kruzik PetscFunctionBegin; 120859a873cSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121859a873cSJakub Kruzik PetscValidLogicalCollectiveInt(pc,red,2); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red))); 123859a873cSJakub Kruzik PetscFunctionReturn(0); 124859a873cSJakub Kruzik } 125859a873cSJakub Kruzik 1268a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact) 1278a71cb68SJakub Kruzik { 1288a71cb68SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1298a71cb68SJakub Kruzik 1308a71cb68SJakub Kruzik PetscFunctionBegin; 1318a71cb68SJakub Kruzik /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 1328a71cb68SJakub Kruzik def->correct = PETSC_TRUE; 1338a71cb68SJakub Kruzik def->correctfact = fact; 134e26ad46dSJakub Kruzik if (def->correct == 0.0) { 1358a71cb68SJakub Kruzik def->correct = PETSC_FALSE; 1368a71cb68SJakub Kruzik } 1378a71cb68SJakub Kruzik PetscFunctionReturn(0); 1388a71cb68SJakub Kruzik } 1398a71cb68SJakub Kruzik 1408a71cb68SJakub Kruzik /*@ 1418a71cb68SJakub Kruzik PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 142e26ad46dSJakub Kruzik The Preconditioner becomes P*M^{-1} + fact*Q. 1438a71cb68SJakub Kruzik 144e26ad46dSJakub Kruzik Logically Collective 1458a71cb68SJakub Kruzik 1468a71cb68SJakub Kruzik Input Parameters: 1478a71cb68SJakub Kruzik + pc - the preconditioner context 148e26ad46dSJakub Kruzik - fact - correction factor 149e26ad46dSJakub Kruzik 150e26ad46dSJakub Kruzik Options Database Keys: 151e26ad46dSJakub Kruzik + -pc_deflation_correction <false> - if true apply coarse problem correction 152e26ad46dSJakub Kruzik - -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 153e26ad46dSJakub Kruzik 154e26ad46dSJakub Kruzik Notes: 155e26ad46dSJakub Kruzik Any non-zero fact enables the coarse problem correction. 1568a71cb68SJakub Kruzik 1578a71cb68SJakub Kruzik Level: intermediate 1588a71cb68SJakub Kruzik 1598a71cb68SJakub Kruzik .seealso: PCDEFLATION 1608a71cb68SJakub Kruzik @*/ 1618a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact) 1628a71cb68SJakub Kruzik { 1638a71cb68SJakub Kruzik PetscFunctionBegin; 1648a71cb68SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1658a71cb68SJakub Kruzik PetscValidLogicalCollectiveScalar(pc,fact,2); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact))); 1678a71cb68SJakub Kruzik PetscFunctionReturn(0); 1688a71cb68SJakub Kruzik } 1698a71cb68SJakub Kruzik 17039417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 17139417d7eSJakub Kruzik { 17239417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 17339417d7eSJakub Kruzik 17439417d7eSJakub Kruzik PetscFunctionBegin; 17539417d7eSJakub Kruzik if (type) def->spacetype = type; 17639417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 17739417d7eSJakub Kruzik PetscFunctionReturn(0); 17839417d7eSJakub Kruzik } 17939417d7eSJakub Kruzik 18039417d7eSJakub Kruzik /*@ 18139417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 18239417d7eSJakub Kruzik 183e26ad46dSJakub Kruzik Logically Collective 18439417d7eSJakub Kruzik 18539417d7eSJakub Kruzik Input Parameters: 18639417d7eSJakub Kruzik + pc - the preconditioner context 18739417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 18839417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 18939417d7eSJakub Kruzik 190e26ad46dSJakub Kruzik Options Database Keys: 191e26ad46dSJakub Kruzik + -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 192e26ad46dSJakub Kruzik - -pc_deflation_compute_space_size <1> - size of the deflation space 193e26ad46dSJakub Kruzik 19439417d7eSJakub Kruzik Notes: 19539417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 196e26ad46dSJakub Kruzik 19734a84908Sprj- The deflation space is computed in PCSetUp(). 19839417d7eSJakub Kruzik 19939417d7eSJakub Kruzik Level: intermediate 20039417d7eSJakub Kruzik 20193b79a5aSJakub Kruzik .seealso: PCDeflationSetLevels(), PCDEFLATION 20239417d7eSJakub Kruzik @*/ 20339417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 20439417d7eSJakub Kruzik { 20539417d7eSJakub Kruzik PetscFunctionBegin; 20639417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 20739417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 20839417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size))); 21039417d7eSJakub Kruzik PetscFunctionReturn(0); 21139417d7eSJakub Kruzik } 2128513960bSJakub Kruzik 213e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 214e662bc50SJakub Kruzik { 215e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 216e662bc50SJakub Kruzik 217e662bc50SJakub Kruzik PetscFunctionBegin; 21870f9219eSJakub Kruzik /* possibly allows W' = Wt (which is valid but not tested) */ 2195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)W)); 220e662bc50SJakub Kruzik if (transpose) { 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->Wt)); 222e662bc50SJakub Kruzik def->Wt = W; 223e662bc50SJakub Kruzik } else { 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->W)); 225e662bc50SJakub Kruzik def->W = W; 226e662bc50SJakub Kruzik } 227e662bc50SJakub Kruzik PetscFunctionReturn(0); 228e662bc50SJakub Kruzik } 229e662bc50SJakub Kruzik 230e662bc50SJakub Kruzik /*@ 2317e9012c5SJakub Kruzik PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose). 232e662bc50SJakub Kruzik 233e26ad46dSJakub Kruzik Logically Collective 234e662bc50SJakub Kruzik 235e662bc50SJakub Kruzik Input Parameters: 236e662bc50SJakub Kruzik + pc - the preconditioner context 237e662bc50SJakub Kruzik . W - deflation matrix 238e26ad46dSJakub Kruzik - transpose - indicates that W is an explicit transpose of the deflation matrix 239e26ad46dSJakub Kruzik 240e26ad46dSJakub Kruzik Notes: 241e26ad46dSJakub Kruzik Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel 242e26ad46dSJakub Kruzik deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and 243e26ad46dSJakub Kruzik the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with 244e26ad46dSJakub Kruzik W1 as the deflation matrix. This repeats until the maximum level set by 24593b79a5aSJakub Kruzik PCDeflationSetLevels() is reached or there are no more matrices available. 246e26ad46dSJakub Kruzik If there are matrices left after reaching the maximum level, 247e26ad46dSJakub Kruzik they are merged into a deflation matrix ...*W{n-1}*Wn. 248e662bc50SJakub Kruzik 249e662bc50SJakub Kruzik Level: intermediate 250e662bc50SJakub Kruzik 25193b79a5aSJakub Kruzik .seealso: PCDeflationSetLevels(), PCDEFLATION 252e662bc50SJakub Kruzik @*/ 253e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 254e662bc50SJakub Kruzik { 255e662bc50SJakub Kruzik PetscFunctionBegin; 256e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 257e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 258e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose))); 260e662bc50SJakub Kruzik PetscFunctionReturn(0); 261e662bc50SJakub Kruzik } 262e662bc50SJakub Kruzik 2633cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 2643cf3a049SJakub Kruzik { 2653cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2663cf3a049SJakub Kruzik 2673cf3a049SJakub Kruzik PetscFunctionBegin; 2685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mat)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->WtA)); 2703cf3a049SJakub Kruzik def->WtA = mat; 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA)); 2723cf3a049SJakub Kruzik PetscFunctionReturn(0); 2733cf3a049SJakub Kruzik } 2743cf3a049SJakub Kruzik 2753cf3a049SJakub Kruzik /*@ 276e26ad46dSJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A). 2773cf3a049SJakub Kruzik 278e26ad46dSJakub Kruzik 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 PetscFunctionBegin; 2913cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2923cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat))); 2943cf3a049SJakub Kruzik PetscFunctionReturn(0); 2953cf3a049SJakub Kruzik } 2963cf3a049SJakub Kruzik 297e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 298e924b002SJakub Kruzik { 299e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 300e924b002SJakub Kruzik 301e924b002SJakub Kruzik PetscFunctionBegin; 3025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mat)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->WtAW)); 304e924b002SJakub Kruzik def->WtAW = mat; 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW)); 306e924b002SJakub Kruzik PetscFunctionReturn(0); 307e924b002SJakub Kruzik } 308e924b002SJakub Kruzik 309e924b002SJakub Kruzik /*@ 310e26ad46dSJakub Kruzik PCDeflationSetCoarseMat - Set the coarse problem Mat. 311e924b002SJakub Kruzik 312e26ad46dSJakub Kruzik Collective 313e924b002SJakub Kruzik 314e924b002SJakub Kruzik Input Parameters: 315e924b002SJakub Kruzik + pc - preconditioner context 316e924b002SJakub Kruzik - mat - coarse problem mat 317e924b002SJakub Kruzik 318e924b002SJakub Kruzik Level: developer 319e924b002SJakub Kruzik 320e924b002SJakub Kruzik .seealso: PCDEFLATION 321e924b002SJakub Kruzik @*/ 322e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 323e924b002SJakub Kruzik { 324e924b002SJakub Kruzik PetscFunctionBegin; 325e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 326e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat))); 328e924b002SJakub Kruzik PetscFunctionReturn(0); 329e924b002SJakub Kruzik } 330e924b002SJakub Kruzik 33198d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 3325c0c31c5SJakub Kruzik { 3335c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3345c0c31c5SJakub Kruzik 3355c0c31c5SJakub Kruzik PetscFunctionBegin; 33698d6e3deSJakub Kruzik *ksp = def->WtAWinv; 3375c0c31c5SJakub Kruzik PetscFunctionReturn(0); 3385c0c31c5SJakub Kruzik } 3395c0c31c5SJakub Kruzik 3405c0c31c5SJakub Kruzik /*@ 3417e9012c5SJakub Kruzik PCDeflationGetCoarseKSP - Returns the coarse problem KSP. 3425c0c31c5SJakub Kruzik 34398d6e3deSJakub Kruzik Not Collective 3445c0c31c5SJakub Kruzik 3455c0c31c5SJakub Kruzik Input Parameters: 34698d6e3deSJakub Kruzik . pc - preconditioner context 3475c0c31c5SJakub Kruzik 348e26ad46dSJakub Kruzik Output Parameters: 34998d6e3deSJakub Kruzik . ksp - coarse problem KSP context 3505c0c31c5SJakub Kruzik 351e26ad46dSJakub Kruzik Level: advanced 35298d6e3deSJakub Kruzik 3538c4db21fSJakub Kruzik .seealso: PCDEFLATION 3545c0c31c5SJakub Kruzik @*/ 35598d6e3deSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3565c0c31c5SJakub Kruzik { 3575c0c31c5SJakub Kruzik PetscFunctionBegin; 35822b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 35998d6e3deSJakub Kruzik PetscValidPointer(ksp,2); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp))); 36198d6e3deSJakub Kruzik PetscFunctionReturn(0); 36298d6e3deSJakub Kruzik } 36398d6e3deSJakub Kruzik 364268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 365268c6673SJakub Kruzik { 366268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 367268c6673SJakub Kruzik 368268c6673SJakub Kruzik PetscFunctionBegin; 369268c6673SJakub Kruzik *apc = def->pc; 370268c6673SJakub Kruzik PetscFunctionReturn(0); 371268c6673SJakub Kruzik } 372268c6673SJakub Kruzik 373268c6673SJakub Kruzik /*@ 3747e9012c5SJakub Kruzik PCDeflationGetPC - Returns the additional preconditioner M^{-1}. 375268c6673SJakub Kruzik 376268c6673SJakub Kruzik Not Collective 377268c6673SJakub Kruzik 378268c6673SJakub Kruzik Input Parameters: 379268c6673SJakub Kruzik . pc - the preconditioner context 380268c6673SJakub Kruzik 381e26ad46dSJakub Kruzik Output Parameters: 382268c6673SJakub Kruzik . apc - additional preconditioner 383268c6673SJakub Kruzik 384268c6673SJakub Kruzik Level: advanced 385268c6673SJakub Kruzik 3868c4db21fSJakub Kruzik .seealso: PCDEFLATION 387268c6673SJakub Kruzik @*/ 388268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 389268c6673SJakub Kruzik { 390268c6673SJakub Kruzik PetscFunctionBegin; 391268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 392064a246eSJacob Faibussowitsch PetscValidPointer(pc,1); 3935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc))); 394268c6673SJakub Kruzik PetscFunctionReturn(0); 395268c6673SJakub Kruzik } 396268c6673SJakub Kruzik 39737eeb815SJakub Kruzik /* 398b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 399b27c8092SJakub Kruzik */ 400b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 401b27c8092SJakub Kruzik { 402b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 403b27c8092SJakub Kruzik Mat A; 404b27c8092SJakub Kruzik Vec r,w1,w2; 405b27c8092SJakub Kruzik PetscBool nonzero; 40637eeb815SJakub Kruzik 407b27c8092SJakub Kruzik PetscFunctionBegin; 408b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 409b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 410b27c8092SJakub Kruzik r = def->work; 4115f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOperators(pc,NULL,&A)); 41237eeb815SJakub Kruzik 4135f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetInitialGuessNonzero(ksp,&nonzero)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetInitialGuessNonzero(ksp,PETSC_TRUE)); 415b27c8092SJakub Kruzik if (nonzero) { 4165f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,r)); /* r <- b - Ax */ 4175f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(r,-1.0,b)); 418b27c8092SJakub Kruzik } else { 4195f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(b,r)); /* r <- b (x is 0) */ 420b27c8092SJakub Kruzik } 421b27c8092SJakub Kruzik 422a3177931SJakub Kruzik if (def->Wt) { 4235f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(def->Wt,r,w1)); /* w1 <- W'*r */ 424a3177931SJakub Kruzik } else { 4255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTranspose(def->W,r,w1)); /* w1 <- W'*r */ 426a3177931SJakub Kruzik } 4275f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(def->WtAWinv,w1,w2)); /* w2 <- (W'*A*W)^{-1}*w1 */ 4285f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(def->W,w2,r)); /* r <- W*w2 */ 4295f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(x,1.0,r)); 430b27c8092SJakub Kruzik PetscFunctionReturn(0); 431b27c8092SJakub Kruzik } 43237eeb815SJakub Kruzik 433f8bf229cSJakub Kruzik /* 434f8bf229cSJakub Kruzik if (def->correct) { 435ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 436f8bf229cSJakub Kruzik } else { 437ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 438f8bf229cSJakub Kruzik } 43937eeb815SJakub Kruzik */ 440f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 441f8bf229cSJakub Kruzik { 442f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 443f8bf229cSJakub Kruzik Mat A; 444f8bf229cSJakub Kruzik Vec u,w1,w2; 445f8bf229cSJakub Kruzik 446f8bf229cSJakub Kruzik PetscFunctionBegin; 447f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 448f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 449f8bf229cSJakub Kruzik u = def->work; 4505f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOperators(pc,NULL,&A)); 451f8bf229cSJakub Kruzik 4525f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(def->pc,r,z)); /* z <- M^{-1}*r */ 45322b0793eSJakub Kruzik if (!def->init) { 4545f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(def->WtA,z,w1)); /* w1 <- W'*A*z */ 455f8bf229cSJakub Kruzik if (def->correct) { 456ae029463SJakub Kruzik if (def->Wt) { 4575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(def->Wt,r,w2)); /* w2 <- W'*r */ 458ae029463SJakub Kruzik } else { 4595f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTranspose(def->W,r,w2)); /* w2 <- W'*r */ 460f8bf229cSJakub Kruzik } 4615f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(w1,-1.0*def->correctfact,w2)); /* w1 <- w1 - l*w2 */ 462f8bf229cSJakub Kruzik } 4635f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(def->WtAWinv,w1,w2)); /* w2 <- (W'*A*W)^{-1}*w1 */ 4645f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(def->W,w2,u)); /* u <- W*w2 */ 4655f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,-1.0,u)); /* z <- z - u */ 46622b0793eSJakub Kruzik } 467f8bf229cSJakub Kruzik PetscFunctionReturn(0); 468f8bf229cSJakub Kruzik } 469f8bf229cSJakub Kruzik 47037eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 47137eeb815SJakub Kruzik { 472409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 473409a357bSJakub Kruzik KSP innerksp; 474409a357bSJakub Kruzik PC pcinner; 475409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 4767b3faf33SJakub Kruzik PetscInt i,m,red,size; 4777b3faf33SJakub Kruzik PetscMPIInt commsize; 478409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 479409a357bSJakub Kruzik MatCompositeType ctype; 480409a357bSJakub Kruzik MPI_Comm comm; 4816c93e71cSJakub Kruzik char prefix[128]=""; 48237eeb815SJakub Kruzik 48337eeb815SJakub Kruzik PetscFunctionBegin; 484409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)pc,&comm)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOperators(pc,NULL,&Amat)); 4876c93e71cSJakub Kruzik if (!def->lvl && !def->prefix) { 4885f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc,&def->prefix)); 4896c93e71cSJakub Kruzik } 4906c93e71cSJakub Kruzik if (def->lvl) { 4915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(prefix,sizeof(prefix),"%d_",(int)def->lvl)); 4926c93e71cSJakub Kruzik } 493f8bf229cSJakub Kruzik 494f8bf229cSJakub Kruzik /* compute a deflation space */ 495409a357bSJakub Kruzik if (def->W || def->Wt) { 496409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 497409a357bSJakub Kruzik } else { 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PCDeflationComputeSpace(pc)); 49937eeb815SJakub Kruzik } 50037eeb815SJakub Kruzik 501409a357bSJakub Kruzik /* nested deflation */ 502409a357bSJakub Kruzik if (def->W) { 5035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match)); 504409a357bSJakub Kruzik if (match) { 5055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetType(def->W,&ctype)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetNumberMat(def->W,&size)); 50737eeb815SJakub Kruzik } 508409a357bSJakub Kruzik } else { 5095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHermitianTranspose(def->Wt,&def->W)); 5105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match)); 511409a357bSJakub Kruzik if (match) { 5125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetType(def->Wt,&ctype)); 5135f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetNumberMat(def->Wt,&size)); 514409a357bSJakub Kruzik } 515409a357bSJakub Kruzik transp = PETSC_TRUE; 516409a357bSJakub Kruzik } 517409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 518409a357bSJakub Kruzik if (!transp) { 5196c93e71cSJakub Kruzik if (def->lvl < def->maxlvl) { 5205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&mats)); 521409a357bSJakub Kruzik for (i=0; i<size; i++) { 5225f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetMat(def->W,i,&mats[i])); 523409a357bSJakub Kruzik } 524409a357bSJakub Kruzik size -= 1; 5255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->W)); 526409a357bSJakub Kruzik def->W = mats[size]; 5275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mats[size])); 528409a357bSJakub Kruzik if (size > 1) { 5295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(comm,size,mats,&nextDef)); 5305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE)); 531409a357bSJakub Kruzik } else { 532409a357bSJakub Kruzik nextDef = mats[0]; 5335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mats[0])); 534409a357bSJakub Kruzik } 5355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mats)); 536409a357bSJakub Kruzik } else { 5371fdb00f9SJakub Kruzik /* TODO test merge side performance */ 5385f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */ 5395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeMerge(def->W)); 540409a357bSJakub Kruzik } 54128da0170SJakub Kruzik } else { 5426c93e71cSJakub Kruzik if (def->lvl < def->maxlvl) { 5435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&mats)); 54428da0170SJakub Kruzik for (i=0; i<size; i++) { 5455f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetMat(def->Wt,i,&mats[i])); 54628da0170SJakub Kruzik } 54728da0170SJakub Kruzik size -= 1; 5485f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->Wt)); 54928da0170SJakub Kruzik def->Wt = mats[0]; 5505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mats[0])); 55128da0170SJakub Kruzik if (size > 1) { 5525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(comm,size,&mats[1],&nextDef)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE)); 55428da0170SJakub Kruzik } else { 55528da0170SJakub Kruzik nextDef = mats[1]; 5565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mats[1])); 557409a357bSJakub Kruzik } 5585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mats)); 55928da0170SJakub Kruzik } else { 5605f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */ 5615f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeMerge(def->Wt)); 56228da0170SJakub Kruzik } 56328da0170SJakub Kruzik } 56428da0170SJakub Kruzik } 56528da0170SJakub Kruzik 56628da0170SJakub Kruzik if (transp) { 5675f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->W)); 5685f80ce2aSJacob Faibussowitsch CHKERRQ(MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W)); 569409a357bSJakub Kruzik } 570409a357bSJakub Kruzik 571ae029463SJakub Kruzik /* assemble WtA */ 572ae029463SJakub Kruzik if (!def->WtA) { 573ae029463SJakub Kruzik if (def->Wt) { 5745f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA)); 575ae029463SJakub Kruzik } else { 576a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 5775f80ce2aSJacob Faibussowitsch CHKERRQ(MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt)); 5785f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA)); 579a3177931SJakub Kruzik #else 5805f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA)); 581a3177931SJakub Kruzik #endif 582ae029463SJakub Kruzik } 583ae029463SJakub Kruzik } 584409a357bSJakub Kruzik /* setup coarse problem */ 585409a357bSJakub Kruzik if (!def->WtAWinv) { 5865f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(def->W,NULL,&m)); 587409a357bSJakub Kruzik if (!def->WtAW) { 5885f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW)); 589409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 5905f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOption(Amat,MAT_SPD,&flgspd)); 5915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(def->WtAW,MAT_SPD,flgspd)); 59276bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 593ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 594409a357bSJakub Kruzik PetscReal *norms; 5955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m,&norms)); 5965f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms)); 597409a357bSJakub Kruzik for (i=0; i<m; i++) { 598409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 59998921bdaSJacob Faibussowitsch SETERRQ(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 600409a357bSJakub Kruzik } 601409a357bSJakub Kruzik } 6025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(norms)); 603377006c5SJakub Kruzik } 604409a357bSJakub Kruzik } else { 6055f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOption(def->WtAW,MAT_SPD,&flgspd)); 606409a357bSJakub Kruzik } 6071fdb00f9SJakub Kruzik /* TODO use MATINV ? */ 6085f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(comm,&def->WtAWinv)); 6095f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW)); 6105f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(def->WtAWinv,&pcinner)); 611557a2f7dSJakub Kruzik /* Setup KSP and PC */ 612557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 613557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 614557a2f7dSJakub Kruzik /* set default KSPtype */ 615557a2f7dSJakub Kruzik if (!def->ksptype) { 616557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 617557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 618557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 619557a2f7dSJakub Kruzik } 620557a2f7dSJakub Kruzik } 6215f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(innerksp,def->ksptype)); /* TODO iherit from KSP + tolerances */ 6225f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pcinner,PCDEFLATION)); /* TODO create coarse preconditinoner M_c = WtMW ? */ 6235f80ce2aSJacob Faibussowitsch CHKERRQ(PCDeflationSetSpace(pcinner,nextDef,transp)); 6245f80ce2aSJacob Faibussowitsch CHKERRQ(PCDeflationSetLevels_Deflation(pcinner,def->lvl+1,def->maxlvl)); 625ae029463SJakub Kruzik /* inherit options */ 6266c93e71cSJakub Kruzik if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix; 6279f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->init = def->init; 628557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 629557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 6309f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->correctfact = def->correctfact; 6319f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact; 6325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&nextDef)); 633557a2f7dSJakub Kruzik } else { /* the last level */ 6345f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(def->WtAWinv,KSPPREONLY)); 6355f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pcinner,PCTELESCOPE)); 6366c93e71cSJakub Kruzik /* do not overwrite PCTELESCOPE */ 6376c93e71cSJakub Kruzik if (def->prefix) { 6385f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(def->WtAWinv,def->prefix)); 639ac66f006SJakub Kruzik } 6405f80ce2aSJacob Faibussowitsch CHKERRQ(KSPAppendOptionsPrefix(def->WtAWinv,"deflation_tel_")); 6415f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions(pcinner)); 6425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pcinner,PCTELESCOPE,&match)); 643*28b400f6SJacob Faibussowitsch PetscCheck(match,comm,PETSC_ERR_SUP,"User can not owerwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead."); 644409a357bSJakub Kruzik /* Reduction factor choice */ 645409a357bSJakub Kruzik red = def->reductionfact; 646409a357bSJakub Kruzik if (red < 0) { 6475f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&commsize)); 648faa75363SBarry Smith red = PetscCeilInt(commsize,PetscCeilInt(m,commsize)); 6495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"")); 650409a357bSJakub Kruzik if (match) red = commsize; 6515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(pc,"Auto choosing reduction factor %D\n",red)); 652409a357bSJakub Kruzik } 6535f80ce2aSJacob Faibussowitsch CHKERRQ(PCTelescopeSetReductionFactor(pcinner,red)); 6545f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetUp(pcinner)); 6555f80ce2aSJacob Faibussowitsch CHKERRQ(PCTelescopeGetKSP(pcinner,&innerksp)); 656ac66f006SJakub Kruzik if (innerksp) { 6575f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(innerksp,&pcinner)); 6585f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pcinner,PCLU)); 659481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU) 6605f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU,MAT_FACTOR_LU,&match)); 661481b7641SJakub Kruzik if (match) { 6625f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU)); 663481b7641SJakub Kruzik } 664481b7641SJakub Kruzik #endif 665481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU_DIST) 6665f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&match)); 667481b7641SJakub Kruzik if (match) { 6685f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST)); 669409a357bSJakub Kruzik } 670481b7641SJakub Kruzik #endif 671409a357bSJakub Kruzik } 672557a2f7dSJakub Kruzik } 673557a2f7dSJakub Kruzik 674557a2f7dSJakub Kruzik if (innerksp) { 6756c93e71cSJakub Kruzik if (def->prefix) { 6765f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(innerksp,def->prefix)); 6775f80ce2aSJacob Faibussowitsch CHKERRQ(KSPAppendOptionsPrefix(innerksp,"deflation_")); 6786c93e71cSJakub Kruzik } else { 6795f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(innerksp,"deflation_")); 6806c93e71cSJakub Kruzik } 6815f80ce2aSJacob Faibussowitsch CHKERRQ(KSPAppendOptionsPrefix(innerksp,prefix)); 6825f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(innerksp)); 6835f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetUp(innerksp)); 684ac66f006SJakub Kruzik } 685409a357bSJakub Kruzik } 6865f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(def->WtAWinv)); 6875f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetUp(def->WtAWinv)); 688409a357bSJakub Kruzik 6891fdb00f9SJakub Kruzik /* create preconditioner */ 69022b0793eSJakub Kruzik if (!def->pc) { 6915f80ce2aSJacob Faibussowitsch CHKERRQ(PCCreate(comm,&def->pc)); 6925f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOperators(def->pc,Amat,Amat)); 6935f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(def->pc,PCNONE)); 6946c93e71cSJakub Kruzik if (def->prefix) { 6955f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOptionsPrefix(def->pc,def->prefix)); 69622b0793eSJakub Kruzik } 6975f80ce2aSJacob Faibussowitsch CHKERRQ(PCAppendOptionsPrefix(def->pc,"deflation_")); 6985f80ce2aSJacob Faibussowitsch CHKERRQ(PCAppendOptionsPrefix(def->pc,prefix)); 6995f80ce2aSJacob Faibussowitsch CHKERRQ(PCAppendOptionsPrefix(def->pc,"pc_")); 7005f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions(def->pc)); 7015f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetUp(def->pc)); 70222b0793eSJakub Kruzik } 70322b0793eSJakub Kruzik 704f8bf229cSJakub Kruzik /* create work vecs */ 7055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(Amat,NULL,&def->work)); 7065f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL)); 70737eeb815SJakub Kruzik PetscFunctionReturn(0); 70837eeb815SJakub Kruzik } 709b30d4775SJakub Kruzik 71037eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 71137eeb815SJakub Kruzik { 712b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 71337eeb815SJakub Kruzik 71437eeb815SJakub Kruzik PetscFunctionBegin; 7155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&def->work)); 7165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroyVecs(2,&def->workcoarse)); 7175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->W)); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->Wt)); 7195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->WtA)); 7205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&def->WtAW)); 7215f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&def->WtAWinv)); 7225f80ce2aSJacob Faibussowitsch CHKERRQ(PCDestroy(&def->pc)); 72337eeb815SJakub Kruzik PetscFunctionReturn(0); 72437eeb815SJakub Kruzik } 72537eeb815SJakub Kruzik 72637eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 72737eeb815SJakub Kruzik { 72837eeb815SJakub Kruzik PetscFunctionBegin; 7295f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_Deflation(pc)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc->data)); 73137eeb815SJakub Kruzik PetscFunctionReturn(0); 73237eeb815SJakub Kruzik } 73337eeb815SJakub Kruzik 7348513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 7358513960bSJakub Kruzik { 7368513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 7371ac6250aSJakub Kruzik PetscInt its; 7388513960bSJakub Kruzik PetscBool iascii; 7398513960bSJakub Kruzik PetscErrorCode ierr; 7408513960bSJakub Kruzik 7418513960bSJakub Kruzik PetscFunctionBegin; 7425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 7438513960bSJakub Kruzik if (iascii) { 7448513960bSJakub Kruzik if (def->correct) { 7451ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n", 7461ac6250aSJakub Kruzik (double)PetscRealPart(def->correctfact), 7471ac6250aSJakub Kruzik (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr); 7488513960bSJakub Kruzik } 7496c93e71cSJakub Kruzik if (!def->lvl) { 7505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype])); 7518513960bSJakub Kruzik } 7528513960bSJakub Kruzik 7535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n")); 7545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 7555f80ce2aSJacob Faibussowitsch CHKERRQ(PCView(def->pc,viewer)); 7565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 75722b0793eSJakub Kruzik 7585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n")); 7595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 7605f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetTotalIterations(def->WtAWinv,&its)); 7615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its)); 7625f80ce2aSJacob Faibussowitsch CHKERRQ(KSPView(def->WtAWinv,viewer)); 7635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 7648513960bSJakub Kruzik } 7658513960bSJakub Kruzik PetscFunctionReturn(0); 7668513960bSJakub Kruzik } 7678513960bSJakub Kruzik 76837eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 76937eeb815SJakub Kruzik { 770880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 77137eeb815SJakub Kruzik 77237eeb815SJakub Kruzik PetscFunctionBegin; 7735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Deflation options")); 7745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL)); 7755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_deflation_levels","Maximum of deflation levels","PCDeflationSetLevels",def->maxlvl,&def->maxlvl,NULL)); 7765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL)); 7775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL)); 7785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL)); 7795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL)); 7805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL)); 7815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL)); 7825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 78337eeb815SJakub Kruzik PetscFunctionReturn(0); 78437eeb815SJakub Kruzik } 78537eeb815SJakub Kruzik 78637eeb815SJakub Kruzik /*MC 787e26ad46dSJakub Kruzik PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value. 78837eeb815SJakub Kruzik 789e26ad46dSJakub Kruzik Options Database Keys: 790e26ad46dSJakub Kruzik + -pc_deflation_init_only <false> - if true computes only the special guess 791e26ad46dSJakub Kruzik . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 79271f2caa7Sprj- . -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem) 793e26ad46dSJakub Kruzik . -pc_deflation_correction <false> - if true apply coarse problem correction 794e26ad46dSJakub Kruzik . -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 795e26ad46dSJakub Kruzik . -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 796e26ad46dSJakub Kruzik - -pc_deflation_compute_space_size <1> - size of the deflation space (corresponds to number of levels for wavelet-based deflation) 79737eeb815SJakub Kruzik 79837eeb815SJakub Kruzik Notes: 7997e9012c5SJakub Kruzik Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2]) 800e26ad46dSJakub Kruzik preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat. 801e26ad46dSJakub Kruzik 802e26ad46dSJakub Kruzik The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space. 803e26ad46dSJakub Kruzik If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the 804e26ad46dSJakub Kruzik preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner 8051fdb00f9SJakub Kruzik application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues 806e26ad46dSJakub Kruzik to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some 807e26ad46dSJakub Kruzik eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates 808e26ad46dSJakub Kruzik of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper. 809e26ad46dSJakub Kruzik 810e26ad46dSJakub Kruzik The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can 811e26ad46dSJakub Kruzik be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size. 812e26ad46dSJakub Kruzik User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix 813e26ad46dSJakub Kruzik is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the 8141fdb00f9SJakub Kruzik deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned 815e26ad46dSJakub Kruzik by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum 8161fdb00f9SJakub Kruzik level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged 8171fdb00f9SJakub Kruzik (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by 81893b79a5aSJakub Kruzik PCDeflationSetLevels() or by -pc_deflation_levels. 819e26ad46dSJakub Kruzik 820e63c2c6dSJakub Kruzik The coarse problem KSP can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1] 8216c93e71cSJakub Kruzik from the second level onward. You can also use 8228c4db21fSJakub Kruzik PCDeflationGetCoarseKSP() to control it from code. The bottom level KSP defaults to 8231fdb00f9SJakub Kruzik KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE. 824481b7641SJakub Kruzik For convenience, the reduction factor can be set by PCDeflationSetReductionFactor() 825481b7641SJakub Kruzik or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size. 826e26ad46dSJakub Kruzik 827e63c2c6dSJakub Kruzik The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for 828e63c2c6dSJakub Kruzik coarse problem KSP apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use 8298c4db21fSJakub Kruzik PCDeflationGetPC() to control the additional preconditioner from code. It defaults to PCNONE. 830e26ad46dSJakub Kruzik 831e26ad46dSJakub Kruzik The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can 832e26ad46dSJakub Kruzik be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can 833e26ad46dSJakub Kruzik significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We 834e26ad46dSJakub Kruzik recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create 835e26ad46dSJakub Kruzik an isolated eigenvalue. 836e26ad46dSJakub Kruzik 8371fdb00f9SJakub Kruzik The options are automatically inherited from the previous deflation level. 8389f604af8SJakub Kruzik 839e26ad46dSJakub Kruzik The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also 8401fdb00f9SJakub Kruzik recommend limiting the number of iterations for the coarse problems. 841e26ad46dSJakub Kruzik 842e26ad46dSJakub Kruzik See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients. 843e26ad46dSJakub Kruzik Section 4 describes some possible choices for the deflation space. 844e26ad46dSJakub Kruzik 845e26ad46dSJakub Kruzik Developer Notes: 846e26ad46dSJakub Kruzik Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech 847e26ad46dSJakub Kruzik Academy of Sciences and VSB - TU Ostrava. 848e26ad46dSJakub Kruzik 849e26ad46dSJakub Kruzik Developed from PERMON code used in [4] while on a research stay with 850e26ad46dSJakub Kruzik Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin. 851e26ad46dSJakub Kruzik 852e26ad46dSJakub Kruzik References: 853606c0280SSatish Balay + * - A. Nicolaides. "Deflation of conjugate gradients with applications to boundary value problems", SIAM J. Numer. Anal. 24.2, 1987. 854606c0280SSatish Balay . * - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988. 855606c0280SSatish Balay . * - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008. 856606c0280SSatish Balay - * - J. Kruzik "Implementation of the Deflated Variants of the Conjugate Gradient Method", Master's thesis, VSB-TUO, 2018 - http://dspace5.vsb.cz/bitstream/handle/10084/130303/KRU0097_USP_N2658_2612T078_2018.pdf 857e26ad46dSJakub Kruzik 858e26ad46dSJakub Kruzik Level: intermediate 85937eeb815SJakub Kruzik 86037eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 86193b79a5aSJakub Kruzik PCDeflationSetInitOnly(), PCDeflationSetLevels(), PCDeflationSetReductionFactor(), 8621fdb00f9SJakub Kruzik PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompute(), 863e26ad46dSJakub Kruzik PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(), 8648c4db21fSJakub Kruzik PCDeflationSetCoarseMat(), PCDeflationGetCoarseKSP(), PCDeflationGetPC() 86537eeb815SJakub Kruzik M*/ 86637eeb815SJakub Kruzik 86737eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 86837eeb815SJakub Kruzik { 86937eeb815SJakub Kruzik PC_Deflation *def; 87037eeb815SJakub Kruzik 87137eeb815SJakub Kruzik PetscFunctionBegin; 8725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&def)); 87337eeb815SJakub Kruzik pc->data = (void*)def; 87437eeb815SJakub Kruzik 875e662bc50SJakub Kruzik def->init = PETSC_FALSE; 876e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 877fcb31d99SJakub Kruzik def->correctfact = 1.0; 878e662bc50SJakub Kruzik def->reductionfact = -1; 879e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 880e662bc50SJakub Kruzik def->spacesize = 1; 881e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 8826c93e71cSJakub Kruzik def->lvl = 0; 8836c93e71cSJakub Kruzik def->maxlvl = 0; 88470f9219eSJakub Kruzik def->W = NULL; 88570f9219eSJakub Kruzik def->Wt = NULL; 88637eeb815SJakub Kruzik 88737eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 888b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 88937eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 89037eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 89137eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 89237eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 8938513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 89437eeb815SJakub Kruzik 8955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation)); 8965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLevels_C",PCDeflationSetLevels_Deflation)); 8975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation)); 8985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation)); 8995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation)); 9005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation)); 9015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation)); 9025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation)); 9035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation)); 9045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation)); 90537eeb815SJakub Kruzik PetscFunctionReturn(0); 90637eeb815SJakub Kruzik } 907