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 5337eeb815SJakub Kruzik const char *const PCDeflationTypes[] = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0}; 5437eeb815SJakub Kruzik 55*8513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = { 56*8513960bSJakub Kruzik "haar", 57*8513960bSJakub Kruzik "jacket-haar", 58*8513960bSJakub Kruzik "db2", 59*8513960bSJakub Kruzik "db4", 60*8513960bSJakub Kruzik "db8", 61*8513960bSJakub Kruzik "db16", 62*8513960bSJakub Kruzik "biorth22", 63*8513960bSJakub Kruzik "meyer", 64*8513960bSJakub Kruzik "aggregation", 65*8513960bSJakub Kruzik "slepc", 66*8513960bSJakub Kruzik "slepc-cheap", 67*8513960bSJakub Kruzik "user", 68*8513960bSJakub Kruzik "DdefSpaceType", 69*8513960bSJakub Kruzik "Ddef_SPACE_", 70*8513960bSJakub Kruzik 0 71*8513960bSJakub Kruzik }; 72*8513960bSJakub Kruzik 73*8513960bSJakub Kruzik 74e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose) 75e662bc50SJakub Kruzik { 76e662bc50SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 77e662bc50SJakub Kruzik PetscErrorCode ierr; 78e662bc50SJakub Kruzik 79e662bc50SJakub Kruzik PetscFunctionBegin; 80e662bc50SJakub Kruzik if (transpose) { 81e662bc50SJakub Kruzik def->Wt = W; 82e662bc50SJakub Kruzik def->W = NULL; 83e662bc50SJakub Kruzik } else { 84e662bc50SJakub Kruzik def->W = W; 85e662bc50SJakub Kruzik } 86e662bc50SJakub Kruzik ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr); 87e662bc50SJakub Kruzik PetscFunctionReturn(0); 88e662bc50SJakub Kruzik } 89e662bc50SJakub Kruzik 90e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */ 91e662bc50SJakub Kruzik /*@ 92e662bc50SJakub Kruzik PCDeflationSetSpace - Set deflation space matrix (or its transpose). 93e662bc50SJakub Kruzik 94e662bc50SJakub Kruzik Logically Collective on PC 95e662bc50SJakub Kruzik 96e662bc50SJakub Kruzik Input Parameters: 97e662bc50SJakub Kruzik + pc - the preconditioner context 98e662bc50SJakub Kruzik . W - deflation matrix 99e662bc50SJakub Kruzik - tranpose - indicates that W is an explicit transpose of the deflation matrix 100e662bc50SJakub Kruzik 101e662bc50SJakub Kruzik Level: intermediate 102e662bc50SJakub Kruzik 103e662bc50SJakub Kruzik .seealso: PCDEFLATION 104e662bc50SJakub Kruzik @*/ 105e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose) 106e662bc50SJakub Kruzik { 107e662bc50SJakub Kruzik PetscErrorCode ierr; 108e662bc50SJakub Kruzik 109e662bc50SJakub Kruzik PetscFunctionBegin; 110e662bc50SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 111e662bc50SJakub Kruzik PetscValidHeaderSpecific(W,MAT_CLASSID,2); 112e662bc50SJakub Kruzik PetscValidLogicalCollectiveBool(pc,transpose,3); 113e662bc50SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr); 114e662bc50SJakub Kruzik PetscFunctionReturn(0); 115e662bc50SJakub Kruzik } 116e662bc50SJakub Kruzik 1175c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max) 1185c0c31c5SJakub Kruzik { 1195c0c31c5SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 1205c0c31c5SJakub Kruzik 1215c0c31c5SJakub Kruzik PetscFunctionBegin; 1225c0c31c5SJakub Kruzik def->nestedlvl = current; 1235c0c31c5SJakub Kruzik def->maxnestedlvl = max; 1245c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1255c0c31c5SJakub Kruzik } 1265c0c31c5SJakub Kruzik 1275c0c31c5SJakub Kruzik /*@ 1285c0c31c5SJakub Kruzik PCDeflationSetMaxLvl - Set maximum level of deflation. 1295c0c31c5SJakub Kruzik 1305c0c31c5SJakub Kruzik Logically Collective on PC 1315c0c31c5SJakub Kruzik 1325c0c31c5SJakub Kruzik Input Parameters: 1335c0c31c5SJakub Kruzik + pc - the preconditioner context 1345c0c31c5SJakub Kruzik . max - maximum deflation level 1355c0c31c5SJakub Kruzik 1365c0c31c5SJakub Kruzik Level: intermediate 1375c0c31c5SJakub Kruzik 1385c0c31c5SJakub Kruzik .seealso: PCDEFLATION 1395c0c31c5SJakub Kruzik @*/ 1405c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max) 1415c0c31c5SJakub Kruzik { 1425c0c31c5SJakub Kruzik PetscErrorCode ierr; 1435c0c31c5SJakub Kruzik 1445c0c31c5SJakub Kruzik PetscFunctionBegin; 1455c0c31c5SJakub Kruzik PetscValidLogicalCollectiveInt(pc,max,2); 1465c0c31c5SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr); 1475c0c31c5SJakub Kruzik PetscFunctionReturn(0); 1485c0c31c5SJakub Kruzik } 149e662bc50SJakub Kruzik 15037eeb815SJakub Kruzik /* 151b27c8092SJakub Kruzik TODO CP corection? 152b27c8092SJakub Kruzik x <- x + W*(W'*A*W)^{-1}*W'*r = x + Q*r 153b27c8092SJakub Kruzik */ 154b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x) 155b27c8092SJakub Kruzik { 156b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 157b27c8092SJakub Kruzik Mat A; 158b27c8092SJakub Kruzik Vec r,w1,w2; 159b27c8092SJakub Kruzik PetscBool nonzero; 160b27c8092SJakub Kruzik PetscErrorCode ierr; 16137eeb815SJakub Kruzik 162b27c8092SJakub Kruzik PetscFunctionBegin; 163b27c8092SJakub Kruzik w1 = def->workcoarse[0]; 164b27c8092SJakub Kruzik w2 = def->workcoarse[1]; 165b27c8092SJakub Kruzik r = def->work; 166b27c8092SJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 16737eeb815SJakub Kruzik 168b27c8092SJakub Kruzik ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 169b27c8092SJakub Kruzik ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 170b27c8092SJakub Kruzik if (nonzero) { 171b27c8092SJakub Kruzik ierr = MatMult(A,x,r);CHKERRQ(ierr); /* r <- b - Ax */ 172b27c8092SJakub Kruzik ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 173b27c8092SJakub Kruzik } else { 174b27c8092SJakub Kruzik ierr = VecCopy(b,r);CHKERRQ(ierr); /* r <- b (x is 0) */ 175b27c8092SJakub Kruzik } 176b27c8092SJakub Kruzik 177b27c8092SJakub Kruzik ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr); /* w1 <- W'*r */ 178b27c8092SJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 179b27c8092SJakub Kruzik ierr = MatMult(def->W,w2,r);CHKERRQ(ierr); /* r <- W*w2 */ 180b27c8092SJakub Kruzik ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr); 181b27c8092SJakub Kruzik PetscFunctionReturn(0); 182b27c8092SJakub Kruzik } 18337eeb815SJakub Kruzik 184f8bf229cSJakub Kruzik /* 185f8bf229cSJakub Kruzik if (def->correct) { 186f8bf229cSJakub Kruzik z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*r 187f8bf229cSJakub Kruzik } else { 188f8bf229cSJakub Kruzik z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*r 189f8bf229cSJakub Kruzik } 19037eeb815SJakub Kruzik */ 191f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r, Vec z) 192f8bf229cSJakub Kruzik { 193f8bf229cSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 194f8bf229cSJakub Kruzik Mat A; 195f8bf229cSJakub Kruzik Vec u,w1,w2; 196f8bf229cSJakub Kruzik PetscErrorCode ierr; 197f8bf229cSJakub Kruzik 198f8bf229cSJakub Kruzik PetscFunctionBegin; 199f8bf229cSJakub Kruzik w1 = def->workcoarse[0]; 200f8bf229cSJakub Kruzik w2 = def->workcoarse[1]; 201f8bf229cSJakub Kruzik u = def->work; 202f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 203f8bf229cSJakub Kruzik 204f8bf229cSJakub Kruzik if (!def->AW) { 205f8bf229cSJakub Kruzik ierr = MatMult(A,r,u);CHKERRQ(ierr); /* u <- A*r */ 206f8bf229cSJakub Kruzik if (def->correct) ierr = VecAXPY(u,-1.0,r);CHKERRQ(ierr); /* u <- A*r -r */ 207f8bf229cSJakub Kruzik ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr); /* w1 <- W'*u */ 208f8bf229cSJakub Kruzik } else { 209f8bf229cSJakub Kruzik ierr = MatMultTranspose(def->AW,u,w1);CHKERRQ(ierr); /* u <- A*r */ 210f8bf229cSJakub Kruzik if (def->correct) { 211f8bf229cSJakub Kruzik ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr); /* w2 <- W'*u */ 212f8bf229cSJakub Kruzik ierr = VecAXPY(w1,-1.0,w2);CHKERRQ(ierr); /* w1 <- w1 - w2 */ 213f8bf229cSJakub Kruzik } 214f8bf229cSJakub Kruzik } 215f8bf229cSJakub Kruzik ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr); /* w2 <- (W'*A*W)^{-1}*w1 */ 216f8bf229cSJakub Kruzik ierr = MatMult(def->W,w2,u);CHKERRQ(ierr); /* u <- W*w2 */ 217f8bf229cSJakub Kruzik ierr = VecWAXPY(z,-1.0,u,r);CHKERRQ(ierr); /* z <- r - u */ 218f8bf229cSJakub Kruzik PetscFunctionReturn(0); 219f8bf229cSJakub Kruzik } 220f8bf229cSJakub Kruzik 221b27c8092SJakub Kruzik static PetscErrorCode PCDeflationSetType_Deflation(PC pc,PCDeflationType type) 222b27c8092SJakub Kruzik { 223b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 224b27c8092SJakub Kruzik 225b27c8092SJakub Kruzik PetscFunctionBegin; 226b27c8092SJakub Kruzik def->init = PETSC_FALSE; 227b27c8092SJakub Kruzik def->pre = PETSC_FALSE; 228b27c8092SJakub Kruzik if (type == PC_DEFLATION_POST) { 229b27c8092SJakub Kruzik //pc->ops->postsolve = PCPostSolve_Deflation; 230b27c8092SJakub Kruzik pc->ops->presolve = 0; 231b27c8092SJakub Kruzik } else { 232b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 233b27c8092SJakub Kruzik pc->ops->postsolve = 0; 234b27c8092SJakub Kruzik if (type == PC_DEFLATION_INIT) { 235b27c8092SJakub Kruzik def->init = PETSC_TRUE; 236b27c8092SJakub Kruzik pc->ops->apply = 0; 237b27c8092SJakub Kruzik } else { 238b27c8092SJakub Kruzik def->pre = PETSC_TRUE; 239b27c8092SJakub Kruzik } 240b27c8092SJakub Kruzik } 241b27c8092SJakub Kruzik PetscFunctionReturn(0); 242b27c8092SJakub Kruzik } 243b27c8092SJakub Kruzik 244b27c8092SJakub Kruzik /*@ 245b27c8092SJakub Kruzik PCDeflationSetType - Causes the deflation preconditioner to use only a special 246b27c8092SJakub Kruzik initial gues or pre/post solve solution update 247b27c8092SJakub Kruzik 248b27c8092SJakub Kruzik Logically Collective on PC 249b27c8092SJakub Kruzik 250b27c8092SJakub Kruzik Input Parameters: 251b27c8092SJakub Kruzik + pc - the preconditioner context 252b27c8092SJakub Kruzik - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 253b27c8092SJakub Kruzik 254b27c8092SJakub Kruzik Options Database Key: 255b27c8092SJakub Kruzik . -pc_deflation_type <pre,init,post> 256b27c8092SJakub Kruzik 257b27c8092SJakub Kruzik Level: intermediate 258b27c8092SJakub Kruzik 259b27c8092SJakub Kruzik Concepts: Deflation preconditioner 260b27c8092SJakub Kruzik 261b27c8092SJakub Kruzik .seealso: PCDeflationGetType() 262b27c8092SJakub Kruzik @*/ 263b27c8092SJakub Kruzik PetscErrorCode PCDeflationSetType(PC pc,PCDeflationType type) 264b27c8092SJakub Kruzik { 265b27c8092SJakub Kruzik PetscErrorCode ierr; 266b27c8092SJakub Kruzik 267b27c8092SJakub Kruzik PetscFunctionBegin; 268b27c8092SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 269b27c8092SJakub Kruzik ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr); 270b27c8092SJakub Kruzik PetscFunctionReturn(0); 271b27c8092SJakub Kruzik } 272b27c8092SJakub Kruzik 273b27c8092SJakub Kruzik static PetscErrorCode PCDeflationGetType_Deflation(PC pc,PCDeflationType *type) 274b27c8092SJakub Kruzik { 275b27c8092SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 276b27c8092SJakub Kruzik 277b27c8092SJakub Kruzik PetscFunctionBegin; 278b27c8092SJakub Kruzik if (def->init) { 279b27c8092SJakub Kruzik *type = PC_DEFLATION_INIT; 280b27c8092SJakub Kruzik } else if (def->pre) { 281b27c8092SJakub Kruzik *type = PC_DEFLATION_PRE; 282b27c8092SJakub Kruzik } else { 283b27c8092SJakub Kruzik *type = PC_DEFLATION_POST; 284b27c8092SJakub Kruzik } 285b27c8092SJakub Kruzik PetscFunctionReturn(0); 286b27c8092SJakub Kruzik } 287b27c8092SJakub Kruzik 288b27c8092SJakub Kruzik /*@ 289b27c8092SJakub Kruzik PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner 290b27c8092SJakub Kruzik 291b27c8092SJakub Kruzik Not Collective on PC 292b27c8092SJakub Kruzik 293b27c8092SJakub Kruzik Input Parameter: 294b27c8092SJakub Kruzik . pc - the preconditioner context 295b27c8092SJakub Kruzik 296b27c8092SJakub Kruzik Output Parameter: 297b27c8092SJakub Kruzik - type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST 298b27c8092SJakub Kruzik 299b27c8092SJakub Kruzik Level: intermediate 300b27c8092SJakub Kruzik 301b27c8092SJakub Kruzik Concepts: Deflation preconditioner 302b27c8092SJakub Kruzik 303b27c8092SJakub Kruzik .seealso: PCDeflationSetType() 304b27c8092SJakub Kruzik @*/ 305b27c8092SJakub Kruzik PetscErrorCode PCDeflationGetType(PC pc,PCDeflationType *type) 306b27c8092SJakub Kruzik { 307b27c8092SJakub Kruzik PetscErrorCode ierr; 308b27c8092SJakub Kruzik 309b27c8092SJakub Kruzik PetscFunctionBegin; 310b27c8092SJakub Kruzik PetscValidHeaderSpecific(pc,PC_CLASSID,1); 311b27c8092SJakub Kruzik ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr); 312b27c8092SJakub Kruzik PetscFunctionReturn(0); 313b27c8092SJakub Kruzik } 314b27c8092SJakub Kruzik 31537eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc) 31637eeb815SJakub Kruzik { 317409a357bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 318409a357bSJakub Kruzik KSP innerksp; 319409a357bSJakub Kruzik PC pcinner; 320409a357bSJakub Kruzik Mat Amat,nextDef=NULL,*mats; 321409a357bSJakub Kruzik PetscInt i,m,red,size,commsize; 322409a357bSJakub Kruzik PetscBool match,flgspd,transp=PETSC_FALSE; 323409a357bSJakub Kruzik MatCompositeType ctype; 324409a357bSJakub Kruzik MPI_Comm comm; 325409a357bSJakub Kruzik const char *prefix; 32637eeb815SJakub Kruzik PetscErrorCode ierr; 32737eeb815SJakub Kruzik 32837eeb815SJakub Kruzik PetscFunctionBegin; 329409a357bSJakub Kruzik if (pc->setupcalled) PetscFunctionReturn(0); 330409a357bSJakub Kruzik ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 331f8bf229cSJakub Kruzik ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr); 332f8bf229cSJakub Kruzik 333f8bf229cSJakub Kruzik /* compute a deflation space */ 334409a357bSJakub Kruzik if (def->W || def->Wt) { 335409a357bSJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_USER; 336409a357bSJakub Kruzik } else { 337e53e0a0dSJakub Kruzik ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr); 33837eeb815SJakub Kruzik } 33937eeb815SJakub Kruzik 340409a357bSJakub Kruzik /* nested deflation */ 341409a357bSJakub Kruzik if (def->W) { 342409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr); 343409a357bSJakub Kruzik if (match) { 344409a357bSJakub Kruzik ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr); 345409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr); 34637eeb815SJakub Kruzik } 347409a357bSJakub Kruzik } else { 348409a357bSJakub Kruzik ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr); 349409a357bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr); 350409a357bSJakub Kruzik if (match) { 351409a357bSJakub Kruzik ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr); 352409a357bSJakub Kruzik ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr); 353409a357bSJakub Kruzik } 354409a357bSJakub Kruzik transp = PETSC_TRUE; 355409a357bSJakub Kruzik } 356409a357bSJakub Kruzik if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) { 357409a357bSJakub Kruzik ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr); 358409a357bSJakub Kruzik if (!transp) { 359409a357bSJakub Kruzik for (i=0; i<size; i++) { 360409a357bSJakub Kruzik ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr); 361409a357bSJakub Kruzik //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr); 362409a357bSJakub Kruzik } 363409a357bSJakub Kruzik if (def->nestedlvl < def->maxnestedlvl) { 364409a357bSJakub Kruzik size -= 1; 365409a357bSJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 366409a357bSJakub Kruzik def->W = mats[size]; 367409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr); 368409a357bSJakub Kruzik if (size > 1) { 369409a357bSJakub Kruzik ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr); 370409a357bSJakub Kruzik ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 371409a357bSJakub Kruzik } else { 372409a357bSJakub Kruzik nextDef = mats[0]; 373409a357bSJakub Kruzik ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr); 374409a357bSJakub Kruzik } 375409a357bSJakub Kruzik } else { 376409a357bSJakub Kruzik /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */ 377409a357bSJakub Kruzik ierr = MatCompositeMerge(def->W);CHKERRQ(ierr); 378409a357bSJakub Kruzik } 379409a357bSJakub Kruzik } 380409a357bSJakub Kruzik ierr = PetscFree(mats);CHKERRQ(ierr); 381409a357bSJakub Kruzik } 382409a357bSJakub Kruzik 383409a357bSJakub Kruzik /* setup coarse problem */ 384409a357bSJakub Kruzik if (!def->WtAWinv) { 385409a357bSJakub Kruzik ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */ 386409a357bSJakub Kruzik if (!def->WtAW) { 387409a357bSJakub Kruzik /* TODO add implicit product version ? */ 388409a357bSJakub Kruzik if (!def->AW) { 389409a357bSJakub Kruzik ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 390409a357bSJakub Kruzik } else { 391409a357bSJakub Kruzik ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr); 392409a357bSJakub Kruzik } 393409a357bSJakub Kruzik /* TODO create MatInheritOption(Mat,MatOption) */ 394409a357bSJakub Kruzik ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr); 395409a357bSJakub Kruzik ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr); 396409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG) 397409a357bSJakub Kruzik /* Check WtAW is not sigular */ 398409a357bSJakub Kruzik PetscReal *norms; 399409a357bSJakub Kruzik ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr); 400409a357bSJakub Kruzik ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr); 401409a357bSJakub Kruzik for (i=0; i<m; i++) { 402409a357bSJakub Kruzik if (norms[i] < 100*PETSC_MACHINE_EPSILON) { 403409a357bSJakub Kruzik SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i); 404409a357bSJakub Kruzik } 405409a357bSJakub Kruzik } 406409a357bSJakub Kruzik ierr = PetscFree(norms);CHKERRQ(ierr); 407409a357bSJakub Kruzik #endif 408409a357bSJakub Kruzik } else { 409409a357bSJakub Kruzik ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr); 410409a357bSJakub Kruzik } 411409a357bSJakub Kruzik /* TODO use MATINV */ 412ac66f006SJakub Kruzik ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 413409a357bSJakub Kruzik ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr); 414409a357bSJakub Kruzik ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr); 415409a357bSJakub Kruzik ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr); 416409a357bSJakub Kruzik ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr); 417409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr); 418ac66f006SJakub Kruzik /* ugly hack to not have overwritten PCTELESCOPE */ 419ac66f006SJakub Kruzik if (prefix) { 420ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr); 421ac66f006SJakub Kruzik ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 422ac66f006SJakub Kruzik } else { 423ac66f006SJakub Kruzik ierr = KSPSetOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr); 424ac66f006SJakub Kruzik } 425ac66f006SJakub Kruzik ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr); 426409a357bSJakub Kruzik /* Reduction factor choice */ 427409a357bSJakub Kruzik red = def->reductionfact; 428409a357bSJakub Kruzik if (red < 0) { 429409a357bSJakub Kruzik ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); 430409a357bSJakub Kruzik red = ceil((float)commsize/ceil((float)m/commsize)); 431409a357bSJakub Kruzik ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr); 432409a357bSJakub Kruzik if (match) red = commsize; 433409a357bSJakub Kruzik ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */ 434409a357bSJakub Kruzik } 435409a357bSJakub Kruzik ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr); 436ac66f006SJakub Kruzik ierr = PCSetUp(pcinner);CHKERRQ(ierr); 437409a357bSJakub Kruzik ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr); 438ac66f006SJakub Kruzik if (innerksp) { 439409a357bSJakub Kruzik ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr); 440409a357bSJakub Kruzik /* Setup KSP and PC */ 441409a357bSJakub Kruzik if (nextDef) { /* next level for multilevel deflation */ 442409a357bSJakub Kruzik /* set default KSPtype */ 443409a357bSJakub Kruzik if (!def->ksptype) { 444409a357bSJakub Kruzik def->ksptype = KSPFGMRES; 445409a357bSJakub Kruzik if (flgspd) { /* SPD system */ 446409a357bSJakub Kruzik def->ksptype = KSPFCG; 447409a357bSJakub Kruzik } 448409a357bSJakub Kruzik } 449409a357bSJakub Kruzik ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */ 450409a357bSJakub Kruzik ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */ 451409a357bSJakub Kruzik ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr); 4525c0c31c5SJakub Kruzik ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr); 453409a357bSJakub Kruzik /* inherit options TODO if not set */ 454409a357bSJakub Kruzik ((PC_Deflation*)(pcinner))->ksptype = def->ksptype; 455409a357bSJakub Kruzik ((PC_Deflation*)(pcinner))->correct = def->correct; 456409a357bSJakub Kruzik ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv; 457409a357bSJakub Kruzik ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst; 458409a357bSJakub Kruzik ierr = MatDestroy(&nextDef);CHKERRQ(ierr); 459409a357bSJakub Kruzik } else { /* the last level */ 460ac66f006SJakub Kruzik ierr = KSPSetType(innerksp,KSPGMRES);CHKERRQ(ierr); 461409a357bSJakub Kruzik /* TODO Cholesky if flgspd? */ 462ac66f006SJakub Kruzik ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr); 463409a357bSJakub Kruzik //TODO remove explicit matSolverPackage 464409a357bSJakub Kruzik if (commsize == red) { 465ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr); 466409a357bSJakub Kruzik } else { 467ac66f006SJakub Kruzik ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); 468409a357bSJakub Kruzik } 469409a357bSJakub Kruzik } 470409a357bSJakub Kruzik /* TODO use def_[lvl]_ if lvl > 0? */ 471409a357bSJakub Kruzik if (prefix) { 472409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr); 473409a357bSJakub Kruzik ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 474409a357bSJakub Kruzik } else { 475409a357bSJakub Kruzik ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr); 476409a357bSJakub Kruzik } 477ac66f006SJakub Kruzik } 478409a357bSJakub Kruzik /* TODO: check if WtAWinv is KSP and move following from this if */ 479409a357bSJakub Kruzik ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr); 480409a357bSJakub Kruzik //if (def->adaptiveconv) { 481409a357bSJakub Kruzik // PetscReal *rnorm; 482409a357bSJakub Kruzik // PetscNew(&rnorm); 483409a357bSJakub Kruzik // ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr); 484409a357bSJakub Kruzik //} 485409a357bSJakub Kruzik ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr); 486409a357bSJakub Kruzik } 487409a357bSJakub Kruzik 488f8bf229cSJakub Kruzik /* create work vecs */ 489f8bf229cSJakub Kruzik ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr); 490f8bf229cSJakub Kruzik ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr); 49137eeb815SJakub Kruzik PetscFunctionReturn(0); 49237eeb815SJakub Kruzik } 493b30d4775SJakub Kruzik 49437eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc) 49537eeb815SJakub Kruzik { 496b30d4775SJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 49737eeb815SJakub Kruzik PetscErrorCode ierr; 49837eeb815SJakub Kruzik 49937eeb815SJakub Kruzik PetscFunctionBegin; 500b30d4775SJakub Kruzik ierr = VecDestroy(&def->work);CHKERRQ(ierr); 501b30d4775SJakub Kruzik ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr); 502b30d4775SJakub Kruzik ierr = MatDestroy(&def->W);CHKERRQ(ierr); 503b30d4775SJakub Kruzik ierr = MatDestroy(&def->Wt);CHKERRQ(ierr); 504b30d4775SJakub Kruzik ierr = MatDestroy(&def->AW);CHKERRQ(ierr); 505b30d4775SJakub Kruzik ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr); 506b30d4775SJakub Kruzik ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr); 50737eeb815SJakub Kruzik PetscFunctionReturn(0); 50837eeb815SJakub Kruzik } 50937eeb815SJakub Kruzik 51037eeb815SJakub Kruzik /* 51137eeb815SJakub Kruzik PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner 51237eeb815SJakub Kruzik that was created with PCCreate_Deflation(). 51337eeb815SJakub Kruzik 51437eeb815SJakub Kruzik Input Parameter: 51537eeb815SJakub Kruzik . pc - the preconditioner context 51637eeb815SJakub Kruzik 51737eeb815SJakub Kruzik Application Interface Routine: PCDestroy() 51837eeb815SJakub Kruzik */ 51937eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc) 52037eeb815SJakub Kruzik { 52137eeb815SJakub Kruzik PetscErrorCode ierr; 52237eeb815SJakub Kruzik 52337eeb815SJakub Kruzik PetscFunctionBegin; 52437eeb815SJakub Kruzik ierr = PCReset_Deflation(pc);CHKERRQ(ierr); 52537eeb815SJakub Kruzik ierr = PetscFree(pc->data);CHKERRQ(ierr); 52637eeb815SJakub Kruzik PetscFunctionReturn(0); 52737eeb815SJakub Kruzik } 52837eeb815SJakub Kruzik 529*8513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer) 530*8513960bSJakub Kruzik { 531*8513960bSJakub Kruzik PC_Deflation *def = (PC_Deflation*)pc->data; 532*8513960bSJakub Kruzik PetscBool iascii; 533*8513960bSJakub Kruzik PetscErrorCode ierr; 534*8513960bSJakub Kruzik 535*8513960bSJakub Kruzik PetscFunctionBegin; 536*8513960bSJakub Kruzik ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 537*8513960bSJakub Kruzik if (iascii) { 538*8513960bSJakub Kruzik //if (cg->adaptiveconv) { 539*8513960bSJakub Kruzik // ierr = PetscViewerASCIIPrintf(viewer," DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr); 540*8513960bSJakub Kruzik //} 541*8513960bSJakub Kruzik if (def->correct) { 542*8513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Using CP correction\n");CHKERRQ(ierr); 543*8513960bSJakub Kruzik } 544*8513960bSJakub Kruzik if (!def->nestedlvl) { 545*8513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr); 546*8513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer," DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr); 547*8513960bSJakub Kruzik } 548*8513960bSJakub Kruzik 549*8513960bSJakub Kruzik ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr); 550*8513960bSJakub Kruzik ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 551*8513960bSJakub Kruzik ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr); 552*8513960bSJakub Kruzik ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 553*8513960bSJakub Kruzik } 554*8513960bSJakub Kruzik PetscFunctionReturn(0); 555*8513960bSJakub Kruzik } 556*8513960bSJakub Kruzik 55737eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc) 55837eeb815SJakub Kruzik { 55937eeb815SJakub Kruzik PC_Deflation *jac = (PC_Deflation*)pc->data; 56037eeb815SJakub Kruzik PetscErrorCode ierr; 56137eeb815SJakub Kruzik PetscBool flg; 56237eeb815SJakub Kruzik PCDeflationType deflt,type; 56337eeb815SJakub Kruzik 56437eeb815SJakub Kruzik PetscFunctionBegin; 56537eeb815SJakub Kruzik ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr); 56637eeb815SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr); 56737eeb815SJakub Kruzik ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 56837eeb815SJakub Kruzik if (flg) { 56937eeb815SJakub Kruzik ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr); 57037eeb815SJakub Kruzik } 57137eeb815SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 57237eeb815SJakub Kruzik PetscFunctionReturn(0); 57337eeb815SJakub Kruzik } 57437eeb815SJakub Kruzik 57537eeb815SJakub Kruzik /*MC 576e361f8b5SJakub Kruzik PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates) 577e361f8b5SJakub Kruzik or to a predefined value 57837eeb815SJakub Kruzik 57937eeb815SJakub Kruzik Options Database Key: 580e361f8b5SJakub Kruzik + -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre) 58137eeb815SJakub Kruzik - -pc_jacobi_abs - use the absolute value of the diagonal entry 58237eeb815SJakub Kruzik 58337eeb815SJakub Kruzik Level: beginner 58437eeb815SJakub Kruzik 58537eeb815SJakub Kruzik Notes: 586e361f8b5SJakub Kruzik todo 58737eeb815SJakub Kruzik 58837eeb815SJakub Kruzik .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 589e662bc50SJakub Kruzik PCDeflationSetType(), PCDeflationSetSpace() 59037eeb815SJakub Kruzik M*/ 59137eeb815SJakub Kruzik 59237eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) 59337eeb815SJakub Kruzik { 59437eeb815SJakub Kruzik PC_Deflation *def; 59537eeb815SJakub Kruzik PetscErrorCode ierr; 59637eeb815SJakub Kruzik 59737eeb815SJakub Kruzik PetscFunctionBegin; 59837eeb815SJakub Kruzik ierr = PetscNewLog(pc,&def);CHKERRQ(ierr); 59937eeb815SJakub Kruzik pc->data = (void*)def; 60037eeb815SJakub Kruzik 601e662bc50SJakub Kruzik def->init = PETSC_FALSE; 602e662bc50SJakub Kruzik def->pre = PETSC_TRUE; 603e662bc50SJakub Kruzik def->correct = PETSC_FALSE; 604e662bc50SJakub Kruzik def->truenorm = PETSC_TRUE; 605e662bc50SJakub Kruzik def->reductionfact = -1; 606e662bc50SJakub Kruzik def->spacetype = PC_DEFLATION_SPACE_HAAR; 607e662bc50SJakub Kruzik def->spacesize = 1; 608e662bc50SJakub Kruzik def->extendsp = PETSC_FALSE; 609e662bc50SJakub Kruzik def->nestedlvl = 0; 610e662bc50SJakub Kruzik def->maxnestedlvl = 0; 611e662bc50SJakub Kruzik def->adaptiveconv = PETSC_FALSE; 612e662bc50SJakub Kruzik def->adaptiveconst = 1.0; 61337eeb815SJakub Kruzik 61437eeb815SJakub Kruzik /* 61537eeb815SJakub Kruzik Set the pointers for the functions that are provided above. 61637eeb815SJakub Kruzik Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 61737eeb815SJakub Kruzik are called, they will automatically call these functions. Note we 61837eeb815SJakub Kruzik choose not to provide a couple of these functions since they are 61937eeb815SJakub Kruzik not needed. 62037eeb815SJakub Kruzik */ 62137eeb815SJakub Kruzik pc->ops->apply = PCApply_Deflation; 62237eeb815SJakub Kruzik pc->ops->applytranspose = PCApply_Deflation; 623b27c8092SJakub Kruzik pc->ops->presolve = PCPreSolve_Deflation; 624b27c8092SJakub Kruzik pc->ops->postsolve = 0; 62537eeb815SJakub Kruzik pc->ops->setup = PCSetUp_Deflation; 62637eeb815SJakub Kruzik pc->ops->reset = PCReset_Deflation; 62737eeb815SJakub Kruzik pc->ops->destroy = PCDestroy_Deflation; 62837eeb815SJakub Kruzik pc->ops->setfromoptions = PCSetFromOptions_Deflation; 629*8513960bSJakub Kruzik pc->ops->view = PCView_Deflation; 63037eeb815SJakub Kruzik pc->ops->applyrichardson = 0; 63137eeb815SJakub Kruzik pc->ops->applysymmetricleft = 0; 63237eeb815SJakub Kruzik pc->ops->applysymmetricright = 0; 63337eeb815SJakub Kruzik 63437eeb815SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr); 63537eeb815SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr); 636e662bc50SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr); 6375c0c31c5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr); 63837eeb815SJakub Kruzik PetscFunctionReturn(0); 63937eeb815SJakub Kruzik } 64037eeb815SJakub Kruzik 641