137eeb815SJakub Kruzik 237eeb815SJakub Kruzik /* -------------------------------------------------------------------- 337eeb815SJakub Kruzik 437eeb815SJakub Kruzik This file implements a Deflation preconditioner in PETSc as part of PC. 537eeb815SJakub Kruzik You can use this as a starting point for implementing your own 637eeb815SJakub Kruzik preconditioner that is not provided with PETSc. (You might also consider 737eeb815SJakub Kruzik just using PCSHELL) 837eeb815SJakub Kruzik 937eeb815SJakub Kruzik The following basic routines are required for each preconditioner. 1037eeb815SJakub Kruzik PCCreate_XXX() - Creates a preconditioner context 1137eeb815SJakub Kruzik PCSetFromOptions_XXX() - Sets runtime options 1237eeb815SJakub Kruzik PCApply_XXX() - Applies the preconditioner 1337eeb815SJakub Kruzik PCDestroy_XXX() - Destroys the preconditioner context 1437eeb815SJakub Kruzik where the suffix "_XXX" denotes a particular implementation, in 1537eeb815SJakub Kruzik this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation). 1637eeb815SJakub Kruzik These routines are actually called via the common user interface 1737eeb815SJakub Kruzik routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 1837eeb815SJakub Kruzik so the application code interface remains identical for all 1937eeb815SJakub Kruzik preconditioners. 2037eeb815SJakub Kruzik 2137eeb815SJakub Kruzik Another key routine is: 2237eeb815SJakub Kruzik PCSetUp_XXX() - Prepares for the use of a preconditioner 2337eeb815SJakub Kruzik by setting data structures and options. The interface routine PCSetUp() 2437eeb815SJakub Kruzik is not usually called directly by the user, but instead is called by 2537eeb815SJakub Kruzik PCApply() if necessary. 2637eeb815SJakub Kruzik 2737eeb815SJakub Kruzik Additional basic routines are: 2837eeb815SJakub Kruzik PCView_XXX() - Prints details of runtime options that 2937eeb815SJakub Kruzik have actually been used. 3037eeb815SJakub Kruzik These are called by application codes via the interface routines 3137eeb815SJakub Kruzik PCView(). 3237eeb815SJakub Kruzik 3337eeb815SJakub Kruzik The various types of solvers (preconditioners, Krylov subspace methods, 3437eeb815SJakub Kruzik nonlinear solvers, timesteppers) are all organized similarly, so the 3537eeb815SJakub Kruzik above description applies to these categories also. One exception is 3637eeb815SJakub Kruzik that the analogues of PCApply() for these components are KSPSolve(), 3737eeb815SJakub Kruzik SNESSolve(), and TSSolve(). 3837eeb815SJakub Kruzik 3937eeb815SJakub Kruzik Additional optional functionality unique to preconditioners is left and 4037eeb815SJakub Kruzik right symmetric preconditioner application via PCApplySymmetricLeft() 4137eeb815SJakub Kruzik and PCApplySymmetricRight(). The Deflation implementation is 4237eeb815SJakub Kruzik PCApplySymmetricLeftOrRight_Deflation(). 4337eeb815SJakub Kruzik 4437eeb815SJakub Kruzik -------------------------------------------------------------------- */ 4537eeb815SJakub Kruzik 4637eeb815SJakub Kruzik /* 4737eeb815SJakub Kruzik Include files needed for the Deflation preconditioner: 4837eeb815SJakub Kruzik pcimpl.h - private include file intended for use by all preconditioners 4937eeb815SJakub Kruzik */ 5037eeb815SJakub Kruzik 51292e2e67SJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 5237eeb815SJakub Kruzik 538513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 548513960bSJakub Kruzik "haar", 558513960bSJakub Kruzik "jacket-haar", 568513960bSJakub Kruzik "db2", 578513960bSJakub Kruzik "db4", 588513960bSJakub Kruzik "db8", 598513960bSJakub Kruzik "db16", 608513960bSJakub Kruzik "biorth22", 618513960bSJakub Kruzik "meyer", 628513960bSJakub Kruzik "aggregation", 638513960bSJakub Kruzik "slepc", 648513960bSJakub Kruzik "slepc-cheap", 658513960bSJakub Kruzik "user", 668513960bSJakub Kruzik "DdefSpaceType", 678513960bSJakub Kruzik "Ddef_SPACE_", 688513960bSJakub Kruzik 0 698513960bSJakub Kruzik }; 708513960bSJakub Kruzik 71*39417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size) 72*39417d7eSJakub Kruzik { 73*39417d7eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 74*39417d7eSJakub Kruzik 75*39417d7eSJakub Kruzik PetscFunctionBegin; 76*39417d7eSJakub Kruzik if (type) def->spacetype = type; 77*39417d7eSJakub Kruzik if (size > 0) def->spacesize = size; 78*39417d7eSJakub Kruzik PetscFunctionReturn(0); 79*39417d7eSJakub Kruzik } 80*39417d7eSJakub Kruzik 81*39417d7eSJakub Kruzik /*@ 82*39417d7eSJakub Kruzik PCDeflationSetSpaceToCompute - Set deflation space type and size to compute. 83*39417d7eSJakub Kruzik 84*39417d7eSJakub Kruzik Logically Collective on PC 85*39417d7eSJakub Kruzik 86*39417d7eSJakub Kruzik Input Parameters: 87*39417d7eSJakub Kruzik + pc - the preconditioner context 88*39417d7eSJakub Kruzik . type - deflation space type to compute (or PETSC_IGNORE) 89*39417d7eSJakub Kruzik - size - size of the space to compute (or PETSC_DEFAULT) 90*39417d7eSJakub Kruzik 91*39417d7eSJakub Kruzik Notes: 92*39417d7eSJakub Kruzik For wavelet-based deflation, size represents number of levels. 93*39417d7eSJakub Kruzik The deflation space is computed in PCSetUP(). 94*39417d7eSJakub Kruzik 95*39417d7eSJakub Kruzik Level: intermediate 96*39417d7eSJakub Kruzik 97*39417d7eSJakub Kruzik .seealso: PCDEFLATION 98*39417d7eSJakub Kruzik @*/ 99*39417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size) 100*39417d7eSJakub Kruzik { 101*39417d7eSJakub Kruzik PetscErrorCode ierr; 102*39417d7eSJakub Kruzik 103*39417d7eSJakub Kruzik PetscFunctionBegin; 104*39417d7eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 105*39417d7eSJakub Kruzik if (type) PetscValidLogicalCollectiveEnum(pc,type,2); 106*39417d7eSJakub Kruzik if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3); 107*39417d7eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr); 108*39417d7eSJakub Kruzik PetscFunctionReturn(0); 109*39417d7eSJakub Kruzik } 1108513960bSJakub Kruzik 111e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 112e662bc50SJakub Kruzik { 113e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 114e662bc50SJakub Kruzik PetscErrorCode ierr; 115e662bc50SJakub Kruzik 116e662bc50SJakub Kruzik PetscFunctionBegin; 117e662bc50SJakub Kruzik if (transpose) { 118e662bc50SJakub Kruzik def->Wt = W; 119e662bc50SJakub Kruzik def->W = NULL; 120e662bc50SJakub Kruzik } else { 121e662bc50SJakub Kruzik def->W = W; 122e662bc50SJakub Kruzik } 123e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 124e662bc50SJakub Kruzik PetscFunctionReturn(0); 125e662bc50SJakub Kruzik } 126e662bc50SJakub Kruzik 127e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 128e662bc50SJakub Kruzik /*@ 129e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 130e662bc50SJakub Kruzik 131e662bc50SJakub Kruzik Logically Collective on PC 132e662bc50SJakub Kruzik 133e662bc50SJakub Kruzik Input Parameters: 134e662bc50SJakub Kruzik + pc - the preconditioner context 135e662bc50SJakub Kruzik . W - deflation matrix 136e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 137e662bc50SJakub Kruzik 138e662bc50SJakub Kruzik Level: intermediate 139e662bc50SJakub Kruzik 140e662bc50SJakub Kruzik .seealso: PCDEFLATION 141e662bc50SJakub Kruzik @*/ 142e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 143e662bc50SJakub Kruzik { 144e662bc50SJakub Kruzik PetscErrorCode ierr; 145e662bc50SJakub Kruzik 146e662bc50SJakub Kruzik PetscFunctionBegin; 147e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 148e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 149e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 150e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 151e662bc50SJakub Kruzik PetscFunctionReturn(0); 152e662bc50SJakub Kruzik } 153e662bc50SJakub Kruzik 1543cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat) 1553cf3a049SJakub Kruzik { 1563cf3a049SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1573cf3a049SJakub Kruzik PetscErrorCode ierr; 1583cf3a049SJakub Kruzik 1593cf3a049SJakub Kruzik PetscFunctionBegin; 1603cf3a049SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 1613cf3a049SJakub Kruzik def->WtA = mat; 1623cf3a049SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1633cf3a049SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr); 1643cf3a049SJakub Kruzik PetscFunctionReturn(0); 1653cf3a049SJakub Kruzik } 1663cf3a049SJakub Kruzik 1673cf3a049SJakub Kruzik /*@ 1683cf3a049SJakub Kruzik PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A). 1693cf3a049SJakub Kruzik 1703cf3a049SJakub Kruzik Not Collective 1713cf3a049SJakub Kruzik 1723cf3a049SJakub Kruzik Input Parameters: 1733cf3a049SJakub Kruzik + pc - preconditioner context 1743cf3a049SJakub Kruzik - mat - projection null space matrix 1753cf3a049SJakub Kruzik 1763cf3a049SJakub Kruzik Level: developer 1773cf3a049SJakub Kruzik 1783cf3a049SJakub Kruzik .seealso: PCDEFLATION 1793cf3a049SJakub Kruzik @*/ 1803cf3a049SJakub Kruzik PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat) 1813cf3a049SJakub Kruzik { 1823cf3a049SJakub Kruzik PetscErrorCode ierr; 1833cf3a049SJakub Kruzik 1843cf3a049SJakub Kruzik PetscFunctionBegin; 1853cf3a049SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1863cf3a049SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 1873cf3a049SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 1883cf3a049SJakub Kruzik PetscFunctionReturn(0); 1893cf3a049SJakub Kruzik } 1903cf3a049SJakub Kruzik 191e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat) 192e924b002SJakub Kruzik { 193e924b002SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 194e924b002SJakub Kruzik PetscErrorCode ierr; 195e924b002SJakub Kruzik 196e924b002SJakub Kruzik PetscFunctionBegin; 197e924b002SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 198e924b002SJakub Kruzik def->WtAW = mat; 199e924b002SJakub Kruzik ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 200e924b002SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr); 201e924b002SJakub Kruzik PetscFunctionReturn(0); 202e924b002SJakub Kruzik } 203e924b002SJakub Kruzik 204e924b002SJakub Kruzik /*@ 205e924b002SJakub Kruzik PCDeflationSetCoarseMat - Set coarse problem Mat. 206e924b002SJakub Kruzik 207e924b002SJakub Kruzik Not Collective 208e924b002SJakub Kruzik 209e924b002SJakub Kruzik Input Parameters: 210e924b002SJakub Kruzik + pc - preconditioner context 211e924b002SJakub Kruzik - mat - coarse problem mat 212e924b002SJakub Kruzik 213e924b002SJakub Kruzik Level: developer 214e924b002SJakub Kruzik 215e924b002SJakub Kruzik .seealso: PCDEFLATION 216e924b002SJakub Kruzik @*/ 217e924b002SJakub Kruzik PetscErrorCode PCDeflationSetCoarseMat(PC pc,Mat mat) 218e924b002SJakub Kruzik { 219e924b002SJakub Kruzik PetscErrorCode ierr; 220e924b002SJakub Kruzik 221e924b002SJakub Kruzik PetscFunctionBegin; 222e924b002SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 223e924b002SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 224e924b002SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr); 225e924b002SJakub Kruzik PetscFunctionReturn(0); 226e924b002SJakub Kruzik } 227e924b002SJakub Kruzik 2285c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 2295c0c31c5SJakub Kruzik { 2305c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 2315c0c31c5SJakub Kruzik 2325c0c31c5SJakub Kruzik PetscFunctionBegin; 23381e2347eSJakub Kruzik if (current) def->nestedlvl = current; 2345c0c31c5SJakub Kruzik def->maxnestedlvl = max; 2355c0c31c5SJakub Kruzik PetscFunctionReturn(0); 2365c0c31c5SJakub Kruzik } 2375c0c31c5SJakub Kruzik 2385c0c31c5SJakub Kruzik /*@ 2395c0c31c5SJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 2405c0c31c5SJakub Kruzik 2415c0c31c5SJakub Kruzik Logically Collective on PC 2425c0c31c5SJakub Kruzik 2435c0c31c5SJakub Kruzik Input Parameters: 2445c0c31c5SJakub Kruzik + pc - the preconditioner context 24522b0793eSJakub Kruzik - max - maximum deflation level 2465c0c31c5SJakub Kruzik 2475c0c31c5SJakub Kruzik Level: intermediate 2485c0c31c5SJakub Kruzik 2495c0c31c5SJakub Kruzik .seealso: PCDEFLATION 2505c0c31c5SJakub Kruzik @*/ 2515c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 2525c0c31c5SJakub Kruzik { 2535c0c31c5SJakub Kruzik PetscErrorCode ierr; 2545c0c31c5SJakub Kruzik 2555c0c31c5SJakub Kruzik PetscFunctionBegin; 25622b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2575c0c31c5SJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 2585c0c31c5SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 2595c0c31c5SJakub Kruzik PetscFunctionReturn(0); 2605c0c31c5SJakub Kruzik } 261e662bc50SJakub Kruzik 262268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc) 263268c6673SJakub Kruzik { 264268c6673SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 265268c6673SJakub Kruzik 266268c6673SJakub Kruzik PetscFunctionBegin; 267268c6673SJakub Kruzik *apc = def->pc; 268268c6673SJakub Kruzik PetscFunctionReturn(0); 269268c6673SJakub Kruzik } 270268c6673SJakub Kruzik 271268c6673SJakub Kruzik /*@ 272268c6673SJakub Kruzik PCDeflationGetPC - Returns a pointer to additional preconditioner. 273268c6673SJakub Kruzik 274268c6673SJakub Kruzik Not Collective 275268c6673SJakub Kruzik 276268c6673SJakub Kruzik Input Parameters: 277268c6673SJakub Kruzik . pc - the preconditioner context 278268c6673SJakub Kruzik 279268c6673SJakub Kruzik Output Parameter: 280268c6673SJakub Kruzik . apc - additional preconditioner 281268c6673SJakub Kruzik 282268c6673SJakub Kruzik Level: advanced 283268c6673SJakub Kruzik 284268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION 285268c6673SJakub Kruzik @*/ 286268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc) 287268c6673SJakub Kruzik { 288268c6673SJakub Kruzik PetscErrorCode ierr; 289268c6673SJakub Kruzik 290268c6673SJakub Kruzik PetscFunctionBegin; 291268c6673SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 292268c6673SJakub Kruzik PetscValidPointer(pc,2); 293268c6673SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr); 294268c6673SJakub Kruzik PetscFunctionReturn(0); 295268c6673SJakub Kruzik } 296268c6673SJakub Kruzik 29722b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc) 29822b0793eSJakub Kruzik { 29922b0793eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 30022b0793eSJakub Kruzik PetscErrorCode ierr; 30122b0793eSJakub Kruzik 30222b0793eSJakub Kruzik PetscFunctionBegin; 30322b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 30422b0793eSJakub Kruzik def->pc = apc; 30522b0793eSJakub Kruzik ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr); 30622b0793eSJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr); 30722b0793eSJakub Kruzik PetscFunctionReturn(0); 30822b0793eSJakub Kruzik } 30922b0793eSJakub Kruzik 31022b0793eSJakub Kruzik /*@ 31122b0793eSJakub Kruzik PCDeflationSetPC - Set additional preconditioner. 31222b0793eSJakub Kruzik 313268c6673SJakub Kruzik Collective on PC 31422b0793eSJakub Kruzik 31522b0793eSJakub Kruzik Input Parameters: 31622b0793eSJakub Kruzik + pc - the preconditioner context 31722b0793eSJakub Kruzik - apc - additional preconditioner 31822b0793eSJakub Kruzik 319268c6673SJakub Kruzik Level: developer 32022b0793eSJakub Kruzik 321268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION 32222b0793eSJakub Kruzik @*/ 32322b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc) 32422b0793eSJakub Kruzik { 32522b0793eSJakub Kruzik PetscErrorCode ierr; 32622b0793eSJakub Kruzik 32722b0793eSJakub Kruzik PetscFunctionBegin; 32822b0793eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 32922b0793eSJakub Kruzik PetscValidHeaderSpecific(apc,PC_CLASSID,2); 33022b0793eSJakub Kruzik PetscCheckSameComm(pc,1,apc,2); 33122b0793eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr); 33222b0793eSJakub Kruzik PetscFunctionReturn(0); 33322b0793eSJakub Kruzik } 33422b0793eSJakub Kruzik 3354a99276eSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp) 3364a99276eSJakub Kruzik { 3374a99276eSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 3384a99276eSJakub Kruzik 3394a99276eSJakub Kruzik PetscFunctionBegin; 3404a99276eSJakub Kruzik *ksp = def->WtAWinv; 3414a99276eSJakub Kruzik PetscFunctionReturn(0); 3424a99276eSJakub Kruzik } 3434a99276eSJakub Kruzik 3444a99276eSJakub Kruzik /*@ 3454a99276eSJakub Kruzik PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP. 3464a99276eSJakub Kruzik 3474a99276eSJakub Kruzik Not Collective 3484a99276eSJakub Kruzik 3494a99276eSJakub Kruzik Input Parameters: 3504a99276eSJakub Kruzik . pc - preconditioner context 3514a99276eSJakub Kruzik 3524a99276eSJakub Kruzik Output Parameter: 3534a99276eSJakub Kruzik . ksp - coarse problem KSP context 3544a99276eSJakub Kruzik 3554a99276eSJakub Kruzik Level: developer 3564a99276eSJakub Kruzik 357f3eaa4f0SJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION 3584a99276eSJakub Kruzik @*/ 3594a99276eSJakub Kruzik PetscErrorCode PCDeflationGetCoarseKSP(PC pc,KSP *ksp) 3604a99276eSJakub Kruzik { 3614a99276eSJakub Kruzik PetscErrorCode ierr; 3624a99276eSJakub Kruzik 3634a99276eSJakub Kruzik PetscFunctionBegin; 3644a99276eSJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3654a99276eSJakub Kruzik PetscValidPointer(ksp,2); 3664a99276eSJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 3674a99276eSJakub Kruzik PetscFunctionReturn(0); 3684a99276eSJakub Kruzik } 3694a99276eSJakub Kruzik 370f3eaa4f0SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp) 371f3eaa4f0SJakub Kruzik { 372f3eaa4f0SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 373f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 374f3eaa4f0SJakub Kruzik 375f3eaa4f0SJakub Kruzik PetscFunctionBegin; 376f3eaa4f0SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 377f3eaa4f0SJakub Kruzik def->WtAWinv = ksp; 378f3eaa4f0SJakub Kruzik ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 379f3eaa4f0SJakub Kruzik ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr); 380f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 381f3eaa4f0SJakub Kruzik } 382f3eaa4f0SJakub Kruzik 383f3eaa4f0SJakub Kruzik /*@ 384f3eaa4f0SJakub Kruzik PCDeflationSetCoarseKSP - Set coarse problem KSP. 385f3eaa4f0SJakub Kruzik 386f3eaa4f0SJakub Kruzik Collective on PC 387f3eaa4f0SJakub Kruzik 388f3eaa4f0SJakub Kruzik Input Parameters: 389f3eaa4f0SJakub Kruzik + pc - preconditioner context 390f3eaa4f0SJakub Kruzik - ksp - coarse problem KSP context 391f3eaa4f0SJakub Kruzik 392f3eaa4f0SJakub Kruzik Level: developer 393f3eaa4f0SJakub Kruzik 394f3eaa4f0SJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION 395f3eaa4f0SJakub Kruzik @*/ 396f3eaa4f0SJakub Kruzik PetscErrorCode PCDeflationSetCoarseKSP(PC pc,KSP ksp) 397f3eaa4f0SJakub Kruzik { 398f3eaa4f0SJakub Kruzik PetscErrorCode ierr; 399f3eaa4f0SJakub Kruzik 400f3eaa4f0SJakub Kruzik PetscFunctionBegin; 401f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 402f3eaa4f0SJakub Kruzik PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 403f3eaa4f0SJakub Kruzik PetscCheckSameComm(pc,1,ksp,2); 404f3eaa4f0SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr); 405f3eaa4f0SJakub Kruzik PetscFunctionReturn(0); 406f3eaa4f0SJakub Kruzik } 407f3eaa4f0SJakub Kruzik 40837eeb815SJakub Kruzik /* 409b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 410b27c8092SJakub Kruzik */ 411b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 412b27c8092SJakub Kruzik { 413b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 414b27c8092SJakub Kruzik Mat A; 415b27c8092SJakub Kruzik Vec r,w1,w2; 416b27c8092SJakub Kruzik PetscBool nonzero; 417b27c8092SJakub Kruzik PetscErrorCode ierr; 41837eeb815SJakub Kruzik 419b27c8092SJakub Kruzik PetscFunctionBegin; 420b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 421b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 422b27c8092SJakub Kruzik r = def->work; 423b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 42437eeb815SJakub Kruzik 425b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 426b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 427b27c8092SJakub Kruzik if (nonzero) { 428b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 429b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 430b27c8092SJakub Kruzik } else { 431b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 432b27c8092SJakub Kruzik } 433b27c8092SJakub Kruzik 434a3177931SJakub Kruzik if (def->Wt) { 435a3177931SJakub Kruzik ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 436a3177931SJakub Kruzik } else { 437a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 438a3177931SJakub Kruzik } 439b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 440b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 441b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 442b27c8092SJakub Kruzik PetscFunctionReturn(0); 443b27c8092SJakub Kruzik } 44437eeb815SJakub Kruzik 445f8bf229cSJakub Kruzik /* 446f8bf229cSJakub Kruzik if (def->correct) { 447ae029463SJakub Kruzik z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r 448f8bf229cSJakub Kruzik } else { 449ae029463SJakub Kruzik z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r 450f8bf229cSJakub Kruzik } 45137eeb815SJakub Kruzik */ 452f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z) 453f8bf229cSJakub Kruzik { 454f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 455f8bf229cSJakub Kruzik Mat A; 456f8bf229cSJakub Kruzik Vec u,w1,w2; 457f8bf229cSJakub Kruzik PetscErrorCode ierr; 458f8bf229cSJakub Kruzik 459f8bf229cSJakub Kruzik PetscFunctionBegin; 460f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 461f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 462f8bf229cSJakub Kruzik u = def->work; 463f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 464f8bf229cSJakub Kruzik 465ae029463SJakub Kruzik ierr = PCApply(def->pc,r,z);CHKERRQ(ierr); /* z <- M^{-1}*r */ 46622b0793eSJakub Kruzik if (!def->init) { 467ae029463SJakub Kruzik ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr); /* w1 <- W'*A*z */ 468f8bf229cSJakub Kruzik if (def->correct) { 469ae029463SJakub Kruzik if (def->Wt) { 470ae029463SJakub Kruzik ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 471ae029463SJakub Kruzik } else { 472a3177931SJakub Kruzik ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*r */ 473f8bf229cSJakub Kruzik } 474ae029463SJakub Kruzik ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr); /* w1 <- w1 - l*w2 */ 475f8bf229cSJakub Kruzik } 476f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 477f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 47822b0793eSJakub Kruzik ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr); /* z <- z - u */ 47922b0793eSJakub Kruzik } 480f8bf229cSJakub Kruzik PetscFunctionReturn(0); 481f8bf229cSJakub Kruzik } 482f8bf229cSJakub Kruzik 48337eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 48437eeb815SJakub Kruzik { 485409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 486409a357bSJakub Kruzik KSP innerksp; 487409a357bSJakub Kruzik PC pcinner; 488409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 489409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 490409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 491409a357bSJakub Kruzik MatCompositeType ctype; 492409a357bSJakub Kruzik MPI_Comm comm; 493409a357bSJakub Kruzik const char *prefix; 49437eeb815SJakub Kruzik PetscErrorCode ierr; 49537eeb815SJakub Kruzik 49637eeb815SJakub Kruzik PetscFunctionBegin; 497409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 498409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 499f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 500f8bf229cSJakub Kruzik 501f8bf229cSJakub Kruzik /* compute a deflation space */ 502409a357bSJakub Kruzik if (def->W || def->Wt) { 503409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 504409a357bSJakub Kruzik } else { 505e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 50637eeb815SJakub Kruzik } 50737eeb815SJakub Kruzik 508409a357bSJakub Kruzik /* nested deflation */ 509409a357bSJakub Kruzik if (def->W) { 510409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 511409a357bSJakub Kruzik if (match) { 512409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 513409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 51437eeb815SJakub Kruzik } 515409a357bSJakub Kruzik } else { 516a3177931SJakub Kruzik ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr); 517409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 518409a357bSJakub Kruzik if (match) { 519409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 520409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 521409a357bSJakub Kruzik } 522409a357bSJakub Kruzik transp = PETSC_TRUE; 523409a357bSJakub Kruzik } 524409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 525409a357bSJakub Kruzik if (!transp) { 52628da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 52728da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 528409a357bSJakub Kruzik for (i=0; i<size; i++) { 529409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 530409a357bSJakub Kruzik } 531409a357bSJakub Kruzik size -= 1; 532409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 533409a357bSJakub Kruzik def->W = mats[size]; 534409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 535409a357bSJakub Kruzik if (size > 1) { 536409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 537409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 538409a357bSJakub Kruzik } else { 539409a357bSJakub Kruzik nextDef = mats[0]; 540409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 541409a357bSJakub Kruzik } 54228da0170SJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 543409a357bSJakub Kruzik } else { 544409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 545409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 546409a357bSJakub Kruzik } 54728da0170SJakub Kruzik } else { 54828da0170SJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 54928da0170SJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 55028da0170SJakub Kruzik for (i=0; i<size; i++) { 55128da0170SJakub Kruzik ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr); 55228da0170SJakub Kruzik } 55328da0170SJakub Kruzik size -= 1; 55428da0170SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 55528da0170SJakub Kruzik def->Wt = mats[0]; 55628da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 55728da0170SJakub Kruzik if (size > 1) { 55828da0170SJakub Kruzik ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr); 55928da0170SJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 56028da0170SJakub Kruzik } else { 56128da0170SJakub Kruzik nextDef = mats[1]; 56228da0170SJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr); 563409a357bSJakub Kruzik } 564409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 56528da0170SJakub Kruzik } else { 56628da0170SJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 56728da0170SJakub Kruzik ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr); 56828da0170SJakub Kruzik } 56928da0170SJakub Kruzik } 57028da0170SJakub Kruzik } 57128da0170SJakub Kruzik 57228da0170SJakub Kruzik if (transp) { 57328da0170SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 574a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr); 575409a357bSJakub Kruzik } 576409a357bSJakub Kruzik 57722b0793eSJakub Kruzik 57822b0793eSJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 57922b0793eSJakub Kruzik 580ae029463SJakub Kruzik /* assemble WtA */ 581ae029463SJakub Kruzik if (!def->WtA) { 582ae029463SJakub Kruzik if (def->Wt) { 583ae029463SJakub Kruzik ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 584ae029463SJakub Kruzik } else { 585a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX) 586a3177931SJakub Kruzik ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr); 587a3177931SJakub Kruzik ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 588a3177931SJakub Kruzik #else 589ae029463SJakub Kruzik ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr); 590a3177931SJakub Kruzik #endif 591ae029463SJakub Kruzik } 592ae029463SJakub Kruzik } 593409a357bSJakub Kruzik /* setup coarse problem */ 594409a357bSJakub Kruzik if (!def->WtAWinv) { 59528da0170SJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); 596409a357bSJakub Kruzik if (!def->WtAW) { 597ae029463SJakub Kruzik ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 598409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 599409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 600409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 601409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 602ae029463SJakub Kruzik /* Check columns of W are not in kernel of A */ 603409a357bSJakub Kruzik PetscReal *norms; 604409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 605409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 606409a357bSJakub Kruzik for (i=0; i<m; i++) { 607409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 608409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 609409a357bSJakub Kruzik } 610409a357bSJakub Kruzik } 611409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 612409a357bSJakub Kruzik #endif 613409a357bSJakub Kruzik } else { 614409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 615409a357bSJakub Kruzik } 616409a357bSJakub Kruzik /* TODO use MATINV */ 617409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 618409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 619409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 620557a2f7dSJakub Kruzik /* Setup KSP and PC */ 621557a2f7dSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 622557a2f7dSJakub Kruzik innerksp = def->WtAWinv; 623557a2f7dSJakub Kruzik /* set default KSPtype */ 624557a2f7dSJakub Kruzik if (!def->ksptype) { 625557a2f7dSJakub Kruzik def->ksptype = KSPFGMRES; 626557a2f7dSJakub Kruzik if (flgspd) { /* SPD system */ 627557a2f7dSJakub Kruzik def->ksptype = KSPFCG; 628557a2f7dSJakub Kruzik } 629557a2f7dSJakub Kruzik } 630ae029463SJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */ 631557a2f7dSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 632557a2f7dSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 633557a2f7dSJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 634ae029463SJakub Kruzik /* inherit options */ 635557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype; 636557a2f7dSJakub Kruzik ((PC_Deflation*)(pcinner->data))->correct = def->correct; 637557a2f7dSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 638557a2f7dSJakub Kruzik } else { /* the last level */ 639557a2f7dSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 640409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 641ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 642ac66f006SJakub Kruzik if (prefix) { 643ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 644ac66f006SJakub Kruzik } 64522b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 646ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 647409a357bSJakub Kruzik /* Reduction factor choice */ 648409a357bSJakub Kruzik red = def->reductionfact; 649409a357bSJakub Kruzik if (red < 0) { 650409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 651409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 652409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 653409a357bSJakub Kruzik if (match) red = commsize; 654409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 655409a357bSJakub Kruzik } 656409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 657ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 658409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 659ac66f006SJakub Kruzik if (innerksp) { 660409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 661409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 662ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 663409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 664409a357bSJakub Kruzik if (commsize == red) { 665ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 666409a357bSJakub Kruzik } else { 667ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 668409a357bSJakub Kruzik } 669409a357bSJakub Kruzik } 670557a2f7dSJakub Kruzik } 671557a2f7dSJakub Kruzik 672557a2f7dSJakub Kruzik if (innerksp) { 673409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 674409a357bSJakub Kruzik if (prefix) { 675409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 676409a357bSJakub Kruzik } 67722b0793eSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 678557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr); 679557a2f7dSJakub Kruzik ierr = KSPSetUp(innerksp);CHKERRQ(ierr); 680ac66f006SJakub Kruzik } 681409a357bSJakub Kruzik } 682557a2f7dSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 683557a2f7dSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 684409a357bSJakub Kruzik 68522b0793eSJakub Kruzik if (!def->pc) { 68622b0793eSJakub Kruzik ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr); 68722b0793eSJakub Kruzik ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr); 68822b0793eSJakub Kruzik ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr); 68922b0793eSJakub Kruzik if (prefix) { 69022b0793eSJakub Kruzik ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr); 69122b0793eSJakub Kruzik } 69222b0793eSJakub Kruzik ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr); 69322b0793eSJakub Kruzik ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr); 69422b0793eSJakub Kruzik ierr = PCSetUp(def->pc);CHKERRQ(ierr); 69522b0793eSJakub Kruzik } 69622b0793eSJakub Kruzik 697f8bf229cSJakub Kruzik /* create work vecs */ 698f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 699f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 70037eeb815SJakub Kruzik PetscFunctionReturn(0); 70137eeb815SJakub Kruzik } 702b30d4775SJakub Kruzik 70337eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 70437eeb815SJakub Kruzik { 705b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 70637eeb815SJakub Kruzik PetscErrorCode ierr; 70737eeb815SJakub Kruzik 70837eeb815SJakub Kruzik PetscFunctionBegin; 709b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 710b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 711b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 712b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 713ae029463SJakub Kruzik ierr = MatDestroy(&def->WtA);CHKERRQ(ierr); 714b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 715b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 71622b0793eSJakub Kruzik ierr = PCDestroy(&def->pc);CHKERRQ(ierr); 71737eeb815SJakub Kruzik PetscFunctionReturn(0); 71837eeb815SJakub Kruzik } 71937eeb815SJakub Kruzik 72037eeb815SJakub Kruzik /* 72137eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 72237eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 72337eeb815SJakub Kruzik 72437eeb815SJakub Kruzik Input Parameter: 72537eeb815SJakub Kruzik . pc - the preconditioner context 72637eeb815SJakub Kruzik 72737eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 72837eeb815SJakub Kruzik */ 72937eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 73037eeb815SJakub Kruzik { 73137eeb815SJakub Kruzik PetscErrorCode ierr; 73237eeb815SJakub Kruzik 73337eeb815SJakub Kruzik PetscFunctionBegin; 73437eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 73537eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 73637eeb815SJakub Kruzik PetscFunctionReturn(0); 73737eeb815SJakub Kruzik } 73837eeb815SJakub Kruzik 7398513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 7408513960bSJakub Kruzik { 7418513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 7428513960bSJakub Kruzik PetscBool iascii; 7438513960bSJakub Kruzik PetscErrorCode ierr; 7448513960bSJakub Kruzik 7458513960bSJakub Kruzik PetscFunctionBegin; 7468513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7478513960bSJakub Kruzik if (iascii) { 7488513960bSJakub Kruzik //if (cg->adaptiveconv) { 7498513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 7508513960bSJakub Kruzik //} 7518513960bSJakub Kruzik if (def->correct) { 7528513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 7538513960bSJakub Kruzik } 7548513960bSJakub Kruzik if (!def->nestedlvl) { 7558513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 7568513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 7578513960bSJakub Kruzik } 7588513960bSJakub Kruzik 75922b0793eSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr); 76022b0793eSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 76122b0793eSJakub Kruzik ierr = PCView(def->pc,viewer);CHKERRQ(ierr); 76222b0793eSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 76322b0793eSJakub Kruzik 7648513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 7658513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7668513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 7678513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7688513960bSJakub Kruzik } 7698513960bSJakub Kruzik PetscFunctionReturn(0); 7708513960bSJakub Kruzik } 7718513960bSJakub Kruzik 77237eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 77337eeb815SJakub Kruzik { 774880d05ceSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 77537eeb815SJakub Kruzik PetscErrorCode ierr; 77637eeb815SJakub Kruzik 77737eeb815SJakub Kruzik PetscFunctionBegin; 77837eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 779880d05ceSJakub Kruzik ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr); 780880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr); 781880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr); 782880d05ceSJakub Kruzik //TODO add set function and fix manpages 783880d05ceSJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr); 784fcb31d99SJakub Kruzik ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr); 785fcb31d99SJakub Kruzik ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr); 786880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr); 787880d05ceSJakub Kruzik ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr); 78837eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 78937eeb815SJakub Kruzik PetscFunctionReturn(0); 79037eeb815SJakub Kruzik } 79137eeb815SJakub Kruzik 79237eeb815SJakub Kruzik /*MC 793e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 794e361f8b5SJakub Kruzik or to a predefined value 79537eeb815SJakub Kruzik 79637eeb815SJakub Kruzik Options Database Key: 797e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 79837eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 79937eeb815SJakub Kruzik 80037eeb815SJakub Kruzik Level: beginner 80137eeb815SJakub Kruzik 80237eeb815SJakub Kruzik Notes: 803e361f8b5SJakub Kruzik todo 80437eeb815SJakub Kruzik 80537eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 806e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 80737eeb815SJakub Kruzik M*/ 80837eeb815SJakub Kruzik 80937eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 81037eeb815SJakub Kruzik { 81137eeb815SJakub Kruzik PC_Deflation *def; 81237eeb815SJakub Kruzik PetscErrorCode ierr; 81337eeb815SJakub Kruzik 81437eeb815SJakub Kruzik PetscFunctionBegin; 81537eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 81637eeb815SJakub Kruzik pc->data = (void*)def; 81737eeb815SJakub Kruzik 818e662bc50SJakub Kruzik def->init = PETSC_FALSE; 819e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 820fcb31d99SJakub Kruzik def->correctfact = 1.0; 821e662bc50SJakub Kruzik def->reductionfact = -1; 822e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 823e662bc50SJakub Kruzik def->spacesize = 1; 824e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 825e662bc50SJakub Kruzik def->nestedlvl = 0; 826e662bc50SJakub Kruzik def->maxnestedlvl = 0; 82737eeb815SJakub Kruzik 82837eeb815SJakub Kruzik /* 82937eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 83037eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 83137eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 83237eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 83337eeb815SJakub Kruzik not needed. 83437eeb815SJakub Kruzik */ 83537eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 83637eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 837b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 838b27c8092SJakub Kruzik pc->ops->postsolve = 0; 83937eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 84037eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 84137eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 84237eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 8438513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 84437eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 84537eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 84637eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 84737eeb815SJakub Kruzik 848*39417d7eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr); 849e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 8503cf3a049SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr); 8515c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 852268c6673SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr); 85322b0793eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr); 854e924b002SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr); 8554a99276eSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr); 856f3eaa4f0SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr); 85737eeb815SJakub Kruzik PetscFunctionReturn(0); 85837eeb815SJakub Kruzik } 85937eeb815SJakub Kruzik 860