xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 0a78af2228aa9a5baa3b5a8e3872e699fa117360)
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",
16*0a78af22SJakub Kruzik   "PCDeflationSpaceType",
17*0a78af22SJakub 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.
32a122ebaeSJakub Kruzik     Sets initial guess to the solution on the deflation space but do not apply deflation preconditioner.
33a122ebaeSJakub Kruzik     The additional preconditioner is still applied.
34a122ebaeSJakub Kruzik 
35a122ebaeSJakub Kruzik    Logically Collective on PC
36a122ebaeSJakub Kruzik 
37a122ebaeSJakub Kruzik    Input Parameters:
38a122ebaeSJakub Kruzik +  pc  - the preconditioner context
39a122ebaeSJakub Kruzik -  flg - default PETSC_FALSE
40a122ebaeSJakub Kruzik 
41a122ebaeSJakub Kruzik    Level: intermediate
42a122ebaeSJakub Kruzik 
43a122ebaeSJakub Kruzik .seealso: PCDEFLATION
44a122ebaeSJakub Kruzik @*/
45a122ebaeSJakub Kruzik PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg)
46a122ebaeSJakub Kruzik {
47a122ebaeSJakub Kruzik   PetscErrorCode ierr;
48a122ebaeSJakub Kruzik 
49a122ebaeSJakub Kruzik   PetscFunctionBegin;
50a122ebaeSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
51a122ebaeSJakub Kruzik   PetscValidLogicalCollectiveBool(pc,flg,2);
52a122ebaeSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
53a122ebaeSJakub Kruzik   PetscFunctionReturn(0);
54a122ebaeSJakub Kruzik }
55a122ebaeSJakub Kruzik 
56a122ebaeSJakub Kruzik 
5798d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
5898d6e3deSJakub Kruzik {
5998d6e3deSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
6098d6e3deSJakub Kruzik 
6198d6e3deSJakub Kruzik   PetscFunctionBegin;
6298d6e3deSJakub Kruzik   if (current) def->nestedlvl = current;
6398d6e3deSJakub Kruzik   def->maxnestedlvl = max;
6498d6e3deSJakub Kruzik   PetscFunctionReturn(0);
6598d6e3deSJakub Kruzik }
6698d6e3deSJakub Kruzik 
6798d6e3deSJakub Kruzik /*@
6898d6e3deSJakub Kruzik    PCDeflationSetMaxLvl - Set maximum level of deflation.
6998d6e3deSJakub Kruzik 
7098d6e3deSJakub Kruzik    Logically Collective on PC
7198d6e3deSJakub Kruzik 
7298d6e3deSJakub Kruzik    Input Parameters:
7398d6e3deSJakub Kruzik +  pc  - the preconditioner context
7498d6e3deSJakub Kruzik -  max - maximum deflation level
7598d6e3deSJakub Kruzik 
7698d6e3deSJakub Kruzik    Level: intermediate
7798d6e3deSJakub Kruzik 
7898d6e3deSJakub Kruzik .seealso: PCDEFLATION
7998d6e3deSJakub Kruzik @*/
8098d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
8198d6e3deSJakub Kruzik {
8298d6e3deSJakub Kruzik   PetscErrorCode ierr;
8398d6e3deSJakub Kruzik 
8498d6e3deSJakub Kruzik   PetscFunctionBegin;
8598d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8698d6e3deSJakub Kruzik   PetscValidLogicalCollectiveInt(pc,max,2);
8798d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
8898d6e3deSJakub Kruzik   PetscFunctionReturn(0);
8998d6e3deSJakub Kruzik }
9098d6e3deSJakub Kruzik 
91859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red)
92859a873cSJakub Kruzik {
93859a873cSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
94859a873cSJakub Kruzik 
95859a873cSJakub Kruzik   PetscFunctionBegin;
96859a873cSJakub Kruzik   def->reductionfact = red;
97859a873cSJakub Kruzik   PetscFunctionReturn(0);
98859a873cSJakub Kruzik }
99859a873cSJakub Kruzik 
100859a873cSJakub Kruzik /*@
101859a873cSJakub Kruzik    PCDeflationSetReductionFactor - Set reduction factor for PCTELESCOPE coarse problem solver.
102859a873cSJakub Kruzik 
103859a873cSJakub Kruzik    Logically Collective on PC
104859a873cSJakub Kruzik 
105859a873cSJakub Kruzik    Input Parameters:
106859a873cSJakub Kruzik +  pc  - the preconditioner context
107859a873cSJakub Kruzik -  red - reduction factor (or PETSC_DETERMINE)
108859a873cSJakub Kruzik 
109859a873cSJakub Kruzik    Level: intermediate
110859a873cSJakub Kruzik 
111859a873cSJakub Kruzik .seealso: PCDEFLATION
112859a873cSJakub Kruzik @*/
113859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red)
114859a873cSJakub Kruzik {
115859a873cSJakub Kruzik   PetscErrorCode ierr;
116859a873cSJakub Kruzik 
117859a873cSJakub Kruzik   PetscFunctionBegin;
118859a873cSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119859a873cSJakub Kruzik   PetscValidLogicalCollectiveInt(pc,red,2);
120859a873cSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr);
121859a873cSJakub Kruzik   PetscFunctionReturn(0);
122859a873cSJakub Kruzik }
123859a873cSJakub Kruzik 
1248a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact)
1258a71cb68SJakub Kruzik {
1268a71cb68SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
1278a71cb68SJakub Kruzik 
1288a71cb68SJakub Kruzik   PetscFunctionBegin;
1298a71cb68SJakub Kruzik   /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */
1308a71cb68SJakub Kruzik   def->correct = PETSC_TRUE;
1318a71cb68SJakub Kruzik   def->correctfact = fact;
1328a71cb68SJakub Kruzik   if (fact == PETSC_DEFAULT) {
1338a71cb68SJakub Kruzik     def->correctfact = 1.0;
1348a71cb68SJakub Kruzik   } else if (def->correct == 0.0) {
1358a71cb68SJakub Kruzik     def->correct = PETSC_FALSE;
1368a71cb68SJakub Kruzik   }
1378a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1388a71cb68SJakub Kruzik }
1398a71cb68SJakub Kruzik 
1408a71cb68SJakub Kruzik /*@
1418a71cb68SJakub Kruzik    PCDeflationSetCorrectionFactor - Set coarse problem correction factor.
1428a71cb68SJakub Kruzik     The Preconditioner becomes P*M^{-1} + factor*Q.
1438a71cb68SJakub Kruzik 
1448a71cb68SJakub Kruzik    Logically Collective on PC
1458a71cb68SJakub Kruzik 
1468a71cb68SJakub Kruzik    Input Parameters:
1478a71cb68SJakub Kruzik +  pc   - the preconditioner context
1488a71cb68SJakub Kruzik -  fact - correction factor (set 0.0 to disable, PETSC_DEFAULT = 1.0)
1498a71cb68SJakub Kruzik 
1508a71cb68SJakub Kruzik    Level: intermediate
1518a71cb68SJakub Kruzik 
1528a71cb68SJakub Kruzik .seealso: PCDEFLATION
1538a71cb68SJakub Kruzik @*/
1548a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact)
1558a71cb68SJakub Kruzik {
1568a71cb68SJakub Kruzik   PetscErrorCode ierr;
1578a71cb68SJakub Kruzik 
1588a71cb68SJakub Kruzik   PetscFunctionBegin;
1598a71cb68SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1608a71cb68SJakub Kruzik   PetscValidLogicalCollectiveScalar(pc,fact,2);
1618a71cb68SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr);
1628a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1638a71cb68SJakub Kruzik }
1648a71cb68SJakub Kruzik 
16539417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size)
16639417d7eSJakub Kruzik {
16739417d7eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
16839417d7eSJakub Kruzik 
16939417d7eSJakub Kruzik   PetscFunctionBegin;
17039417d7eSJakub Kruzik   if (type) def->spacetype = type;
17139417d7eSJakub Kruzik   if (size > 0) def->spacesize = size;
17239417d7eSJakub Kruzik   PetscFunctionReturn(0);
17339417d7eSJakub Kruzik }
17439417d7eSJakub Kruzik 
17539417d7eSJakub Kruzik /*@
17639417d7eSJakub Kruzik    PCDeflationSetSpaceToCompute - Set deflation space type and size to compute.
17739417d7eSJakub Kruzik 
17839417d7eSJakub Kruzik    Logically Collective on PC
17939417d7eSJakub Kruzik 
18039417d7eSJakub Kruzik    Input Parameters:
18139417d7eSJakub Kruzik +  pc   - the preconditioner context
18239417d7eSJakub Kruzik .  type - deflation space type to compute (or PETSC_IGNORE)
18339417d7eSJakub Kruzik -  size - size of the space to compute (or PETSC_DEFAULT)
18439417d7eSJakub Kruzik 
18539417d7eSJakub Kruzik    Notes:
18639417d7eSJakub Kruzik     For wavelet-based deflation, size represents number of levels.
18739417d7eSJakub Kruzik     The deflation space is computed in PCSetUP().
18839417d7eSJakub Kruzik 
18939417d7eSJakub Kruzik    Level: intermediate
19039417d7eSJakub Kruzik 
19139417d7eSJakub Kruzik .seealso: PCDEFLATION
19239417d7eSJakub Kruzik @*/
19339417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size)
19439417d7eSJakub Kruzik {
19539417d7eSJakub Kruzik   PetscErrorCode ierr;
19639417d7eSJakub Kruzik 
19739417d7eSJakub Kruzik   PetscFunctionBegin;
19839417d7eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
19939417d7eSJakub Kruzik   if (type) PetscValidLogicalCollectiveEnum(pc,type,2);
20039417d7eSJakub Kruzik   if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3);
20139417d7eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr);
20239417d7eSJakub Kruzik   PetscFunctionReturn(0);
20339417d7eSJakub Kruzik }
2048513960bSJakub Kruzik 
205e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
206e662bc50SJakub Kruzik {
207e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
208e662bc50SJakub Kruzik   PetscErrorCode ierr;
209e662bc50SJakub Kruzik 
210e662bc50SJakub Kruzik   PetscFunctionBegin;
211e662bc50SJakub Kruzik   if (transpose) {
212e662bc50SJakub Kruzik     def->Wt = W;
213e662bc50SJakub Kruzik     def->W = NULL;
214e662bc50SJakub Kruzik   } else {
215e662bc50SJakub Kruzik     def->W = W;
216e662bc50SJakub Kruzik   }
217e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
218e662bc50SJakub Kruzik   PetscFunctionReturn(0);
219e662bc50SJakub Kruzik }
220e662bc50SJakub Kruzik 
221e662bc50SJakub Kruzik /*@
222e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
223e662bc50SJakub Kruzik 
224e662bc50SJakub Kruzik    Logically Collective on PC
225e662bc50SJakub Kruzik 
226e662bc50SJakub Kruzik    Input Parameters:
227e662bc50SJakub Kruzik +  pc - the preconditioner context
228e662bc50SJakub Kruzik .  W  - deflation matrix
229e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
230e662bc50SJakub Kruzik 
231e662bc50SJakub Kruzik    Level: intermediate
232e662bc50SJakub Kruzik 
233e662bc50SJakub Kruzik .seealso: PCDEFLATION
234e662bc50SJakub Kruzik @*/
235e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
236e662bc50SJakub Kruzik {
237e662bc50SJakub Kruzik   PetscErrorCode ierr;
238e662bc50SJakub Kruzik 
239e662bc50SJakub Kruzik   PetscFunctionBegin;
240e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
241e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
242e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
243e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
244e662bc50SJakub Kruzik   PetscFunctionReturn(0);
245e662bc50SJakub Kruzik }
246e662bc50SJakub Kruzik 
2473cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat)
2483cf3a049SJakub Kruzik {
2493cf3a049SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
2503cf3a049SJakub Kruzik   PetscErrorCode   ierr;
2513cf3a049SJakub Kruzik 
2523cf3a049SJakub Kruzik   PetscFunctionBegin;
2533cf3a049SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
2543cf3a049SJakub Kruzik   def->WtA = mat;
2553cf3a049SJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2563cf3a049SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr);
2573cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2583cf3a049SJakub Kruzik }
2593cf3a049SJakub Kruzik 
2603cf3a049SJakub Kruzik /*@
2613cf3a049SJakub Kruzik    PCDeflationSetProjectionNullSpaceMat - Set projection null space matrix (W'*A).
2623cf3a049SJakub Kruzik 
2633cf3a049SJakub Kruzik    Not Collective
2643cf3a049SJakub Kruzik 
2653cf3a049SJakub Kruzik    Input Parameters:
2663cf3a049SJakub Kruzik +  pc  - preconditioner context
2673cf3a049SJakub Kruzik -  mat - projection null space matrix
2683cf3a049SJakub Kruzik 
2693cf3a049SJakub Kruzik    Level: developer
2703cf3a049SJakub Kruzik 
2713cf3a049SJakub Kruzik .seealso: PCDEFLATION
2723cf3a049SJakub Kruzik @*/
2733cf3a049SJakub Kruzik PetscErrorCode  PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat)
2743cf3a049SJakub Kruzik {
2753cf3a049SJakub Kruzik   PetscErrorCode ierr;
2763cf3a049SJakub Kruzik 
2773cf3a049SJakub Kruzik   PetscFunctionBegin;
2783cf3a049SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2793cf3a049SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
2803cf3a049SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
2813cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2823cf3a049SJakub Kruzik }
2833cf3a049SJakub Kruzik 
284e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
285e924b002SJakub Kruzik {
286e924b002SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
287e924b002SJakub Kruzik   PetscErrorCode   ierr;
288e924b002SJakub Kruzik 
289e924b002SJakub Kruzik   PetscFunctionBegin;
290e924b002SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
291e924b002SJakub Kruzik   def->WtAW = mat;
292e924b002SJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
293e924b002SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
294e924b002SJakub Kruzik   PetscFunctionReturn(0);
295e924b002SJakub Kruzik }
296e924b002SJakub Kruzik 
297e924b002SJakub Kruzik /*@
298e924b002SJakub Kruzik    PCDeflationSetCoarseMat - Set coarse problem Mat.
299e924b002SJakub Kruzik 
300e924b002SJakub Kruzik    Not Collective
301e924b002SJakub Kruzik 
302e924b002SJakub Kruzik    Input Parameters:
303e924b002SJakub Kruzik +  pc - preconditioner context
304e924b002SJakub Kruzik -  mat - coarse problem mat
305e924b002SJakub Kruzik 
306e924b002SJakub Kruzik    Level: developer
307e924b002SJakub Kruzik 
308e924b002SJakub Kruzik .seealso: PCDEFLATION
309e924b002SJakub Kruzik @*/
310e924b002SJakub Kruzik PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
311e924b002SJakub Kruzik {
312e924b002SJakub Kruzik   PetscErrorCode ierr;
313e924b002SJakub Kruzik 
314e924b002SJakub Kruzik   PetscFunctionBegin;
315e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
316e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
317e924b002SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
318e924b002SJakub Kruzik   PetscFunctionReturn(0);
319e924b002SJakub Kruzik }
320e924b002SJakub Kruzik 
32198d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
3225c0c31c5SJakub Kruzik {
3235c0c31c5SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
3245c0c31c5SJakub Kruzik 
3255c0c31c5SJakub Kruzik   PetscFunctionBegin;
32698d6e3deSJakub Kruzik   *ksp = def->WtAWinv;
3275c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
3285c0c31c5SJakub Kruzik }
3295c0c31c5SJakub Kruzik 
3305c0c31c5SJakub Kruzik /*@
33198d6e3deSJakub Kruzik    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
3325c0c31c5SJakub Kruzik 
33398d6e3deSJakub Kruzik    Not Collective
3345c0c31c5SJakub Kruzik 
3355c0c31c5SJakub Kruzik    Input Parameters:
33698d6e3deSJakub Kruzik .  pc - preconditioner context
3375c0c31c5SJakub Kruzik 
33898d6e3deSJakub Kruzik    Output Parameter:
33998d6e3deSJakub Kruzik .  ksp - coarse problem KSP context
3405c0c31c5SJakub Kruzik 
34198d6e3deSJakub Kruzik    Level: developer
34298d6e3deSJakub Kruzik 
34398d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION
3445c0c31c5SJakub Kruzik @*/
34598d6e3deSJakub Kruzik PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
3465c0c31c5SJakub Kruzik {
3475c0c31c5SJakub Kruzik   PetscErrorCode ierr;
3485c0c31c5SJakub Kruzik 
3495c0c31c5SJakub Kruzik   PetscFunctionBegin;
35022b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
35198d6e3deSJakub Kruzik   PetscValidPointer(ksp,2);
35298d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
35398d6e3deSJakub Kruzik   PetscFunctionReturn(0);
35498d6e3deSJakub Kruzik }
35598d6e3deSJakub Kruzik 
35698d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp)
35798d6e3deSJakub Kruzik {
35898d6e3deSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
35998d6e3deSJakub Kruzik   PetscErrorCode   ierr;
36098d6e3deSJakub Kruzik 
36198d6e3deSJakub Kruzik   PetscFunctionBegin;
36298d6e3deSJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
36398d6e3deSJakub Kruzik   def->WtAWinv = ksp;
36498d6e3deSJakub Kruzik   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
36598d6e3deSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr);
36698d6e3deSJakub Kruzik   PetscFunctionReturn(0);
36798d6e3deSJakub Kruzik }
36898d6e3deSJakub Kruzik 
36998d6e3deSJakub Kruzik /*@
37098d6e3deSJakub Kruzik    PCDeflationSetCoarseKSP - Set coarse problem KSP.
37198d6e3deSJakub Kruzik 
37298d6e3deSJakub Kruzik    Collective on PC
37398d6e3deSJakub Kruzik 
37498d6e3deSJakub Kruzik    Input Parameters:
37598d6e3deSJakub Kruzik +  pc - preconditioner context
37698d6e3deSJakub Kruzik -  ksp - coarse problem KSP context
37798d6e3deSJakub Kruzik 
37898d6e3deSJakub Kruzik    Level: developer
37998d6e3deSJakub Kruzik 
38098d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION
38198d6e3deSJakub Kruzik @*/
38298d6e3deSJakub Kruzik PetscErrorCode  PCDeflationSetCoarseKSP(PC pc,KSP ksp)
38398d6e3deSJakub Kruzik {
38498d6e3deSJakub Kruzik   PetscErrorCode ierr;
38598d6e3deSJakub Kruzik 
38698d6e3deSJakub Kruzik   PetscFunctionBegin;
38798d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
38898d6e3deSJakub Kruzik   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
38998d6e3deSJakub Kruzik   PetscCheckSameComm(pc,1,ksp,2);
39098d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr);
3915c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
3925c0c31c5SJakub Kruzik }
393e662bc50SJakub Kruzik 
394268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
395268c6673SJakub Kruzik {
396268c6673SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
397268c6673SJakub Kruzik 
398268c6673SJakub Kruzik   PetscFunctionBegin;
399268c6673SJakub Kruzik   *apc = def->pc;
400268c6673SJakub Kruzik   PetscFunctionReturn(0);
401268c6673SJakub Kruzik }
402268c6673SJakub Kruzik 
403268c6673SJakub Kruzik /*@
404268c6673SJakub Kruzik    PCDeflationGetPC - Returns a pointer to additional preconditioner.
405268c6673SJakub Kruzik 
406268c6673SJakub Kruzik    Not Collective
407268c6673SJakub Kruzik 
408268c6673SJakub Kruzik    Input Parameters:
409268c6673SJakub Kruzik .  pc  - the preconditioner context
410268c6673SJakub Kruzik 
411268c6673SJakub Kruzik    Output Parameter:
412268c6673SJakub Kruzik .  apc - additional preconditioner
413268c6673SJakub Kruzik 
414268c6673SJakub Kruzik    Level: advanced
415268c6673SJakub Kruzik 
416268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION
417268c6673SJakub Kruzik @*/
418268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
419268c6673SJakub Kruzik {
420268c6673SJakub Kruzik   PetscErrorCode ierr;
421268c6673SJakub Kruzik 
422268c6673SJakub Kruzik   PetscFunctionBegin;
423268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
424268c6673SJakub Kruzik   PetscValidPointer(pc,2);
425268c6673SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
426268c6673SJakub Kruzik   PetscFunctionReturn(0);
427268c6673SJakub Kruzik }
428268c6673SJakub Kruzik 
42922b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
43022b0793eSJakub Kruzik {
43122b0793eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
43222b0793eSJakub Kruzik   PetscErrorCode ierr;
43322b0793eSJakub Kruzik 
43422b0793eSJakub Kruzik   PetscFunctionBegin;
43522b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
43622b0793eSJakub Kruzik   def->pc = apc;
43722b0793eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
43822b0793eSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
43922b0793eSJakub Kruzik   PetscFunctionReturn(0);
44022b0793eSJakub Kruzik }
44122b0793eSJakub Kruzik 
44222b0793eSJakub Kruzik /*@
44322b0793eSJakub Kruzik    PCDeflationSetPC - Set additional preconditioner.
44422b0793eSJakub Kruzik 
445268c6673SJakub Kruzik    Collective on PC
44622b0793eSJakub Kruzik 
44722b0793eSJakub Kruzik    Input Parameters:
44822b0793eSJakub Kruzik +  pc  - the preconditioner context
44922b0793eSJakub Kruzik -  apc - additional preconditioner
45022b0793eSJakub Kruzik 
451268c6673SJakub Kruzik    Level: developer
45222b0793eSJakub Kruzik 
453268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION
45422b0793eSJakub Kruzik @*/
45522b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
45622b0793eSJakub Kruzik {
45722b0793eSJakub Kruzik   PetscErrorCode ierr;
45822b0793eSJakub Kruzik 
45922b0793eSJakub Kruzik   PetscFunctionBegin;
46022b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
46122b0793eSJakub Kruzik   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
46222b0793eSJakub Kruzik   PetscCheckSameComm(pc,1,apc,2);
46322b0793eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
46422b0793eSJakub Kruzik   PetscFunctionReturn(0);
46522b0793eSJakub Kruzik }
46622b0793eSJakub Kruzik 
46737eeb815SJakub Kruzik /*
468b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
469b27c8092SJakub Kruzik */
470b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
471b27c8092SJakub Kruzik {
472b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
473b27c8092SJakub Kruzik   Mat              A;
474b27c8092SJakub Kruzik   Vec              r,w1,w2;
475b27c8092SJakub Kruzik   PetscBool        nonzero;
476b27c8092SJakub Kruzik   PetscErrorCode   ierr;
47737eeb815SJakub Kruzik 
478b27c8092SJakub Kruzik   PetscFunctionBegin;
479b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
480b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
481b27c8092SJakub Kruzik   r  = def->work;
482b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
48337eeb815SJakub Kruzik 
484b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
485b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
486b27c8092SJakub Kruzik   if (nonzero) {
487b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);                          /*    r  <- b - Ax              */
488b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
489b27c8092SJakub Kruzik   } else {
490b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                            /*    r  <- b (x is 0)          */
491b27c8092SJakub Kruzik   }
492b27c8092SJakub Kruzik 
493a3177931SJakub Kruzik   if (def->Wt) {
494a3177931SJakub Kruzik     ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr);                   /*    w1 <- W'*r                */
495a3177931SJakub Kruzik   } else {
496a3177931SJakub Kruzik     ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
497a3177931SJakub Kruzik   }
498b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);              /*    w2 <- (W'*A*W)^{-1}*w1    */
499b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);                      /*    r  <- W*w2                */
500b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
501b27c8092SJakub Kruzik   PetscFunctionReturn(0);
502b27c8092SJakub Kruzik }
50337eeb815SJakub Kruzik 
504f8bf229cSJakub Kruzik /*
505f8bf229cSJakub Kruzik   if (def->correct) {
506ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
507f8bf229cSJakub Kruzik   } else {
508ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
509f8bf229cSJakub Kruzik   }
51037eeb815SJakub Kruzik */
511f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
512f8bf229cSJakub Kruzik {
513f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
514f8bf229cSJakub Kruzik   Mat              A;
515f8bf229cSJakub Kruzik   Vec              u,w1,w2;
516f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
517f8bf229cSJakub Kruzik 
518f8bf229cSJakub Kruzik   PetscFunctionBegin;
519f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
520f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
521f8bf229cSJakub Kruzik   u  = def->work;
522f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
523f8bf229cSJakub Kruzik 
524ae029463SJakub Kruzik   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                          /*    z <- M^{-1}*r             */
52522b0793eSJakub Kruzik   if (!def->init) {
526ae029463SJakub Kruzik     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                      /*    w1 <- W'*A*z              */
527f8bf229cSJakub Kruzik     if (def->correct) {
528ae029463SJakub Kruzik       if (def->Wt) {
529ae029463SJakub Kruzik         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);                   /*    w2 <- W'*r                */
530ae029463SJakub Kruzik       } else {
531a3177931SJakub Kruzik         ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr);  /*    w2 <- W'*r                */
532f8bf229cSJakub Kruzik       }
533ae029463SJakub Kruzik       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);      /*    w1 <- w1 - l*w2           */
534f8bf229cSJakub Kruzik     }
535f8bf229cSJakub Kruzik     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);                /*    w2 <- (W'*A*W)^{-1}*w1    */
536f8bf229cSJakub Kruzik     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                        /*    u  <- W*w2                */
53722b0793eSJakub Kruzik     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                           /*    z  <- z - u               */
53822b0793eSJakub Kruzik   }
539f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
540f8bf229cSJakub Kruzik }
541f8bf229cSJakub Kruzik 
54237eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
54337eeb815SJakub Kruzik {
544409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
545409a357bSJakub Kruzik   KSP              innerksp;
546409a357bSJakub Kruzik   PC               pcinner;
547409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
548409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
549409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
550409a357bSJakub Kruzik   MatCompositeType ctype;
551409a357bSJakub Kruzik   MPI_Comm         comm;
552409a357bSJakub Kruzik   const char       *prefix;
55337eeb815SJakub Kruzik   PetscErrorCode   ierr;
55437eeb815SJakub Kruzik 
55537eeb815SJakub Kruzik   PetscFunctionBegin;
556409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
557409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
558f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
559*0a78af22SJakub Kruzik   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
560f8bf229cSJakub Kruzik 
561f8bf229cSJakub Kruzik   /* compute a deflation space */
562409a357bSJakub Kruzik   if (def->W || def->Wt) {
563409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
564409a357bSJakub Kruzik   } else {
565e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
56637eeb815SJakub Kruzik   }
56737eeb815SJakub Kruzik 
568409a357bSJakub Kruzik   /* nested deflation */
569409a357bSJakub Kruzik   if (def->W) {
570409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
571409a357bSJakub Kruzik     if (match) {
572409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
573409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
57437eeb815SJakub Kruzik     }
575409a357bSJakub Kruzik   } else {
576a3177931SJakub Kruzik     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
577409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
578409a357bSJakub Kruzik     if (match) {
579409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
580409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
581409a357bSJakub Kruzik     }
582409a357bSJakub Kruzik     transp = PETSC_TRUE;
583409a357bSJakub Kruzik   }
584409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
585409a357bSJakub Kruzik     if (!transp) {
58628da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
58728da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
588409a357bSJakub Kruzik         for (i=0; i<size; i++) {
589409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
590409a357bSJakub Kruzik         }
591409a357bSJakub Kruzik         size -= 1;
592409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
593409a357bSJakub Kruzik         def->W = mats[size];
594409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
595409a357bSJakub Kruzik         if (size > 1) {
596409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
597409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
598409a357bSJakub Kruzik         } else {
599409a357bSJakub Kruzik           nextDef = mats[0];
600409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
601409a357bSJakub Kruzik         }
60228da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
603409a357bSJakub Kruzik       } else {
604409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
605409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
606409a357bSJakub Kruzik       }
60728da0170SJakub Kruzik     } else {
60828da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
60928da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
61028da0170SJakub Kruzik         for (i=0; i<size; i++) {
61128da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
61228da0170SJakub Kruzik         }
61328da0170SJakub Kruzik         size -= 1;
61428da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
61528da0170SJakub Kruzik         def->Wt = mats[0];
61628da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
61728da0170SJakub Kruzik         if (size > 1) {
61828da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
61928da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
62028da0170SJakub Kruzik         } else {
62128da0170SJakub Kruzik           nextDef = mats[1];
62228da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
623409a357bSJakub Kruzik         }
624409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
62528da0170SJakub Kruzik       } else {
62628da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
62728da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
62828da0170SJakub Kruzik       }
62928da0170SJakub Kruzik     }
63028da0170SJakub Kruzik   }
63128da0170SJakub Kruzik 
63228da0170SJakub Kruzik   if (transp) {
63328da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
634a3177931SJakub Kruzik     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
635409a357bSJakub Kruzik   }
636409a357bSJakub Kruzik 
637ae029463SJakub Kruzik   /* assemble WtA */
638ae029463SJakub Kruzik   if (!def->WtA) {
639ae029463SJakub Kruzik     if (def->Wt) {
640ae029463SJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
641ae029463SJakub Kruzik     } else {
642a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
643a3177931SJakub Kruzik       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
644a3177931SJakub Kruzik       ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
645a3177931SJakub Kruzik #else
646ae029463SJakub Kruzik       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
647a3177931SJakub Kruzik #endif
648ae029463SJakub Kruzik     }
649ae029463SJakub Kruzik   }
650409a357bSJakub Kruzik   /* setup coarse problem */
651409a357bSJakub Kruzik   if (!def->WtAWinv) {
65228da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
653409a357bSJakub Kruzik     if (!def->WtAW) {
654ae029463SJakub Kruzik       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
655409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
656409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
657409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
658409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
659ae029463SJakub Kruzik       /* Check columns of W are not in kernel of A */
660409a357bSJakub Kruzik       PetscReal *norms;
661409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
662409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
663409a357bSJakub Kruzik       for (i=0; i<m; i++) {
664409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
665409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
666409a357bSJakub Kruzik         }
667409a357bSJakub Kruzik       }
668409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
669409a357bSJakub Kruzik #endif
670409a357bSJakub Kruzik     } else {
671409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
672409a357bSJakub Kruzik     }
673409a357bSJakub Kruzik     /* TODO use MATINV */
674409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
675409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
676409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
677557a2f7dSJakub Kruzik     /* Setup KSP and PC */
678557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
679557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
680557a2f7dSJakub Kruzik       /* set default KSPtype */
681557a2f7dSJakub Kruzik       if (!def->ksptype) {
682557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
683557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
684557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
685557a2f7dSJakub Kruzik         }
686557a2f7dSJakub Kruzik       }
687ae029463SJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
688557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
689557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
690557a2f7dSJakub Kruzik       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
691ae029463SJakub Kruzik       /* inherit options */
692557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
693557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
694557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
695557a2f7dSJakub Kruzik     } else { /* the last level */
696557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
697409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
698ac66f006SJakub Kruzik       /* ugly hack to not have overwritten PCTELESCOPE */
699ac66f006SJakub Kruzik       if (prefix) {
700ac66f006SJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
701ac66f006SJakub Kruzik       }
70222b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
703ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
704409a357bSJakub Kruzik       /* Reduction factor choice */
705409a357bSJakub Kruzik       red = def->reductionfact;
706409a357bSJakub Kruzik       if (red < 0) {
707409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
708409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
709409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
710409a357bSJakub Kruzik         if (match) red = commsize;
711409a357bSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
712409a357bSJakub Kruzik       }
713409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
714ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
715409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
716ac66f006SJakub Kruzik       if (innerksp) {
717409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
718409a357bSJakub Kruzik         /* TODO Cholesky if flgspd? */
719ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
720409a357bSJakub Kruzik         //TODO remove explicit matSolverPackage
721409a357bSJakub Kruzik         if (commsize == red) {
722ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
723409a357bSJakub Kruzik         } else {
724ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
725409a357bSJakub Kruzik         }
726409a357bSJakub Kruzik       }
727557a2f7dSJakub Kruzik     }
728557a2f7dSJakub Kruzik 
729557a2f7dSJakub Kruzik     if (innerksp) {
730409a357bSJakub Kruzik       /* TODO use def_[lvl]_ if lvl > 0? */
731409a357bSJakub Kruzik       if (prefix) {
732409a357bSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
733409a357bSJakub Kruzik       }
73422b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
735557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
736557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
737ac66f006SJakub Kruzik     }
738409a357bSJakub Kruzik   }
739557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
740557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
741409a357bSJakub Kruzik 
74222b0793eSJakub Kruzik   if (!def->pc) {
74322b0793eSJakub Kruzik     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
74422b0793eSJakub Kruzik     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
74522b0793eSJakub Kruzik     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
74622b0793eSJakub Kruzik     if (prefix) {
74722b0793eSJakub Kruzik       ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
74822b0793eSJakub Kruzik     }
74922b0793eSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr);
75022b0793eSJakub Kruzik     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
75122b0793eSJakub Kruzik     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
75222b0793eSJakub Kruzik   }
75322b0793eSJakub Kruzik 
754f8bf229cSJakub Kruzik   /* create work vecs */
755f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
756f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
75737eeb815SJakub Kruzik   PetscFunctionReturn(0);
75837eeb815SJakub Kruzik }
759b30d4775SJakub Kruzik 
76037eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
76137eeb815SJakub Kruzik {
762b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
76337eeb815SJakub Kruzik   PetscErrorCode    ierr;
76437eeb815SJakub Kruzik 
76537eeb815SJakub Kruzik   PetscFunctionBegin;
766b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
767b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
768b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
769b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
770ae029463SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
771b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
772b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
77322b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
77437eeb815SJakub Kruzik   PetscFunctionReturn(0);
77537eeb815SJakub Kruzik }
77637eeb815SJakub Kruzik 
77737eeb815SJakub Kruzik /*
77837eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
77937eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
78037eeb815SJakub Kruzik 
78137eeb815SJakub Kruzik    Input Parameter:
78237eeb815SJakub Kruzik .  pc - the preconditioner context
78337eeb815SJakub Kruzik 
78437eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
78537eeb815SJakub Kruzik */
78637eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
78737eeb815SJakub Kruzik {
78837eeb815SJakub Kruzik   PetscErrorCode ierr;
78937eeb815SJakub Kruzik 
79037eeb815SJakub Kruzik   PetscFunctionBegin;
79137eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
79237eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
79337eeb815SJakub Kruzik   PetscFunctionReturn(0);
79437eeb815SJakub Kruzik }
79537eeb815SJakub Kruzik 
7968513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
7978513960bSJakub Kruzik {
7988513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
7991ac6250aSJakub Kruzik   PetscInt          its;
8008513960bSJakub Kruzik   PetscBool         iascii;
8018513960bSJakub Kruzik   PetscErrorCode    ierr;
8028513960bSJakub Kruzik 
8038513960bSJakub Kruzik   PetscFunctionBegin;
8048513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8058513960bSJakub Kruzik   if (iascii) {
8068513960bSJakub Kruzik     if (def->correct) {
8071ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n",
8081ac6250aSJakub Kruzik                                     (double)PetscRealPart(def->correctfact),
8091ac6250aSJakub Kruzik                                     (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr);
8108513960bSJakub Kruzik     }
8118513960bSJakub Kruzik     if (!def->nestedlvl) {
8121ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
8138513960bSJakub Kruzik     }
8148513960bSJakub Kruzik 
8151ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr);
81622b0793eSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
81722b0793eSJakub Kruzik     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
81822b0793eSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
81922b0793eSJakub Kruzik 
8201ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr);
8218513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
8221ac6250aSJakub Kruzik     ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr);
8231ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr);
8248513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
8258513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
8268513960bSJakub Kruzik   }
8278513960bSJakub Kruzik   PetscFunctionReturn(0);
8288513960bSJakub Kruzik }
8298513960bSJakub Kruzik 
83037eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
83137eeb815SJakub Kruzik {
832880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
83337eeb815SJakub Kruzik   PetscErrorCode    ierr;
83437eeb815SJakub Kruzik 
83537eeb815SJakub Kruzik   PetscFunctionBegin;
83637eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
837a122ebaeSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr);
838859a873cSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
839859a873cSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
8408a71cb68SJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr);
8418a71cb68SJakub 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);
842880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
843880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
844880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
845880d05ceSJakub Kruzik //TODO add set function and fix manpages
84637eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
84737eeb815SJakub Kruzik   PetscFunctionReturn(0);
84837eeb815SJakub Kruzik }
84937eeb815SJakub Kruzik 
85037eeb815SJakub Kruzik /*MC
851e361f8b5SJakub Kruzik      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
852e361f8b5SJakub Kruzik      or to a predefined value
85337eeb815SJakub Kruzik 
85437eeb815SJakub Kruzik    Options Database Key:
855e361f8b5SJakub Kruzik +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
85637eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
85737eeb815SJakub Kruzik 
85837eeb815SJakub Kruzik    Level: beginner
85937eeb815SJakub Kruzik 
86037eeb815SJakub Kruzik   Notes:
861e361f8b5SJakub Kruzik     todo
86237eeb815SJakub Kruzik 
86337eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
864e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
86537eeb815SJakub Kruzik M*/
86637eeb815SJakub Kruzik 
86737eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
86837eeb815SJakub Kruzik {
86937eeb815SJakub Kruzik   PC_Deflation   *def;
87037eeb815SJakub Kruzik   PetscErrorCode ierr;
87137eeb815SJakub Kruzik 
87237eeb815SJakub Kruzik   PetscFunctionBegin;
87337eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
87437eeb815SJakub Kruzik   pc->data = (void*)def;
87537eeb815SJakub Kruzik 
876e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
877e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
878fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
879e662bc50SJakub Kruzik   def->reductionfact = -1;
880e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
881e662bc50SJakub Kruzik   def->spacesize     = 1;
882e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
883e662bc50SJakub Kruzik   def->nestedlvl     = 0;
884e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
88537eeb815SJakub Kruzik 
88637eeb815SJakub Kruzik   pc->ops->apply          = PCApply_Deflation;
887b27c8092SJakub Kruzik   pc->ops->presolve       = PCPreSolve_Deflation;
88837eeb815SJakub Kruzik   pc->ops->setup          = PCSetUp_Deflation;
88937eeb815SJakub Kruzik   pc->ops->reset          = PCReset_Deflation;
89037eeb815SJakub Kruzik   pc->ops->destroy        = PCDestroy_Deflation;
89137eeb815SJakub Kruzik   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
8928513960bSJakub Kruzik   pc->ops->view           = PCView_Deflation;
89337eeb815SJakub Kruzik 
894a122ebaeSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr);
89598d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
896859a873cSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr);
8978a71cb68SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr);
89839417d7eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr);
899e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
9003cf3a049SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr);
901e924b002SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
9024a99276eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
903f3eaa4f0SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
90498d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
90598d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
90637eeb815SJakub Kruzik   PetscFunctionReturn(0);
90737eeb815SJakub Kruzik }
90837eeb815SJakub Kruzik 
909