1292e2e67SJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 237eeb815SJakub Kruzik 38513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 48513960bSJakub Kruzik "haar", 58513960bSJakub Kruzik "jacket-haar", 68513960bSJakub Kruzik "db2", 78513960bSJakub Kruzik "db4", 88513960bSJakub Kruzik "db8", 98513960bSJakub Kruzik "db16", 108513960bSJakub Kruzik "biorth22", 118513960bSJakub Kruzik "meyer", 128513960bSJakub Kruzik "aggregation", 138513960bSJakub Kruzik "slepc", 148513960bSJakub Kruzik "slepc-cheap", 158513960bSJakub Kruzik "user", 160a78af22SJakub Kruzik "PCDeflationSpaceType", 170a78af22SJakub Kruzik "PC_DEFLATION_SPACE_", 188513960bSJakub Kruzik 0 198513960bSJakub Kruzik }; 208513960bSJakub Kruzik 21a122ebaeSJakub Kruzik static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc,PetscBool flg) 22a122ebaeSJakub Kruzik { 23a122ebaeSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 24a122ebaeSJakub Kruzik 25a122ebaeSJakub Kruzik PetscFunctionBegin; 26a122ebaeSJakub Kruzik def->init = flg; 27a122ebaeSJakub Kruzik PetscFunctionReturn(0); 28a122ebaeSJakub Kruzik } 29a122ebaeSJakub Kruzik 30a122ebaeSJakub Kruzik /*@ 31a122ebaeSJakub Kruzik PCDeflationSetInitOnly - Do only initialization step. 32e26ad46dSJakub Kruzik Sets initial guess to the solution on the deflation space but does not apply 33e26ad46dSJakub Kruzik the deflation preconditioner. The additional preconditioner is still applied. 34a122ebaeSJakub Kruzik 35e26ad46dSJakub Kruzik Logically Collective 36a122ebaeSJakub Kruzik 37a122ebaeSJakub Kruzik Input Parameters: 38a122ebaeSJakub Kruzik + pc - the preconditioner context 39a122ebaeSJakub Kruzik - flg - default PETSC_FALSE 40a122ebaeSJakub Kruzik 41e26ad46dSJakub Kruzik Options Database Keys: 42e26ad46dSJakub Kruzik . -pc_deflation_init_only <false> - if true computes only the special guess 43e26ad46dSJakub Kruzik 44a122ebaeSJakub Kruzik Level: intermediate 45a122ebaeSJakub Kruzik 46a122ebaeSJakub Kruzik .seealso: PCDEFLATION 47a122ebaeSJakub Kruzik @*/ 48a122ebaeSJakub Kruzik PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg) 49a122ebaeSJakub Kruzik { 50a122ebaeSJakub Kruzik PetscErrorCode ierr; 51a122ebaeSJakub Kruzik 52a122ebaeSJakub Kruzik PetscFunctionBegin; 53a122ebaeSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 54a122ebaeSJakub Kruzik PetscValidLogicalCollectiveBool(pc,flg,2); 55a122ebaeSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 56a122ebaeSJakub Kruzik PetscFunctionReturn(0); 57a122ebaeSJakub Kruzik } 58a122ebaeSJakub Kruzik 59a122ebaeSJakub Kruzik 6098d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 6198d6e3deSJakub Kruzik { 6298d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 6398d6e3deSJakub Kruzik 6498d6e3deSJakub Kruzik PetscFunctionBegin; 656c93e71cSJakub Kruzik if (current) def->lvl = current; 666c93e71cSJakub Kruzik def->maxlvl = max; 6798d6e3deSJakub Kruzik PetscFunctionReturn(0); 6898d6e3deSJakub Kruzik } 6998d6e3deSJakub Kruzik 7098d6e3deSJakub Kruzik /*@ 71e26ad46dSJakub Kruzik PCDeflationSetMaxLvl - Set the maximum level of deflation nesting. 7298d6e3deSJakub Kruzik 73e26ad46dSJakub Kruzik Logically Collective 7498d6e3deSJakub Kruzik 7598d6e3deSJakub Kruzik Input Parameters: 7698d6e3deSJakub Kruzik + pc - the preconditioner context 7798d6e3deSJakub Kruzik - max - maximum deflation level 7898d6e3deSJakub Kruzik 79e26ad46dSJakub Kruzik Options Database Keys: 80e26ad46dSJakub Kruzik . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 81e26ad46dSJakub Kruzik 8298d6e3deSJakub Kruzik Level: intermediate 8398d6e3deSJakub Kruzik 84e26ad46dSJakub Kruzik .seealso: PCDeflationSetSpaceToCompute(), PCDeflationSetSpace(), PCDEFLATION 8598d6e3deSJakub Kruzik @*/ 8698d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 8798d6e3deSJakub Kruzik { 8898d6e3deSJakub Kruzik PetscErrorCode ierr; 8998d6e3deSJakub Kruzik 9098d6e3deSJakub Kruzik PetscFunctionBegin; 9198d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9298d6e3deSJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 9398d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 9498d6e3deSJakub Kruzik PetscFunctionReturn(0); 9598d6e3deSJakub Kruzik } 9698d6e3deSJakub Kruzik 97859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red) 98859a873cSJakub Kruzik { 99859a873cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 100859a873cSJakub Kruzik 101859a873cSJakub Kruzik PetscFunctionBegin; 102859a873cSJakub Kruzik def->reductionfact = red; 103859a873cSJakub Kruzik PetscFunctionReturn(0); 104859a873cSJakub Kruzik } 105859a873cSJakub Kruzik 106859a873cSJakub Kruzik /*@ 107e26ad46dSJakub Kruzik PCDeflationSetReductionFactor - Set reduction factor for the bottom PCTELESCOPE coarse problem solver. 108859a873cSJakub Kruzik 109e26ad46dSJakub Kruzik Logically Collective 110859a873cSJakub Kruzik 111859a873cSJakub Kruzik Input Parameters: 112859a873cSJakub Kruzik + pc - the preconditioner context 113859a873cSJakub Kruzik - red - reduction factor (or PETSC_DETERMINE) 114859a873cSJakub Kruzik 115e26ad46dSJakub Kruzik Options Database Keys: 116e26ad46dSJakub Kruzik . -pc_deflation_reduction_factor <-1> - reduction factor on bottom level coarse problem for PCTELESCOPE 117e26ad46dSJakub Kruzik 118e26ad46dSJakub Kruzik Notes: 119e26ad46dSJakub Kruzik Default is computed based on the size of the coarse problem. 120e26ad46dSJakub Kruzik 121859a873cSJakub Kruzik Level: intermediate 122859a873cSJakub Kruzik 123e26ad46dSJakub Kruzik .seealso: PCTELESCOPE, PCDEFLATION 124859a873cSJakub Kruzik @*/ 125859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red) 126859a873cSJakub Kruzik { 127859a873cSJakub Kruzik PetscErrorCode ierr; 128859a873cSJakub Kruzik 129859a873cSJakub Kruzik PetscFunctionBegin; 130859a873cSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 131859a873cSJakub Kruzik PetscValidLogicalCollectiveInt(pc,red,2); 132859a873cSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr); 133859a873cSJakub Kruzik PetscFunctionReturn(0); 134859a873cSJakub Kruzik } 135859a873cSJakub Kruzik 1368a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact) 1378a71cb68SJakub Kruzik { 1388a71cb68SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1398a71cb68SJakub Kruzik 1408a71cb68SJakub Kruzik PetscFunctionBegin; 1418a71cb68SJakub Kruzik /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */ 1428a71cb68SJakub Kruzik def->correct = PETSC_TRUE; 1438a71cb68SJakub Kruzik def->correctfact = fact; 144e26ad46dSJakub Kruzik if (def->correct == 0.0) { 1458a71cb68SJakub Kruzik def->correct = PETSC_FALSE; 1468a71cb68SJakub Kruzik } 1478a71cb68SJakub Kruzik PetscFunctionReturn(0); 1488a71cb68SJakub Kruzik } 1498a71cb68SJakub Kruzik 1508a71cb68SJakub Kruzik /*@ 1518a71cb68SJakub Kruzik PCDeflationSetCorrectionFactor - Set coarse problem correction factor. 152e26ad46dSJakub Kruzik The Preconditioner becomes P*M^{-1} + fact*Q. 1538a71cb68SJakub Kruzik 154e26ad46dSJakub Kruzik Logically Collective 1558a71cb68SJakub Kruzik 1568a71cb68SJakub Kruzik Input Parameters: 1578a71cb68SJakub Kruzik + pc - the preconditioner context 158e26ad46dSJakub Kruzik - fact - correction factor 159e26ad46dSJakub Kruzik 160e26ad46dSJakub Kruzik Options Database Keys: 161e26ad46dSJakub Kruzik + -pc_deflation_correction <false> - if true apply coarse problem correction 162e26ad46dSJakub Kruzik - -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 163e26ad46dSJakub Kruzik 164e26ad46dSJakub Kruzik Notes: 165e26ad46dSJakub Kruzik Any non-zero fact enables the coarse problem correction. 1668a71cb68SJakub Kruzik 1678a71cb68SJakub Kruzik Level: intermediate 1688a71cb68SJakub Kruzik 1698a71cb68SJakub Kruzik .seealso: PCDEFLATION 1708a71cb68SJakub Kruzik @*/ 1718a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact) 1728a71cb68SJakub Kruzik { 1738a71cb68SJakub Kruzik PetscErrorCode ierr; 1748a71cb68SJakub Kruzik 1758a71cb68SJakub Kruzik PetscFunctionBegin; 1768a71cb68SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1778a71cb68SJakub Kruzik PetscValidLogicalCollectiveScalar(pc,fact,2); 1788a71cb68SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr); 1798a71cb68SJakub Kruzik PetscFunctionReturn(0); 1808a71cb68SJakub Kruzik } 1818a71cb68SJakub Kruzik 18239417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 18339417d7eSJakub Kruzik { 18439417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 18539417d7eSJakub Kruzik 18639417d7eSJakub Kruzik PetscFunctionBegin; 18739417d7eSJakub Kruzik if (type) def->spacetype = type; 18839417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 18939417d7eSJakub Kruzik PetscFunctionReturn(0); 19039417d7eSJakub Kruzik } 19139417d7eSJakub Kruzik 19239417d7eSJakub Kruzik /*@ 19339417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 19439417d7eSJakub Kruzik 195e26ad46dSJakub Kruzik Logically Collective 19639417d7eSJakub Kruzik 19739417d7eSJakub Kruzik Input Parameters: 19839417d7eSJakub Kruzik + pc - the preconditioner context 19939417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 20039417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 20139417d7eSJakub Kruzik 202e26ad46dSJakub Kruzik Options Database Keys: 203e26ad46dSJakub Kruzik + -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 204e26ad46dSJakub Kruzik - -pc_deflation_compute_space_size <1> - size of the deflation space 205e26ad46dSJakub Kruzik 20639417d7eSJakub Kruzik Notes: 20739417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 208e26ad46dSJakub Kruzik 20939417d7eSJakub Kruzik The deflation space is computed in PCSetUP(). 21039417d7eSJakub Kruzik 21139417d7eSJakub Kruzik Level: intermediate 21239417d7eSJakub Kruzik 213e26ad46dSJakub Kruzik .seealso: PCDeflationSetMaxLvl(), PCDEFLATION 21439417d7eSJakub Kruzik @*/ 21539417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 21639417d7eSJakub Kruzik { 21739417d7eSJakub Kruzik PetscErrorCode ierr; 21839417d7eSJakub Kruzik 21939417d7eSJakub Kruzik PetscFunctionBegin; 22039417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 22139417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 22239417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 22339417d7eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 22439417d7eSJakub Kruzik PetscFunctionReturn(0); 22539417d7eSJakub Kruzik } 2268513960bSJakub Kruzik 227e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 228e662bc50SJakub Kruzik { 229e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 230e662bc50SJakub Kruzik PetscErrorCode ierr; 231e662bc50SJakub Kruzik 232e662bc50SJakub Kruzik PetscFunctionBegin; 233e662bc50SJakub Kruzik if (transpose) { 234e662bc50SJakub Kruzik def->Wt = W; 235e662bc50SJakub Kruzik def->W = NULL; 236e662bc50SJakub Kruzik } else { 237e662bc50SJakub Kruzik def->W = W; 238e662bc50SJakub Kruzik } 239e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 240e662bc50SJakub Kruzik PetscFunctionReturn(0); 241e662bc50SJakub Kruzik } 242e662bc50SJakub Kruzik 243e662bc50SJakub Kruzik /*@ 244e26ad46dSJakub Kruzik PCDeflationSetSpace - Set the deflation space matrix (or its (hermitian) transpose). 245e662bc50SJakub Kruzik 246e26ad46dSJakub Kruzik Logically Collective 247e662bc50SJakub Kruzik 248e662bc50SJakub Kruzik Input Parameters: 249e662bc50SJakub Kruzik + pc - the preconditioner context 250e662bc50SJakub Kruzik . W - deflation matrix 251e26ad46dSJakub Kruzik - transpose - indicates that W is an explicit transpose of the deflation matrix 252e26ad46dSJakub Kruzik 253e26ad46dSJakub Kruzik Notes: 254e26ad46dSJakub Kruzik Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel 255e26ad46dSJakub Kruzik deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and 256e26ad46dSJakub Kruzik the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with 257e26ad46dSJakub Kruzik W1 as the deflation matrix. This repeats until the maximum level set by 258e26ad46dSJakub Kruzik PCDeflationSetMaxLvl is reached or there are no more matrices available. 259e26ad46dSJakub Kruzik If there are matrices left after reaching the maximum level, 260e26ad46dSJakub Kruzik they are merged into a deflation matrix ...*W{n-1}*Wn. 261e662bc50SJakub Kruzik 262e662bc50SJakub Kruzik Level: intermediate 263e662bc50SJakub Kruzik 264e26ad46dSJakub Kruzik .seealso: PCDeflationSetMaxLvl(), PCDEFLATION 265e662bc50SJakub Kruzik @*/ 266e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 267e662bc50SJakub Kruzik { 268e662bc50SJakub Kruzik PetscErrorCode ierr; 269e662bc50SJakub Kruzik 270e662bc50SJakub Kruzik PetscFunctionBegin; 271e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 272e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 273e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 274e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 275e662bc50SJakub Kruzik PetscFunctionReturn(0); 276e662bc50SJakub Kruzik } 277e662bc50SJakub Kruzik 2783cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 2793cf3a049SJakub Kruzik { 2803cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2813cf3a049SJakub Kruzik PetscErrorCode ierr; 2823cf3a049SJakub Kruzik 2833cf3a049SJakub Kruzik PetscFunctionBegin; 2843cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 2853cf3a049SJakub Kruzik def->WtA = mat; 2863cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2873cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 2883cf3a049SJakub Kruzik PetscFunctionReturn(0); 2893cf3a049SJakub Kruzik } 2903cf3a049SJakub Kruzik 2913cf3a049SJakub Kruzik /*@ 292e26ad46dSJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A). 2933cf3a049SJakub Kruzik 294e26ad46dSJakub Kruzik Collective 2953cf3a049SJakub Kruzik 2963cf3a049SJakub Kruzik Input Parameters: 2973cf3a049SJakub Kruzik + pc - preconditioner context 2983cf3a049SJakub Kruzik - mat - projection null space matrix 2993cf3a049SJakub Kruzik 3003cf3a049SJakub Kruzik Level: developer 3013cf3a049SJakub Kruzik 3023cf3a049SJakub Kruzik .seealso: PCDEFLATION 3033cf3a049SJakub Kruzik @*/ 3043cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 3053cf3a049SJakub Kruzik { 3063cf3a049SJakub Kruzik PetscErrorCode ierr; 3073cf3a049SJakub Kruzik 3083cf3a049SJakub Kruzik PetscFunctionBegin; 3093cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3103cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 3113cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 3123cf3a049SJakub Kruzik PetscFunctionReturn(0); 3133cf3a049SJakub Kruzik } 3143cf3a049SJakub Kruzik 315e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 316e924b002SJakub Kruzik { 317e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 318e924b002SJakub Kruzik PetscErrorCode ierr; 319e924b002SJakub Kruzik 320e924b002SJakub Kruzik PetscFunctionBegin; 321e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 322e924b002SJakub Kruzik def->WtAW = mat; 323e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 324e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 325e924b002SJakub Kruzik PetscFunctionReturn(0); 326e924b002SJakub Kruzik } 327e924b002SJakub Kruzik 328e924b002SJakub Kruzik /*@ 329e26ad46dSJakub Kruzik PCDeflationSetCoarseMat - Set the coarse problem Mat. 330e924b002SJakub Kruzik 331e26ad46dSJakub Kruzik Collective 332e924b002SJakub Kruzik 333e924b002SJakub Kruzik Input Parameters: 334e924b002SJakub Kruzik + pc - preconditioner context 335e924b002SJakub Kruzik - mat - coarse problem mat 336e924b002SJakub Kruzik 337e924b002SJakub Kruzik Level: developer 338e924b002SJakub Kruzik 339e924b002SJakub Kruzik .seealso: PCDEFLATION 340e924b002SJakub Kruzik @*/ 341e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 342e924b002SJakub Kruzik { 343e924b002SJakub Kruzik PetscErrorCode ierr; 344e924b002SJakub Kruzik 345e924b002SJakub Kruzik PetscFunctionBegin; 346e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 347e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 348e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 349e924b002SJakub Kruzik PetscFunctionReturn(0); 350e924b002SJakub Kruzik } 351e924b002SJakub Kruzik 35298d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 3535c0c31c5SJakub Kruzik { 3545c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3555c0c31c5SJakub Kruzik 3565c0c31c5SJakub Kruzik PetscFunctionBegin; 35798d6e3deSJakub Kruzik *ksp = def->WtAWinv; 3585c0c31c5SJakub Kruzik PetscFunctionReturn(0); 3595c0c31c5SJakub Kruzik } 3605c0c31c5SJakub Kruzik 3615c0c31c5SJakub Kruzik /*@ 36298d6e3deSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 3635c0c31c5SJakub Kruzik 36498d6e3deSJakub Kruzik Not Collective 3655c0c31c5SJakub Kruzik 3665c0c31c5SJakub Kruzik Input Parameters: 36798d6e3deSJakub Kruzik . pc - preconditioner context 3685c0c31c5SJakub Kruzik 369e26ad46dSJakub Kruzik Output Parameters: 37098d6e3deSJakub Kruzik . ksp - coarse problem KSP context 3715c0c31c5SJakub Kruzik 372e26ad46dSJakub Kruzik Level: advanced 37398d6e3deSJakub Kruzik 37498d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 3755c0c31c5SJakub Kruzik @*/ 37698d6e3deSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3775c0c31c5SJakub Kruzik { 3785c0c31c5SJakub Kruzik PetscErrorCode ierr; 3795c0c31c5SJakub Kruzik 3805c0c31c5SJakub Kruzik PetscFunctionBegin; 38122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 38298d6e3deSJakub Kruzik PetscValidPointer(ksp,2); 38398d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 38498d6e3deSJakub Kruzik PetscFunctionReturn(0); 38598d6e3deSJakub Kruzik } 38698d6e3deSJakub Kruzik 38798d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 38898d6e3deSJakub Kruzik { 38998d6e3deSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 39098d6e3deSJakub Kruzik PetscErrorCode ierr; 39198d6e3deSJakub Kruzik 39298d6e3deSJakub Kruzik PetscFunctionBegin; 39398d6e3deSJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 39498d6e3deSJakub Kruzik def->WtAWinv = ksp; 39598d6e3deSJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 39698d6e3deSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 39798d6e3deSJakub Kruzik PetscFunctionReturn(0); 39898d6e3deSJakub Kruzik } 39998d6e3deSJakub Kruzik 40098d6e3deSJakub Kruzik /*@ 401e26ad46dSJakub Kruzik PCDeflationSetCoarseKSP - Set the coarse problem KSP. 40298d6e3deSJakub Kruzik 403e26ad46dSJakub Kruzik Collective 40498d6e3deSJakub Kruzik 40598d6e3deSJakub Kruzik Input Parameters: 40698d6e3deSJakub Kruzik + pc - preconditioner context 40798d6e3deSJakub Kruzik - ksp - coarse problem KSP context 40898d6e3deSJakub Kruzik 40998d6e3deSJakub Kruzik Level: developer 41098d6e3deSJakub Kruzik 41198d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 41298d6e3deSJakub Kruzik @*/ 41398d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 41498d6e3deSJakub Kruzik { 41598d6e3deSJakub Kruzik PetscErrorCode ierr; 41698d6e3deSJakub Kruzik 41798d6e3deSJakub Kruzik PetscFunctionBegin; 41898d6e3deSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 41998d6e3deSJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 42098d6e3deSJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 42198d6e3deSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 4225c0c31c5SJakub Kruzik PetscFunctionReturn(0); 4235c0c31c5SJakub Kruzik } 424e662bc50SJakub Kruzik 425268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 426268c6673SJakub Kruzik { 427268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 428268c6673SJakub Kruzik 429268c6673SJakub Kruzik PetscFunctionBegin; 430268c6673SJakub Kruzik *apc = def->pc; 431268c6673SJakub Kruzik PetscFunctionReturn(0); 432268c6673SJakub Kruzik } 433268c6673SJakub Kruzik 434268c6673SJakub Kruzik /*@ 435e26ad46dSJakub Kruzik PCDeflationGetPC - Returns a pointer to the additional preconditioner. 436268c6673SJakub Kruzik 437268c6673SJakub Kruzik Not Collective 438268c6673SJakub Kruzik 439268c6673SJakub Kruzik Input Parameters: 440268c6673SJakub Kruzik . pc - the preconditioner context 441268c6673SJakub Kruzik 442e26ad46dSJakub Kruzik Output Parameters: 443268c6673SJakub Kruzik . apc - additional preconditioner 444268c6673SJakub Kruzik 445268c6673SJakub Kruzik Level: advanced 446268c6673SJakub Kruzik 447268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 448268c6673SJakub Kruzik @*/ 449268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 450268c6673SJakub Kruzik { 451268c6673SJakub Kruzik PetscErrorCode ierr; 452268c6673SJakub Kruzik 453268c6673SJakub Kruzik PetscFunctionBegin; 454268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 455268c6673SJakub Kruzik PetscValidPointer(pc,2); 456268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 457268c6673SJakub Kruzik PetscFunctionReturn(0); 458268c6673SJakub Kruzik } 459268c6673SJakub Kruzik 46022b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 46122b0793eSJakub Kruzik { 46222b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 46322b0793eSJakub Kruzik PetscErrorCode ierr; 46422b0793eSJakub Kruzik 46522b0793eSJakub Kruzik PetscFunctionBegin; 46622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 46722b0793eSJakub Kruzik def->pc = apc; 46822b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 46922b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 47022b0793eSJakub Kruzik PetscFunctionReturn(0); 47122b0793eSJakub Kruzik } 47222b0793eSJakub Kruzik 47322b0793eSJakub Kruzik /*@ 474e26ad46dSJakub Kruzik PCDeflationSetPC - Set the additional preconditioner. 47522b0793eSJakub Kruzik 476e26ad46dSJakub Kruzik Collective 47722b0793eSJakub Kruzik 47822b0793eSJakub Kruzik Input Parameters: 47922b0793eSJakub Kruzik + pc - the preconditioner context 48022b0793eSJakub Kruzik - apc - additional preconditioner 48122b0793eSJakub Kruzik 482268c6673SJakub Kruzik Level: developer 48322b0793eSJakub Kruzik 484268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 48522b0793eSJakub Kruzik @*/ 48622b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 48722b0793eSJakub Kruzik { 48822b0793eSJakub Kruzik PetscErrorCode ierr; 48922b0793eSJakub Kruzik 49022b0793eSJakub Kruzik PetscFunctionBegin; 49122b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 49222b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 49322b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 49422b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 49522b0793eSJakub Kruzik PetscFunctionReturn(0); 49622b0793eSJakub Kruzik } 49722b0793eSJakub Kruzik 49837eeb815SJakub Kruzik /* 499b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 500b27c8092SJakub Kruzik */ 501b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 502b27c8092SJakub Kruzik { 503b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 504b27c8092SJakub Kruzik Mat A; 505b27c8092SJakub Kruzik Vec r,w1,w2; 506b27c8092SJakub Kruzik PetscBool nonzero; 507b27c8092SJakub Kruzik PetscErrorCode ierr; 50837eeb815SJakub Kruzik 509b27c8092SJakub Kruzik PetscFunctionBegin; 510b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 511b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 512b27c8092SJakub Kruzik r = def->work; 513b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 51437eeb815SJakub Kruzik 515b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 516b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 517b27c8092SJakub Kruzik if (nonzero) { 518b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 519b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 520b27c8092SJakub Kruzik } else { 521b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 522b27c8092SJakub Kruzik } 523b27c8092SJakub Kruzik 524a3177931SJakub Kruzik if (def->Wt) { 525a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 526a3177931SJakub Kruzik } else { 527a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 528a3177931SJakub Kruzik } 529b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 530b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 531b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 532b27c8092SJakub Kruzik PetscFunctionReturn(0); 533b27c8092SJakub Kruzik } 53437eeb815SJakub Kruzik 535f8bf229cSJakub Kruzik /* 536f8bf229cSJakub Kruzik if (def->correct) { 537ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 538f8bf229cSJakub Kruzik } else { 539ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 540f8bf229cSJakub Kruzik } 54137eeb815SJakub Kruzik */ 542f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 543f8bf229cSJakub Kruzik { 544f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 545f8bf229cSJakub Kruzik Mat A; 546f8bf229cSJakub Kruzik Vec u,w1,w2; 547f8bf229cSJakub Kruzik PetscErrorCode ierr; 548f8bf229cSJakub Kruzik 549f8bf229cSJakub Kruzik PetscFunctionBegin; 550f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 551f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 552f8bf229cSJakub Kruzik u = def->work; 553f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 554f8bf229cSJakub Kruzik 555ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 55622b0793eSJakub Kruzik if (!def->init) { 557ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 558f8bf229cSJakub Kruzik if (def->correct) { 559ae029463SJakub Kruzik if (def->Wt) { 560ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 561ae029463SJakub Kruzik } else { 562a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 563f8bf229cSJakub Kruzik } 564ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 565f8bf229cSJakub Kruzik } 566f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 567f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 56822b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 56922b0793eSJakub Kruzik } 570f8bf229cSJakub Kruzik PetscFunctionReturn(0); 571f8bf229cSJakub Kruzik } 572f8bf229cSJakub Kruzik 57337eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 57437eeb815SJakub Kruzik { 575409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 576409a357bSJakub Kruzik KSP innerksp; 577409a357bSJakub Kruzik PC pcinner; 578409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 579409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 580409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 581409a357bSJakub Kruzik MatCompositeType ctype; 582409a357bSJakub Kruzik MPI_Comm comm; 5836c93e71cSJakub Kruzik char prefix[128]=""; 58437eeb815SJakub Kruzik PetscErrorCode ierr; 58537eeb815SJakub Kruzik 58637eeb815SJakub Kruzik PetscFunctionBegin; 587409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 588409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 589f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 5906c93e71cSJakub Kruzik if (!def->lvl && !def->prefix) { 5916c93e71cSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&def->prefix);CHKERRQ(ierr); 5926c93e71cSJakub Kruzik } 5936c93e71cSJakub Kruzik if (def->lvl) { 5946c93e71cSJakub Kruzik sprintf(prefix,"%d_",(int)def->lvl); 5956c93e71cSJakub Kruzik } 596f8bf229cSJakub Kruzik 597f8bf229cSJakub Kruzik /* compute a deflation space */ 598409a357bSJakub Kruzik if (def->W || def->Wt) { 599409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 600409a357bSJakub Kruzik } else { 601e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 60237eeb815SJakub Kruzik } 60337eeb815SJakub Kruzik 604409a357bSJakub Kruzik /* nested deflation */ 605409a357bSJakub Kruzik if (def->W) { 606409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 607409a357bSJakub Kruzik if (match) { 608409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 609409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 61037eeb815SJakub Kruzik } 611409a357bSJakub Kruzik } else { 612a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 613409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 614409a357bSJakub Kruzik if (match) { 615409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 616409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 617409a357bSJakub Kruzik } 618409a357bSJakub Kruzik transp = PETSC_TRUE; 619409a357bSJakub Kruzik } 620409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 621409a357bSJakub Kruzik if (!transp) { 6226c93e71cSJakub Kruzik if (def->lvl < def->maxlvl) { 62328da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 624409a357bSJakub Kruzik for (i=0; i<size; i++) { 625409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 626409a357bSJakub Kruzik } 627409a357bSJakub Kruzik size -= 1; 628409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 629409a357bSJakub Kruzik def->W = mats[size]; 630409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 631409a357bSJakub Kruzik if (size > 1) { 632409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 633409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 634409a357bSJakub Kruzik } else { 635409a357bSJakub Kruzik nextDef = mats[0]; 636409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 637409a357bSJakub Kruzik } 63828da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 639409a357bSJakub Kruzik } else { 640409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 641409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 642409a357bSJakub Kruzik } 64328da0170SJakub Kruzik } else { 6446c93e71cSJakub Kruzik if (def->lvl < def->maxlvl) { 64528da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 64628da0170SJakub Kruzik for (i=0; i<size; i++) { 64728da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 64828da0170SJakub Kruzik } 64928da0170SJakub Kruzik size -= 1; 65028da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 65128da0170SJakub Kruzik def->Wt = mats[0]; 65228da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 65328da0170SJakub Kruzik if (size > 1) { 65428da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 65528da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 65628da0170SJakub Kruzik } else { 65728da0170SJakub Kruzik nextDef = mats[1]; 65828da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 659409a357bSJakub Kruzik } 660409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 66128da0170SJakub Kruzik } else { 66228da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 66328da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 66428da0170SJakub Kruzik } 66528da0170SJakub Kruzik } 66628da0170SJakub Kruzik } 66728da0170SJakub Kruzik 66828da0170SJakub Kruzik if (transp) { 66928da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 670a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 671409a357bSJakub Kruzik } 672409a357bSJakub Kruzik 673ae029463SJakub Kruzik /* assemble WtA */ 674ae029463SJakub Kruzik if (!def->WtA) { 675ae029463SJakub Kruzik if (def->Wt) { 676ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 677ae029463SJakub Kruzik } else { 678a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 679a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 680a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 681a3177931SJakub Kruzik #else 682ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 683a3177931SJakub Kruzik #endif 684ae029463SJakub Kruzik } 685ae029463SJakub Kruzik } 686409a357bSJakub Kruzik /* setup coarse problem */ 687409a357bSJakub Kruzik if (!def->WtAWinv) { 68828da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 689409a357bSJakub Kruzik if (!def->WtAW) { 690ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 691409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 692409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 693409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 694409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 695ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 696409a357bSJakub Kruzik PetscReal *norms; 697409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 698409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 699409a357bSJakub Kruzik for (i=0; i<m; i++) { 700409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 701409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 702409a357bSJakub Kruzik } 703409a357bSJakub Kruzik } 704409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 705409a357bSJakub Kruzik #endif 706409a357bSJakub Kruzik } else { 707409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 708409a357bSJakub Kruzik } 709409a357bSJakub Kruzik /* TODO use MATINV */ 710409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 711409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 712409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 713557a2f7dSJakub Kruzik /* Setup KSP and PC */ 714557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 715557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 716557a2f7dSJakub Kruzik /* set default KSPtype */ 717557a2f7dSJakub Kruzik if (!def->ksptype) { 718557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 719557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 720557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 721557a2f7dSJakub Kruzik } 722557a2f7dSJakub Kruzik } 723ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 724557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 725557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 7266c93e71cSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->lvl+1,def->maxlvl);CHKERRQ(ierr); 727ae029463SJakub Kruzik /* inherit options */ 7286c93e71cSJakub Kruzik if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix; 729*9f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->init = def->init; 730557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 731557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 732*9f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->correctfact = def->correctfact; 733*9f604af8SJakub Kruzik ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact; 734557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 735557a2f7dSJakub Kruzik } else { /* the last level */ 736557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 737409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 7386c93e71cSJakub Kruzik /* do not overwrite PCTELESCOPE */ 7396c93e71cSJakub Kruzik if (def->prefix) { 7406c93e71cSJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,def->prefix);CHKERRQ(ierr); 741ac66f006SJakub Kruzik } 7426c93e71cSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"def_tel_");CHKERRQ(ierr); 743ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 744409a357bSJakub Kruzik /* Reduction factor choice */ 745409a357bSJakub Kruzik red = def->reductionfact; 746409a357bSJakub Kruzik if (red < 0) { 747409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 748409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 749409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 750409a357bSJakub Kruzik if (match) red = commsize; 7516c93e71cSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); 752409a357bSJakub Kruzik } 753409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 754ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 755409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 756ac66f006SJakub Kruzik if (innerksp) { 757409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 758409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 759ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 760409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 761409a357bSJakub Kruzik if (commsize == red) { 762ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 763409a357bSJakub Kruzik } else { 764ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 765409a357bSJakub Kruzik } 766409a357bSJakub Kruzik } 767557a2f7dSJakub Kruzik } 768557a2f7dSJakub Kruzik 769557a2f7dSJakub Kruzik if (innerksp) { 7706c93e71cSJakub Kruzik if (def->prefix) { 7716c93e71cSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,def->prefix);CHKERRQ(ierr); 77222b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 7736c93e71cSJakub Kruzik } else { 7746c93e71cSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 7756c93e71cSJakub Kruzik } 7766c93e71cSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 777557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 778557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 779ac66f006SJakub Kruzik } 780409a357bSJakub Kruzik } 781557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 782557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 783409a357bSJakub Kruzik 78422b0793eSJakub Kruzik if (!def->pc) { 78522b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 78622b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 78722b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 7886c93e71cSJakub Kruzik if (def->prefix) { 7896c93e71cSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,def->prefix);CHKERRQ(ierr); 79022b0793eSJakub Kruzik } 7916c93e71cSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_");CHKERRQ(ierr); 7926c93e71cSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 7936c93e71cSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"pc_");CHKERRQ(ierr); 79422b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 79522b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 79622b0793eSJakub Kruzik } 79722b0793eSJakub Kruzik 798f8bf229cSJakub Kruzik /* create work vecs */ 799f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 800f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 80137eeb815SJakub Kruzik PetscFunctionReturn(0); 80237eeb815SJakub Kruzik } 803b30d4775SJakub Kruzik 80437eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 80537eeb815SJakub Kruzik { 806b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 80737eeb815SJakub Kruzik PetscErrorCode ierr; 80837eeb815SJakub Kruzik 80937eeb815SJakub Kruzik PetscFunctionBegin; 810b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 811b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 812b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 813b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 814ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 815b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 816b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 81722b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 81837eeb815SJakub Kruzik PetscFunctionReturn(0); 81937eeb815SJakub Kruzik } 82037eeb815SJakub Kruzik 82137eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 82237eeb815SJakub Kruzik { 82337eeb815SJakub Kruzik PetscErrorCode ierr; 82437eeb815SJakub Kruzik 82537eeb815SJakub Kruzik PetscFunctionBegin; 82637eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 82737eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 82837eeb815SJakub Kruzik PetscFunctionReturn(0); 82937eeb815SJakub Kruzik } 83037eeb815SJakub Kruzik 8318513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 8328513960bSJakub Kruzik { 8338513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 8341ac6250aSJakub Kruzik PetscInt its; 8358513960bSJakub Kruzik PetscBool iascii; 8368513960bSJakub Kruzik PetscErrorCode ierr; 8378513960bSJakub Kruzik 8388513960bSJakub Kruzik PetscFunctionBegin; 8398513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8408513960bSJakub Kruzik if (iascii) { 8418513960bSJakub Kruzik if (def->correct) { 8421ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n", 8431ac6250aSJakub Kruzik (double)PetscRealPart(def->correctfact), 8441ac6250aSJakub Kruzik (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr); 8458513960bSJakub Kruzik } 8466c93e71cSJakub Kruzik if (!def->lvl) { 8471ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 8488513960bSJakub Kruzik } 8498513960bSJakub Kruzik 8501ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr); 85122b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 85222b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 85322b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 85422b0793eSJakub Kruzik 8551ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr); 8568513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 8571ac6250aSJakub Kruzik ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr); 8581ac6250aSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr); 8598513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 8608513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 8618513960bSJakub Kruzik } 8628513960bSJakub Kruzik PetscFunctionReturn(0); 8638513960bSJakub Kruzik } 8648513960bSJakub Kruzik 86537eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 86637eeb815SJakub Kruzik { 867880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 86837eeb815SJakub Kruzik PetscErrorCode ierr; 86937eeb815SJakub Kruzik 87037eeb815SJakub Kruzik PetscFunctionBegin; 87137eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 872a122ebaeSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr); 8736c93e71cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxlvl,&def->maxlvl,NULL);CHKERRQ(ierr); 874859a873cSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 8758a71cb68SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr); 8768a71cb68SJakub 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); 877880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 878880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 879880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 880880d05ceSJakub Kruzik //TODO add set function and fix manpages 88137eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 88237eeb815SJakub Kruzik PetscFunctionReturn(0); 88337eeb815SJakub Kruzik } 88437eeb815SJakub Kruzik 88537eeb815SJakub Kruzik /*MC 886e26ad46dSJakub Kruzik PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value. 88737eeb815SJakub Kruzik 888e26ad46dSJakub Kruzik Options Database Keys: 889e26ad46dSJakub Kruzik + -pc_deflation_init_only <false> - if true computes only the special guess 890e26ad46dSJakub Kruzik . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation 891e26ad46dSJakub Kruzik . -pc_deflation_reduction_factor <-1> - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem) 892e26ad46dSJakub Kruzik . -pc_deflation_correction <false> - if true apply coarse problem correction 893e26ad46dSJakub Kruzik . -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor 894e26ad46dSJakub Kruzik . -pc_deflation_compute_space <haar> - compute PCDeflationSpaceType deflation space 895e26ad46dSJakub Kruzik - -pc_deflation_compute_space_size <1> - size of the deflation space (corresponds to number of levels for wavelet-based deflation) 89637eeb815SJakub Kruzik 89737eeb815SJakub Kruzik Notes: 898e26ad46dSJakub Kruzik Given a (complex - transpose is always hermitian) full rank deflation matrix W, the deflation (introduced in [1,2]) 899e26ad46dSJakub Kruzik preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat. 900e26ad46dSJakub Kruzik 901e26ad46dSJakub Kruzik The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space. 902e26ad46dSJakub Kruzik If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the 903e26ad46dSJakub Kruzik preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner 904e26ad46dSJakub Kruzik application consists of P*M{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues 905e26ad46dSJakub Kruzik to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some 906e26ad46dSJakub Kruzik eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates 907e26ad46dSJakub Kruzik of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper. 908e26ad46dSJakub Kruzik 909e26ad46dSJakub Kruzik The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can 910e26ad46dSJakub Kruzik be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size. 911e26ad46dSJakub Kruzik User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix 912e26ad46dSJakub Kruzik is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the 913e26ad46dSJakub Kruzik deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by FCG (if A is MAT_SPD) or FGMRES preconditioned 914e26ad46dSJakub Kruzik by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum 915e26ad46dSJakub Kruzik level is reached or there are no more matrices. If the maximal level is reached, the remaining matrices are merged 916e26ad46dSJakub Kruzik (multiplied) to create the last deflation matrix. The maximal level defaults to 0 and can be set by 917e26ad46dSJakub Kruzik PCDeflationSetMaxLvl() or by -pc_deflation_max_lvl. 918e26ad46dSJakub Kruzik 9196c93e71cSJakub Kruzik The coarse problem KSP can be controlled from the command line with prefix -def_ for the first level and -def_[lvl-1] 9206c93e71cSJakub Kruzik from the second level onward. You can also use 921e26ad46dSJakub Kruzik PCDeflationGetCoarseKSP() or PCDeflationSetCoarseKSP() to control it from code. The bottom level KSP defaults to 922e26ad46dSJakub Kruzik KSPREONLY with PCLU direct solver wrapped into PCTELESCOPE. For convenience, the reduction factor can be set by 923e26ad46dSJakub Kruzik PCDeflationSetReductionFactor() or -pc_deflation_recduction_factor. The default is chosen heuristically based on the 924e26ad46dSJakub Kruzik coarse problem size. 925e26ad46dSJakub Kruzik 9266c93e71cSJakub Kruzik The additional preconditioner can be controlled from command line with prefix -def_[lvl]_pc (same rules used for 9276c93e71cSJakub Kruzik coarse problem KSP apply for [lvl]_ part of prefix), e.g., -def_1_pc_pc_type bjacobi. You can also use 9286c93e71cSJakub Kruzik PCDeflationGetPC() or PCDeflationSetPC() to control the additional preconditioner from code. It defaults to PCNONE. 929e26ad46dSJakub Kruzik 930e26ad46dSJakub Kruzik The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can 931e26ad46dSJakub Kruzik be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can 932e26ad46dSJakub Kruzik significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We 933e26ad46dSJakub Kruzik recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create 934e26ad46dSJakub Kruzik an isolated eigenvalue. 935e26ad46dSJakub Kruzik 936*9f604af8SJakub Kruzik The options are automatically inherited from previous deflation level. 937*9f604af8SJakub Kruzik 938e26ad46dSJakub Kruzik The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also 939e26ad46dSJakub Kruzik recommend limiting the number of iterations for the coarse problem. 940e26ad46dSJakub Kruzik 941e26ad46dSJakub Kruzik See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients. 942e26ad46dSJakub Kruzik Section 4 describes some possible choices for the deflation space. 943e26ad46dSJakub Kruzik 944e26ad46dSJakub Kruzik Developer Notes: 945e26ad46dSJakub Kruzik Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech 946e26ad46dSJakub Kruzik Academy of Sciences and VSB - TU Ostrava. 947e26ad46dSJakub Kruzik 948e26ad46dSJakub Kruzik Developed from PERMON code used in [4] while on a research stay with 949e26ad46dSJakub Kruzik Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin. 950e26ad46dSJakub Kruzik 951e26ad46dSJakub Kruzik References: 952e26ad46dSJakub Kruzik + [1] - A. Nicolaides. “Deflation of conjugate gradients with applications to boundary valueproblems”, SIAM J. Numer. Anal. 24.2, 1987. 953e26ad46dSJakub Kruzik . [2] - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988. 954e26ad46dSJakub Kruzik . [3] - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008. 955e26ad46dSJakub Kruzik . [4] - 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 956e26ad46dSJakub Kruzik 957e26ad46dSJakub Kruzik Level: intermediate 95837eeb815SJakub Kruzik 95937eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 960e26ad46dSJakub Kruzik PCDeflationSetInitOnly(), PCDeflationSetMaxLvl(), PCDeflationSetReductionFactor(), 961e26ad46dSJakub Kruzik PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompoute(), 962e26ad46dSJakub Kruzik PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(), 963e26ad46dSJakub Kruzik PCDeflationSetCoarseMat(), PCDeflationSetCoarseKSP(), PCDeflationGetCoarseKSP(), 964e26ad46dSJakub Kruzik PCDeflationGetPC(), PCDeflationSetPC() 96537eeb815SJakub Kruzik M*/ 96637eeb815SJakub Kruzik 96737eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 96837eeb815SJakub Kruzik { 96937eeb815SJakub Kruzik PC_Deflation *def; 97037eeb815SJakub Kruzik PetscErrorCode ierr; 97137eeb815SJakub Kruzik 97237eeb815SJakub Kruzik PetscFunctionBegin; 97337eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 97437eeb815SJakub Kruzik pc->data = (void*)def; 97537eeb815SJakub Kruzik 976e662bc50SJakub Kruzik def->init = PETSC_FALSE; 977e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 978fcb31d99SJakub Kruzik def->correctfact = 1.0; 979e662bc50SJakub Kruzik def->reductionfact = -1; 980e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 981e662bc50SJakub Kruzik def->spacesize = 1; 982e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 9836c93e71cSJakub Kruzik def->lvl = 0; 9846c93e71cSJakub Kruzik def->maxlvl = 0; 98537eeb815SJakub Kruzik 98637eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 987b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 98837eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 98937eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 99037eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 99137eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 9928513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 99337eeb815SJakub Kruzik 994a122ebaeSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr); 99598d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 996859a873cSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr); 9978a71cb68SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr); 99839417d7eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 999e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 10003cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 1001e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 10024a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 1003f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 100498d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 100598d6e3deSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 100637eeb815SJakub Kruzik PetscFunctionReturn(0); 100737eeb815SJakub Kruzik } 100837eeb815SJakub Kruzik 1009