xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
15170378fSJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/ /* includes for fortran wrappers */
237eeb815SJakub Kruzik 
39371c9d4SSatish Balay const char *const PCDeflationSpaceTypes[] = {"haar", "db2", "db4", "db8", "db16", "biorth22", "meyer", "aggregation", "user", "PCDeflationSpaceType", "PC_DEFLATION_SPACE_", NULL};
48513960bSJakub Kruzik 
59371c9d4SSatish Balay static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc, PetscBool flg) {
6a122ebaeSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
7a122ebaeSJakub Kruzik 
8a122ebaeSJakub Kruzik   PetscFunctionBegin;
9a122ebaeSJakub Kruzik   def->init = flg;
10a122ebaeSJakub Kruzik   PetscFunctionReturn(0);
11a122ebaeSJakub Kruzik }
12a122ebaeSJakub Kruzik 
13a122ebaeSJakub Kruzik /*@
14a122ebaeSJakub Kruzik    PCDeflationSetInitOnly - Do only initialization step.
15e26ad46dSJakub Kruzik     Sets initial guess to the solution on the deflation space but does not apply
16e26ad46dSJakub Kruzik     the deflation preconditioner. The additional preconditioner is still applied.
17a122ebaeSJakub Kruzik 
18e26ad46dSJakub Kruzik    Logically Collective
19a122ebaeSJakub Kruzik 
20a122ebaeSJakub Kruzik    Input Parameters:
21a122ebaeSJakub Kruzik +  pc  - the preconditioner context
22*f1580f4eSBarry Smith -  flg - default `PETSC_FALSE`
23a122ebaeSJakub Kruzik 
24*f1580f4eSBarry Smith    Options Database Key:
25e26ad46dSJakub Kruzik .    -pc_deflation_init_only <false> - if true computes only the special guess
26e26ad46dSJakub Kruzik 
27a122ebaeSJakub Kruzik    Level: intermediate
28a122ebaeSJakub Kruzik 
29db781477SPatrick Sanan .seealso: `PCDEFLATION`
30a122ebaeSJakub Kruzik @*/
319371c9d4SSatish Balay PetscErrorCode PCDeflationSetInitOnly(PC pc, PetscBool flg) {
32a122ebaeSJakub Kruzik   PetscFunctionBegin;
33a122ebaeSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
34a122ebaeSJakub Kruzik   PetscValidLogicalCollectiveBool(pc, flg, 2);
35cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetInitOnly_C", (PC, PetscBool), (pc, flg));
36a122ebaeSJakub Kruzik   PetscFunctionReturn(0);
37a122ebaeSJakub Kruzik }
38a122ebaeSJakub Kruzik 
399371c9d4SSatish Balay static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc, PetscInt current, PetscInt max) {
4098d6e3deSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
4198d6e3deSJakub Kruzik 
4298d6e3deSJakub Kruzik   PetscFunctionBegin;
436c93e71cSJakub Kruzik   if (current) def->lvl = current;
446c93e71cSJakub Kruzik   def->maxlvl = max;
4598d6e3deSJakub Kruzik   PetscFunctionReturn(0);
4698d6e3deSJakub Kruzik }
4798d6e3deSJakub Kruzik 
4898d6e3deSJakub Kruzik /*@
4993b79a5aSJakub Kruzik    PCDeflationSetLevels - Set the maximum level of deflation nesting.
5098d6e3deSJakub Kruzik 
51e26ad46dSJakub Kruzik    Logically Collective
5298d6e3deSJakub Kruzik 
5398d6e3deSJakub Kruzik    Input Parameters:
5498d6e3deSJakub Kruzik +  pc  - the preconditioner context
5598d6e3deSJakub Kruzik -  max - maximum deflation level
5698d6e3deSJakub Kruzik 
57*f1580f4eSBarry Smith    Options Database Key:
58e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation
59e26ad46dSJakub Kruzik 
6098d6e3deSJakub Kruzik    Level: intermediate
6198d6e3deSJakub Kruzik 
62db781477SPatrick Sanan .seealso: `PCDeflationSetSpaceToCompute()`, `PCDeflationSetSpace()`, `PCDEFLATION`
6398d6e3deSJakub Kruzik @*/
649371c9d4SSatish Balay PetscErrorCode PCDeflationSetLevels(PC pc, PetscInt max) {
6598d6e3deSJakub Kruzik   PetscFunctionBegin;
6698d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6798d6e3deSJakub Kruzik   PetscValidLogicalCollectiveInt(pc, max, 2);
68cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetLevels_C", (PC, PetscInt, PetscInt), (pc, 0, max));
6998d6e3deSJakub Kruzik   PetscFunctionReturn(0);
7098d6e3deSJakub Kruzik }
7198d6e3deSJakub Kruzik 
729371c9d4SSatish Balay static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc, PetscInt red) {
73859a873cSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
74859a873cSJakub Kruzik 
75859a873cSJakub Kruzik   PetscFunctionBegin;
76859a873cSJakub Kruzik   def->reductionfact = red;
77859a873cSJakub Kruzik   PetscFunctionReturn(0);
78859a873cSJakub Kruzik }
79859a873cSJakub Kruzik 
80859a873cSJakub Kruzik /*@
81*f1580f4eSBarry Smith    PCDeflationSetReductionFactor - Set reduction factor for the `PCDEFLATION`
82859a873cSJakub Kruzik 
83e26ad46dSJakub Kruzik    Logically Collective
84859a873cSJakub Kruzik 
85859a873cSJakub Kruzik    Input Parameters:
86859a873cSJakub Kruzik +  pc  - the preconditioner context
87*f1580f4eSBarry Smith -  red - reduction factor (or `PETSC_DETERMINE`)
88859a873cSJakub Kruzik 
89*f1580f4eSBarry Smith    Options Database Key:
90*f1580f4eSBarry Smith .    -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for `PCDEFLATION`
91e26ad46dSJakub Kruzik 
92*f1580f4eSBarry Smith    Note:
93e26ad46dSJakub Kruzik    Default is computed based on the size of the coarse problem.
94e26ad46dSJakub Kruzik 
95859a873cSJakub Kruzik    Level: intermediate
96859a873cSJakub Kruzik 
97*f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCDEFLATION`, `PCDeflationSetLevels()`
98859a873cSJakub Kruzik @*/
999371c9d4SSatish Balay PetscErrorCode PCDeflationSetReductionFactor(PC pc, PetscInt red) {
100859a873cSJakub Kruzik   PetscFunctionBegin;
101859a873cSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
102859a873cSJakub Kruzik   PetscValidLogicalCollectiveInt(pc, red, 2);
103cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetReductionFactor_C", (PC, PetscInt), (pc, red));
104859a873cSJakub Kruzik   PetscFunctionReturn(0);
105859a873cSJakub Kruzik }
106859a873cSJakub Kruzik 
1079371c9d4SSatish Balay static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc, PetscScalar fact) {
1088a71cb68SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
1098a71cb68SJakub Kruzik 
1108a71cb68SJakub Kruzik   PetscFunctionBegin;
1118a71cb68SJakub Kruzik   /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */
1128a71cb68SJakub Kruzik   def->correct     = PETSC_TRUE;
1138a71cb68SJakub Kruzik   def->correctfact = fact;
114ad540459SPierre Jolivet   if (def->correct == 0.0) def->correct = PETSC_FALSE;
1158a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1168a71cb68SJakub Kruzik }
1178a71cb68SJakub Kruzik 
1188a71cb68SJakub Kruzik /*@
1198a71cb68SJakub Kruzik    PCDeflationSetCorrectionFactor - Set coarse problem correction factor.
120*f1580f4eSBarry Smith     The preconditioner becomes P*M^{-1} + fact*Q.
1218a71cb68SJakub Kruzik 
122e26ad46dSJakub Kruzik    Logically Collective
1238a71cb68SJakub Kruzik 
1248a71cb68SJakub Kruzik    Input Parameters:
1258a71cb68SJakub Kruzik +  pc   - the preconditioner context
126e26ad46dSJakub Kruzik -  fact - correction factor
127e26ad46dSJakub Kruzik 
128e26ad46dSJakub Kruzik    Options Database Keys:
129e26ad46dSJakub Kruzik +    -pc_deflation_correction        <false> - if true apply coarse problem correction
130e26ad46dSJakub Kruzik -    -pc_deflation_correction_factor <1.0>   - sets coarse problem correction factor
131e26ad46dSJakub Kruzik 
132*f1580f4eSBarry Smith    Note:
133e26ad46dSJakub Kruzik     Any non-zero fact enables the coarse problem correction.
1348a71cb68SJakub Kruzik 
1358a71cb68SJakub Kruzik    Level: intermediate
1368a71cb68SJakub Kruzik 
137*f1580f4eSBarry Smith .seealso: `PCDEFLATION`, `PCDeflationSetLevels()`, `PCDeflationSetReductionFactor()`
1388a71cb68SJakub Kruzik @*/
1399371c9d4SSatish Balay PetscErrorCode PCDeflationSetCorrectionFactor(PC pc, PetscScalar fact) {
1408a71cb68SJakub Kruzik   PetscFunctionBegin;
1418a71cb68SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1428a71cb68SJakub Kruzik   PetscValidLogicalCollectiveScalar(pc, fact, 2);
143cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetCorrectionFactor_C", (PC, PetscScalar), (pc, fact));
1448a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1458a71cb68SJakub Kruzik }
1468a71cb68SJakub Kruzik 
1479371c9d4SSatish Balay static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc, PCDeflationSpaceType type, PetscInt size) {
14839417d7eSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
14939417d7eSJakub Kruzik 
15039417d7eSJakub Kruzik   PetscFunctionBegin;
15139417d7eSJakub Kruzik   if (type) def->spacetype = type;
15239417d7eSJakub Kruzik   if (size > 0) def->spacesize = size;
15339417d7eSJakub Kruzik   PetscFunctionReturn(0);
15439417d7eSJakub Kruzik }
15539417d7eSJakub Kruzik 
15639417d7eSJakub Kruzik /*@
15739417d7eSJakub Kruzik    PCDeflationSetSpaceToCompute - Set deflation space type and size to compute.
15839417d7eSJakub Kruzik 
159e26ad46dSJakub Kruzik    Logically Collective
16039417d7eSJakub Kruzik 
16139417d7eSJakub Kruzik    Input Parameters:
16239417d7eSJakub Kruzik +  pc   - the preconditioner context
163*f1580f4eSBarry Smith .  type - deflation space type to compute (or `PETSC_IGNORE`)
164*f1580f4eSBarry Smith -  size - size of the space to compute (or `PETSC_DEFAULT`)
16539417d7eSJakub Kruzik 
166e26ad46dSJakub Kruzik    Options Database Keys:
167*f1580f4eSBarry Smith +    -pc_deflation_compute_space      <haar> - compute `PCDeflationSpaceType` deflation space
168e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>    - size of the deflation space
169e26ad46dSJakub Kruzik 
17039417d7eSJakub Kruzik    Notes:
17139417d7eSJakub Kruzik     For wavelet-based deflation, size represents number of levels.
172e26ad46dSJakub Kruzik 
173*f1580f4eSBarry Smith     The deflation space is computed in `PCSetUp()`.
17439417d7eSJakub Kruzik 
17539417d7eSJakub Kruzik    Level: intermediate
17639417d7eSJakub Kruzik 
177db781477SPatrick Sanan .seealso: `PCDeflationSetLevels()`, `PCDEFLATION`
17839417d7eSJakub Kruzik @*/
1799371c9d4SSatish Balay PetscErrorCode PCDeflationSetSpaceToCompute(PC pc, PCDeflationSpaceType type, PetscInt size) {
18039417d7eSJakub Kruzik   PetscFunctionBegin;
18139417d7eSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
18239417d7eSJakub Kruzik   if (type) PetscValidLogicalCollectiveEnum(pc, type, 2);
18339417d7eSJakub Kruzik   if (size > 0) PetscValidLogicalCollectiveInt(pc, size, 3);
184cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetSpaceToCompute_C", (PC, PCDeflationSpaceType, PetscInt), (pc, type, size));
18539417d7eSJakub Kruzik   PetscFunctionReturn(0);
18639417d7eSJakub Kruzik }
1878513960bSJakub Kruzik 
1889371c9d4SSatish Balay static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc, Mat W, PetscBool transpose) {
189e662bc50SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
190e662bc50SJakub Kruzik 
191e662bc50SJakub Kruzik   PetscFunctionBegin;
19270f9219eSJakub Kruzik   /* possibly allows W' = Wt (which is valid but not tested) */
1939566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)W));
194e662bc50SJakub Kruzik   if (transpose) {
1959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->Wt));
196e662bc50SJakub Kruzik     def->Wt = W;
197e662bc50SJakub Kruzik   } else {
1989566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->W));
199e662bc50SJakub Kruzik     def->W = W;
200e662bc50SJakub Kruzik   }
201e662bc50SJakub Kruzik   PetscFunctionReturn(0);
202e662bc50SJakub Kruzik }
203e662bc50SJakub Kruzik 
204e662bc50SJakub Kruzik /*@
2057e9012c5SJakub Kruzik    PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose).
206e662bc50SJakub Kruzik 
207e26ad46dSJakub Kruzik    Logically Collective
208e662bc50SJakub Kruzik 
209e662bc50SJakub Kruzik    Input Parameters:
210e662bc50SJakub Kruzik +  pc        - the preconditioner context
211e662bc50SJakub Kruzik .  W         - deflation matrix
212e26ad46dSJakub Kruzik -  transpose - indicates that W is an explicit transpose of the deflation matrix
213e26ad46dSJakub Kruzik 
214e26ad46dSJakub Kruzik    Notes:
215*f1580f4eSBarry Smith     Setting W as a multipliplicative `MATCOMPOSITE` enables use of the multilevel
216e26ad46dSJakub Kruzik     deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and
217e26ad46dSJakub Kruzik     the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with
218e26ad46dSJakub Kruzik     W1 as the deflation matrix. This repeats until the maximum level set by
21993b79a5aSJakub Kruzik     PCDeflationSetLevels() is reached or there are no more matrices available.
220e26ad46dSJakub Kruzik     If there are matrices left after reaching the maximum level,
221e26ad46dSJakub Kruzik     they are merged into a deflation matrix ...*W{n-1}*Wn.
222e662bc50SJakub Kruzik 
223e662bc50SJakub Kruzik    Level: intermediate
224e662bc50SJakub Kruzik 
225*f1580f4eSBarry Smith .seealso: `PCDeflationSetLevels()`, `PCDEFLATION`, `PCDeflationSetProjectionNullSpaceMat()`
226e662bc50SJakub Kruzik @*/
2279371c9d4SSatish Balay PetscErrorCode PCDeflationSetSpace(PC pc, Mat W, PetscBool transpose) {
228e662bc50SJakub Kruzik   PetscFunctionBegin;
229e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
230e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W, MAT_CLASSID, 2);
231e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc, transpose, 3);
232cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetSpace_C", (PC, Mat, PetscBool), (pc, W, transpose));
233e662bc50SJakub Kruzik   PetscFunctionReturn(0);
234e662bc50SJakub Kruzik }
235e662bc50SJakub Kruzik 
2369371c9d4SSatish Balay static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc, Mat mat) {
2373cf3a049SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
2383cf3a049SJakub Kruzik 
2393cf3a049SJakub Kruzik   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtA));
2423cf3a049SJakub Kruzik   def->WtA = mat;
2439566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)def->WtA));
2443cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2453cf3a049SJakub Kruzik }
2463cf3a049SJakub Kruzik 
2473cf3a049SJakub Kruzik /*@
248e26ad46dSJakub Kruzik    PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A).
2493cf3a049SJakub Kruzik 
250e26ad46dSJakub Kruzik    Collective
2513cf3a049SJakub Kruzik 
2523cf3a049SJakub Kruzik    Input Parameters:
2533cf3a049SJakub Kruzik +  pc  - preconditioner context
2543cf3a049SJakub Kruzik -  mat - projection null space matrix
2553cf3a049SJakub Kruzik 
2563cf3a049SJakub Kruzik    Level: developer
2573cf3a049SJakub Kruzik 
258*f1580f4eSBarry Smith .seealso: `PCDEFLATION`, `PCDeflationSetSpace()`
2593cf3a049SJakub Kruzik @*/
2609371c9d4SSatish Balay PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc, Mat mat) {
2613cf3a049SJakub Kruzik   PetscFunctionBegin;
2623cf3a049SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2633cf3a049SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 2);
264cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetProjectionNullSpaceMat_C", (PC, Mat), (pc, mat));
2653cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2663cf3a049SJakub Kruzik }
2673cf3a049SJakub Kruzik 
2689371c9d4SSatish Balay static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc, Mat mat) {
269e924b002SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
270e924b002SJakub Kruzik 
271e924b002SJakub Kruzik   PetscFunctionBegin;
2729566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtAW));
274e924b002SJakub Kruzik   def->WtAW = mat;
2759566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)def->WtAW));
276e924b002SJakub Kruzik   PetscFunctionReturn(0);
277e924b002SJakub Kruzik }
278e924b002SJakub Kruzik 
279e924b002SJakub Kruzik /*@
280*f1580f4eSBarry Smith    PCDeflationSetCoarseMat - Set the coarse problem `Mat`.
281e924b002SJakub Kruzik 
282e26ad46dSJakub Kruzik    Collective
283e924b002SJakub Kruzik 
284e924b002SJakub Kruzik    Input Parameters:
285e924b002SJakub Kruzik +  pc  - preconditioner context
286e924b002SJakub Kruzik -  mat - coarse problem mat
287e924b002SJakub Kruzik 
288e924b002SJakub Kruzik    Level: developer
289e924b002SJakub Kruzik 
290*f1580f4eSBarry Smith .seealso: `PCDEFLATION`, `PCDeflationGetCoarseKSP()`
291e924b002SJakub Kruzik @*/
2929371c9d4SSatish Balay PetscErrorCode PCDeflationSetCoarseMat(PC pc, Mat mat) {
293e924b002SJakub Kruzik   PetscFunctionBegin;
294e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
295e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 2);
296cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetCoarseMat_C", (PC, Mat), (pc, mat));
297e924b002SJakub Kruzik   PetscFunctionReturn(0);
298e924b002SJakub Kruzik }
299e924b002SJakub Kruzik 
3009371c9d4SSatish Balay static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc, KSP *ksp) {
3015c0c31c5SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
3025c0c31c5SJakub Kruzik 
3035c0c31c5SJakub Kruzik   PetscFunctionBegin;
30498d6e3deSJakub Kruzik   *ksp = def->WtAWinv;
3055c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
3065c0c31c5SJakub Kruzik }
3075c0c31c5SJakub Kruzik 
3085c0c31c5SJakub Kruzik /*@
309*f1580f4eSBarry Smith    PCDeflationGetCoarseKSP - Returns the coarse problem `KSP`.
3105c0c31c5SJakub Kruzik 
31198d6e3deSJakub Kruzik    Not Collective
3125c0c31c5SJakub Kruzik 
313*f1580f4eSBarry Smith    Input Parameter:
31498d6e3deSJakub Kruzik .  pc - preconditioner context
3155c0c31c5SJakub Kruzik 
316*f1580f4eSBarry Smith    Output Parameter:
317*f1580f4eSBarry Smith .  ksp - coarse problem `KSP` context
3185c0c31c5SJakub Kruzik 
319e26ad46dSJakub Kruzik    Level: advanced
32098d6e3deSJakub Kruzik 
321*f1580f4eSBarry Smith .seealso: `PCDEFLATION`, `PCDeflationSetCoarseMat()`
3225c0c31c5SJakub Kruzik @*/
3239371c9d4SSatish Balay PetscErrorCode PCDeflationGetCoarseKSP(PC pc, KSP *ksp) {
3245c0c31c5SJakub Kruzik   PetscFunctionBegin;
32522b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
32698d6e3deSJakub Kruzik   PetscValidPointer(ksp, 2);
327cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationGetCoarseKSP_C", (PC, KSP *), (pc, ksp));
32898d6e3deSJakub Kruzik   PetscFunctionReturn(0);
32998d6e3deSJakub Kruzik }
33098d6e3deSJakub Kruzik 
3319371c9d4SSatish Balay static PetscErrorCode PCDeflationGetPC_Deflation(PC pc, PC *apc) {
332268c6673SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
333268c6673SJakub Kruzik 
334268c6673SJakub Kruzik   PetscFunctionBegin;
335268c6673SJakub Kruzik   *apc = def->pc;
336268c6673SJakub Kruzik   PetscFunctionReturn(0);
337268c6673SJakub Kruzik }
338268c6673SJakub Kruzik 
339268c6673SJakub Kruzik /*@
3407e9012c5SJakub Kruzik    PCDeflationGetPC - Returns the additional preconditioner M^{-1}.
341268c6673SJakub Kruzik 
342268c6673SJakub Kruzik    Not Collective
343268c6673SJakub Kruzik 
344*f1580f4eSBarry Smith    Input Parameter:
345268c6673SJakub Kruzik .  pc  - the preconditioner context
346268c6673SJakub Kruzik 
347*f1580f4eSBarry Smith    Output Parameter:
348268c6673SJakub Kruzik .  apc - additional preconditioner
349268c6673SJakub Kruzik 
350268c6673SJakub Kruzik    Level: advanced
351268c6673SJakub Kruzik 
352*f1580f4eSBarry Smith .seealso: `PCDEFLATION`, `PCDeflationGetCoarseKSP()`
353268c6673SJakub Kruzik @*/
3549371c9d4SSatish Balay PetscErrorCode PCDeflationGetPC(PC pc, PC *apc) {
355268c6673SJakub Kruzik   PetscFunctionBegin;
356268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
357064a246eSJacob Faibussowitsch   PetscValidPointer(pc, 1);
358cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationGetPC_C", (PC, PC *), (pc, apc));
359268c6673SJakub Kruzik   PetscFunctionReturn(0);
360268c6673SJakub Kruzik }
361268c6673SJakub Kruzik 
36237eeb815SJakub Kruzik /*
363b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
364b27c8092SJakub Kruzik */
3659371c9d4SSatish Balay static PetscErrorCode PCPreSolve_Deflation(PC pc, KSP ksp, Vec b, Vec x) {
366b27c8092SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
367b27c8092SJakub Kruzik   Mat           A;
368b27c8092SJakub Kruzik   Vec           r, w1, w2;
369b27c8092SJakub Kruzik   PetscBool     nonzero;
37037eeb815SJakub Kruzik 
371b27c8092SJakub Kruzik   PetscFunctionBegin;
372b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
373b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
374b27c8092SJakub Kruzik   r  = def->work;
3759566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
37637eeb815SJakub Kruzik 
3779566063dSJacob Faibussowitsch   PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero));
3789566063dSJacob Faibussowitsch   PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
379b27c8092SJakub Kruzik   if (nonzero) {
3809566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, r)); /*    r  <- b - Ax              */
3819566063dSJacob Faibussowitsch     PetscCall(VecAYPX(r, -1.0, b));
382b27c8092SJakub Kruzik   } else {
3839566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, r)); /*    r  <- b (x is 0)          */
384b27c8092SJakub Kruzik   }
385b27c8092SJakub Kruzik 
386a3177931SJakub Kruzik   if (def->Wt) {
3879566063dSJacob Faibussowitsch     PetscCall(MatMult(def->Wt, r, w1)); /*    w1 <- W'*r                */
388a3177931SJakub Kruzik   } else {
3899566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(def->W, r, w1)); /*    w1 <- W'*r                */
390a3177931SJakub Kruzik   }
3919566063dSJacob Faibussowitsch   PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /*    w2 <- (W'*A*W)^{-1}*w1    */
3929566063dSJacob Faibussowitsch   PetscCall(MatMult(def->W, w2, r));         /*    r  <- W*w2                */
3939566063dSJacob Faibussowitsch   PetscCall(VecAYPX(x, 1.0, r));
394b27c8092SJakub Kruzik   PetscFunctionReturn(0);
395b27c8092SJakub Kruzik }
39637eeb815SJakub Kruzik 
397f8bf229cSJakub Kruzik /*
398f8bf229cSJakub Kruzik   if (def->correct) {
399ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
400f8bf229cSJakub Kruzik   } else {
401ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
402f8bf229cSJakub Kruzik   }
40337eeb815SJakub Kruzik */
4049371c9d4SSatish Balay static PetscErrorCode PCApply_Deflation(PC pc, Vec r, Vec z) {
405f8bf229cSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
406f8bf229cSJakub Kruzik   Mat           A;
407f8bf229cSJakub Kruzik   Vec           u, w1, w2;
408f8bf229cSJakub Kruzik 
409f8bf229cSJakub Kruzik   PetscFunctionBegin;
410f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
411f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
412f8bf229cSJakub Kruzik   u  = def->work;
4139566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
414f8bf229cSJakub Kruzik 
4159566063dSJacob Faibussowitsch   PetscCall(PCApply(def->pc, r, z)); /*    z <- M^{-1}*r             */
41622b0793eSJakub Kruzik   if (!def->init) {
4179566063dSJacob Faibussowitsch     PetscCall(MatMult(def->WtA, z, w1)); /*    w1 <- W'*A*z              */
418f8bf229cSJakub Kruzik     if (def->correct) {
419ae029463SJakub Kruzik       if (def->Wt) {
4209566063dSJacob Faibussowitsch         PetscCall(MatMult(def->Wt, r, w2)); /*    w2 <- W'*r                */
421ae029463SJakub Kruzik       } else {
4229566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(def->W, r, w2)); /*    w2 <- W'*r                */
423f8bf229cSJakub Kruzik       }
4249566063dSJacob Faibussowitsch       PetscCall(VecAXPY(w1, -1.0 * def->correctfact, w2)); /*    w1 <- w1 - l*w2           */
425f8bf229cSJakub Kruzik     }
4269566063dSJacob Faibussowitsch     PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /*    w2 <- (W'*A*W)^{-1}*w1    */
4279566063dSJacob Faibussowitsch     PetscCall(MatMult(def->W, w2, u));         /*    u  <- W*w2                */
4289566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, -1.0, u));            /*    z  <- z - u               */
42922b0793eSJakub Kruzik   }
430f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
431f8bf229cSJakub Kruzik }
432f8bf229cSJakub Kruzik 
4339371c9d4SSatish Balay static PetscErrorCode PCSetUp_Deflation(PC pc) {
434409a357bSJakub Kruzik   PC_Deflation    *def = (PC_Deflation *)pc->data;
435409a357bSJakub Kruzik   KSP              innerksp;
436409a357bSJakub Kruzik   PC               pcinner;
437409a357bSJakub Kruzik   Mat              Amat, nextDef = NULL, *mats;
4387b3faf33SJakub Kruzik   PetscInt         i, m, red, size;
4397b3faf33SJakub Kruzik   PetscMPIInt      commsize;
440b94d7dedSBarry Smith   PetscBool        match, flgspd, isset, transp = PETSC_FALSE;
441409a357bSJakub Kruzik   MatCompositeType ctype;
442409a357bSJakub Kruzik   MPI_Comm         comm;
4436c93e71cSJakub Kruzik   char             prefix[128] = "";
44437eeb815SJakub Kruzik 
44537eeb815SJakub Kruzik   PetscFunctionBegin;
446409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
4489566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &Amat));
44948a46eb9SPierre Jolivet   if (!def->lvl && !def->prefix) PetscCall(PCGetOptionsPrefix(pc, &def->prefix));
45048a46eb9SPierre Jolivet   if (def->lvl) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%d_", (int)def->lvl));
451f8bf229cSJakub Kruzik 
452f8bf229cSJakub Kruzik   /* compute a deflation space */
453409a357bSJakub Kruzik   if (def->W || def->Wt) {
454409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
455409a357bSJakub Kruzik   } else {
4569566063dSJacob Faibussowitsch     PetscCall(PCDeflationComputeSpace(pc));
45737eeb815SJakub Kruzik   }
45837eeb815SJakub Kruzik 
459409a357bSJakub Kruzik   /* nested deflation */
460409a357bSJakub Kruzik   if (def->W) {
4619566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)def->W, MATCOMPOSITE, &match));
462409a357bSJakub Kruzik     if (match) {
4639566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetType(def->W, &ctype));
4649566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetNumberMat(def->W, &size));
46537eeb815SJakub Kruzik     }
466409a357bSJakub Kruzik   } else {
4679566063dSJacob Faibussowitsch     PetscCall(MatCreateHermitianTranspose(def->Wt, &def->W));
4689566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)def->Wt, MATCOMPOSITE, &match));
469409a357bSJakub Kruzik     if (match) {
4709566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetType(def->Wt, &ctype));
4719566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetNumberMat(def->Wt, &size));
472409a357bSJakub Kruzik     }
473409a357bSJakub Kruzik     transp = PETSC_TRUE;
474409a357bSJakub Kruzik   }
475409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
476409a357bSJakub Kruzik     if (!transp) {
4776c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
4789566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(size, &mats));
47948a46eb9SPierre Jolivet         for (i = 0; i < size; i++) PetscCall(MatCompositeGetMat(def->W, i, &mats[i]));
480409a357bSJakub Kruzik         size -= 1;
4819566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&def->W));
482409a357bSJakub Kruzik         def->W = mats[size];
4839566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)mats[size]));
484409a357bSJakub Kruzik         if (size > 1) {
4859566063dSJacob Faibussowitsch           PetscCall(MatCreateComposite(comm, size, mats, &nextDef));
4869566063dSJacob Faibussowitsch           PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE));
487409a357bSJakub Kruzik         } else {
488409a357bSJakub Kruzik           nextDef = mats[0];
4899566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)mats[0]));
490409a357bSJakub Kruzik         }
4919566063dSJacob Faibussowitsch         PetscCall(PetscFree(mats));
492409a357bSJakub Kruzik       } else {
4931fdb00f9SJakub Kruzik         /* TODO test merge side performance */
4949566063dSJacob Faibussowitsch         /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */
4959566063dSJacob Faibussowitsch         PetscCall(MatCompositeMerge(def->W));
496409a357bSJakub Kruzik       }
49728da0170SJakub Kruzik     } else {
4986c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
4999566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(size, &mats));
50048a46eb9SPierre Jolivet         for (i = 0; i < size; i++) PetscCall(MatCompositeGetMat(def->Wt, i, &mats[i]));
50128da0170SJakub Kruzik         size -= 1;
5029566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&def->Wt));
50328da0170SJakub Kruzik         def->Wt = mats[0];
5049566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)mats[0]));
50528da0170SJakub Kruzik         if (size > 1) {
5069566063dSJacob Faibussowitsch           PetscCall(MatCreateComposite(comm, size, &mats[1], &nextDef));
5079566063dSJacob Faibussowitsch           PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE));
50828da0170SJakub Kruzik         } else {
50928da0170SJakub Kruzik           nextDef = mats[1];
5109566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)mats[1]));
511409a357bSJakub Kruzik         }
5129566063dSJacob Faibussowitsch         PetscCall(PetscFree(mats));
51328da0170SJakub Kruzik       } else {
5149566063dSJacob Faibussowitsch         /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */
5159566063dSJacob Faibussowitsch         PetscCall(MatCompositeMerge(def->Wt));
51628da0170SJakub Kruzik       }
51728da0170SJakub Kruzik     }
51828da0170SJakub Kruzik   }
51928da0170SJakub Kruzik 
52028da0170SJakub Kruzik   if (transp) {
5219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->W));
5229566063dSJacob Faibussowitsch     PetscCall(MatHermitianTranspose(def->Wt, MAT_INITIAL_MATRIX, &def->W));
523409a357bSJakub Kruzik   }
524409a357bSJakub Kruzik 
525ae029463SJakub Kruzik   /* assemble WtA */
526ae029463SJakub Kruzik   if (!def->WtA) {
527ae029463SJakub Kruzik     if (def->Wt) {
5289566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
529ae029463SJakub Kruzik     } else {
530a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
5319566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(def->W, MAT_INITIAL_MATRIX, &def->Wt));
5329566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
533a3177931SJakub Kruzik #else
5349566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(def->W, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
535a3177931SJakub Kruzik #endif
536ae029463SJakub Kruzik     }
537ae029463SJakub Kruzik   }
538409a357bSJakub Kruzik   /* setup coarse problem */
539409a357bSJakub Kruzik   if (!def->WtAWinv) {
5409566063dSJacob Faibussowitsch     PetscCall(MatGetSize(def->W, NULL, &m));
541409a357bSJakub Kruzik     if (!def->WtAW) {
5429566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->WtA, def->W, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtAW));
543409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
544b94d7dedSBarry Smith       PetscCall(MatIsSPDKnown(Amat, &isset, &flgspd));
545b94d7dedSBarry Smith       if (isset) PetscCall(MatSetOption(def->WtAW, MAT_SPD, flgspd));
54676bd3646SJed Brown       if (PetscDefined(USE_DEBUG)) {
547ae029463SJakub Kruzik         /* Check columns of W are not in kernel of A */
548409a357bSJakub Kruzik         PetscReal *norms;
549b94d7dedSBarry Smith 
5509566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m, &norms));
5519566063dSJacob Faibussowitsch         PetscCall(MatGetColumnNorms(def->WtAW, NORM_INFINITY, norms));
552ad540459SPierre Jolivet         for (i = 0; i < m; i++) PetscCheck(norms[i] > 100 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)def->WtAW), PETSC_ERR_SUP, "Column %" PetscInt_FMT " of W is in kernel of A.", i);
5539566063dSJacob Faibussowitsch         PetscCall(PetscFree(norms));
554377006c5SJakub Kruzik       }
555b94d7dedSBarry Smith     } else PetscCall(MatIsSPDKnown(def->WtAW, &isset, &flgspd));
556b94d7dedSBarry Smith 
5571fdb00f9SJakub Kruzik     /* TODO use MATINV ? */
5589566063dSJacob Faibussowitsch     PetscCall(KSPCreate(comm, &def->WtAWinv));
5599566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(def->WtAWinv, def->WtAW, def->WtAW));
5609566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(def->WtAWinv, &pcinner));
561557a2f7dSJakub Kruzik     /* Setup KSP and PC */
562557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
563557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
564557a2f7dSJakub Kruzik       /* set default KSPtype */
565557a2f7dSJakub Kruzik       if (!def->ksptype) {
566557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
567b94d7dedSBarry Smith         if (isset && flgspd) { /* SPD system */
568557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
569557a2f7dSJakub Kruzik         }
570557a2f7dSJakub Kruzik       }
5719566063dSJacob Faibussowitsch       PetscCall(KSPSetType(innerksp, def->ksptype)); /* TODO iherit from KSP + tolerances */
5729566063dSJacob Faibussowitsch       PetscCall(PCSetType(pcinner, PCDEFLATION));    /* TODO create coarse preconditinoner M_c = WtMW ? */
5739566063dSJacob Faibussowitsch       PetscCall(PCDeflationSetSpace(pcinner, nextDef, transp));
5749566063dSJacob Faibussowitsch       PetscCall(PCDeflationSetLevels_Deflation(pcinner, def->lvl + 1, def->maxlvl));
575ae029463SJakub Kruzik       /* inherit options */
5766c93e71cSJakub Kruzik       if (def->prefix) ((PC_Deflation *)(pcinner->data))->prefix = def->prefix;
5779f604af8SJakub Kruzik       ((PC_Deflation *)(pcinner->data))->init          = def->init;
578557a2f7dSJakub Kruzik       ((PC_Deflation *)(pcinner->data))->ksptype       = def->ksptype;
579557a2f7dSJakub Kruzik       ((PC_Deflation *)(pcinner->data))->correct       = def->correct;
5809f604af8SJakub Kruzik       ((PC_Deflation *)(pcinner->data))->correctfact   = def->correctfact;
5819f604af8SJakub Kruzik       ((PC_Deflation *)(pcinner->data))->reductionfact = def->reductionfact;
5829566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&nextDef));
583557a2f7dSJakub Kruzik     } else { /* the last level */
5849566063dSJacob Faibussowitsch       PetscCall(KSPSetType(def->WtAWinv, KSPPREONLY));
5859566063dSJacob Faibussowitsch       PetscCall(PCSetType(pcinner, PCTELESCOPE));
5866c93e71cSJakub Kruzik       /* do not overwrite PCTELESCOPE */
5871baa6e33SBarry Smith       if (def->prefix) PetscCall(KSPSetOptionsPrefix(def->WtAWinv, def->prefix));
5889566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(def->WtAWinv, "deflation_tel_"));
5899566063dSJacob Faibussowitsch       PetscCall(PCSetFromOptions(pcinner));
5909566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pcinner, PCTELESCOPE, &match));
59128b400f6SJacob Faibussowitsch       PetscCheck(match, comm, PETSC_ERR_SUP, "User can not owerwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead.");
592409a357bSJakub Kruzik       /* Reduction factor choice */
593409a357bSJakub Kruzik       red = def->reductionfact;
594409a357bSJakub Kruzik       if (red < 0) {
5959566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm, &commsize));
596faa75363SBarry Smith         red = PetscCeilInt(commsize, PetscCeilInt(m, commsize));
5979566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompareAny((PetscObject)(def->WtAW), &match, MATSEQDENSE, MATMPIDENSE, MATDENSE, ""));
598409a357bSJakub Kruzik         if (match) red = commsize;
59963a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "Auto choosing reduction factor %" PetscInt_FMT "\n", red));
600409a357bSJakub Kruzik       }
6019566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSetReductionFactor(pcinner, red));
6029566063dSJacob Faibussowitsch       PetscCall(PCSetUp(pcinner));
6039566063dSJacob Faibussowitsch       PetscCall(PCTelescopeGetKSP(pcinner, &innerksp));
604ac66f006SJakub Kruzik       if (innerksp) {
6059566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(innerksp, &pcinner));
6069566063dSJacob Faibussowitsch         PetscCall(PCSetType(pcinner, PCLU));
607481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU)
6089566063dSJacob Faibussowitsch         PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU, MAT_FACTOR_LU, &match));
6091baa6e33SBarry Smith         if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU));
610481b7641SJakub Kruzik #endif
611481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU_DIST)
6129566063dSJacob Faibussowitsch         PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &match));
6131baa6e33SBarry Smith         if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU_DIST));
614481b7641SJakub Kruzik #endif
615409a357bSJakub Kruzik       }
616557a2f7dSJakub Kruzik     }
617557a2f7dSJakub Kruzik 
618557a2f7dSJakub Kruzik     if (innerksp) {
6196c93e71cSJakub Kruzik       if (def->prefix) {
6209566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(innerksp, def->prefix));
6219566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(innerksp, "deflation_"));
6226c93e71cSJakub Kruzik       } else {
6239566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(innerksp, "deflation_"));
6246c93e71cSJakub Kruzik       }
6259566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(innerksp, prefix));
6269566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(innerksp));
6279566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(innerksp));
628ac66f006SJakub Kruzik     }
629409a357bSJakub Kruzik   }
6309566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(def->WtAWinv));
6319566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(def->WtAWinv));
632409a357bSJakub Kruzik 
6331fdb00f9SJakub Kruzik   /* create preconditioner */
63422b0793eSJakub Kruzik   if (!def->pc) {
6359566063dSJacob Faibussowitsch     PetscCall(PCCreate(comm, &def->pc));
6369566063dSJacob Faibussowitsch     PetscCall(PCSetOperators(def->pc, Amat, Amat));
6379566063dSJacob Faibussowitsch     PetscCall(PCSetType(def->pc, PCNONE));
6381baa6e33SBarry Smith     if (def->prefix) PetscCall(PCSetOptionsPrefix(def->pc, def->prefix));
6399566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, "deflation_"));
6409566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, prefix));
6419566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, "pc_"));
6429566063dSJacob Faibussowitsch     PetscCall(PCSetFromOptions(def->pc));
6439566063dSJacob Faibussowitsch     PetscCall(PCSetUp(def->pc));
64422b0793eSJakub Kruzik   }
64522b0793eSJakub Kruzik 
646f8bf229cSJakub Kruzik   /* create work vecs */
6479566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Amat, NULL, &def->work));
6489566063dSJacob Faibussowitsch   PetscCall(KSPCreateVecs(def->WtAWinv, 2, &def->workcoarse, 0, NULL));
64937eeb815SJakub Kruzik   PetscFunctionReturn(0);
65037eeb815SJakub Kruzik }
651b30d4775SJakub Kruzik 
6529371c9d4SSatish Balay static PetscErrorCode PCReset_Deflation(PC pc) {
653b30d4775SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
65437eeb815SJakub Kruzik 
65537eeb815SJakub Kruzik   PetscFunctionBegin;
6569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&def->work));
6579566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2, &def->workcoarse));
6589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->W));
6599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->Wt));
6609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtA));
6619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtAW));
6629566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&def->WtAWinv));
6639566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&def->pc));
66437eeb815SJakub Kruzik   PetscFunctionReturn(0);
66537eeb815SJakub Kruzik }
66637eeb815SJakub Kruzik 
6679371c9d4SSatish Balay static PetscErrorCode PCDestroy_Deflation(PC pc) {
66837eeb815SJakub Kruzik   PetscFunctionBegin;
6699566063dSJacob Faibussowitsch   PetscCall(PCReset_Deflation(pc));
6702e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", NULL));
6712e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", NULL));
6722e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", NULL));
6732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", NULL));
6742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", NULL));
6752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", NULL));
6762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", NULL));
6772e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", NULL));
6782e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", NULL));
6792e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", NULL));
6809566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
68137eeb815SJakub Kruzik   PetscFunctionReturn(0);
68237eeb815SJakub Kruzik }
68337eeb815SJakub Kruzik 
6849371c9d4SSatish Balay static PetscErrorCode PCView_Deflation(PC pc, PetscViewer viewer) {
6858513960bSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
6861ac6250aSJakub Kruzik   PetscInt      its;
6878513960bSJakub Kruzik   PetscBool     iascii;
6888513960bSJakub Kruzik 
6898513960bSJakub Kruzik   PetscFunctionBegin;
6909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
6918513960bSJakub Kruzik   if (iascii) {
69248a46eb9SPierre Jolivet     if (def->correct) PetscCall(PetscViewerASCIIPrintf(viewer, "using CP correction, factor = %g+%gi\n", (double)PetscRealPart(def->correctfact), (double)PetscImaginaryPart(def->correctfact)));
69348a46eb9SPierre Jolivet     if (!def->lvl) PetscCall(PetscViewerASCIIPrintf(viewer, "deflation space type: %s\n", PCDeflationSpaceTypes[def->spacetype]));
6948513960bSJakub Kruzik 
6959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "--- Additional PC:\n"));
6969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
6979566063dSJacob Faibussowitsch     PetscCall(PCView(def->pc, viewer));
6989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
69922b0793eSJakub Kruzik 
7009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "--- Coarse problem solver:\n"));
7019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
7029566063dSJacob Faibussowitsch     PetscCall(KSPGetTotalIterations(def->WtAWinv, &its));
70363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations: %" PetscInt_FMT "\n", its));
7049566063dSJacob Faibussowitsch     PetscCall(KSPView(def->WtAWinv, viewer));
7059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
7068513960bSJakub Kruzik   }
7078513960bSJakub Kruzik   PetscFunctionReturn(0);
7088513960bSJakub Kruzik }
7098513960bSJakub Kruzik 
7109371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Deflation(PC pc, PetscOptionItems *PetscOptionsObject) {
711880d05ceSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
71237eeb815SJakub Kruzik 
71337eeb815SJakub Kruzik   PetscFunctionBegin;
714d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Deflation options");
7159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_init_only", "Use only initialization step - Initdef", "PCDeflationSetInitOnly", def->init, &def->init, NULL));
7169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_levels", "Maximum of deflation levels", "PCDeflationSetLevels", def->maxlvl, &def->maxlvl, NULL));
7179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_reduction_factor", "Reduction factor for coarse problem solution using PCTELESCOPE", "PCDeflationSetReductionFactor", def->reductionfact, &def->reductionfact, NULL));
7189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_correction", "Add coarse problem correction Q to P", "PCDeflationSetCorrectionFactor", def->correct, &def->correct, NULL));
7199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsScalar("-pc_deflation_correction_factor", "Set multiple of Q to use as coarse problem correction", "PCDeflationSetCorrectionFactor", def->correctfact, &def->correctfact, NULL));
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_deflation_compute_space", "Compute deflation space", "PCDeflationSetSpace", PCDeflationSpaceTypes, (PetscEnum)def->spacetype, (PetscEnum *)&def->spacetype, NULL));
7219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_compute_space_size", "Set size of the deflation space to compute", "PCDeflationSetSpace", def->spacesize, &def->spacesize, NULL));
7229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_space_extend", "Extend deflation space instead of truncating (wavelets)", "PCDeflation", def->extendsp, &def->extendsp, NULL));
723d0609cedSBarry Smith   PetscOptionsHeadEnd();
72437eeb815SJakub Kruzik   PetscFunctionReturn(0);
72537eeb815SJakub Kruzik }
72637eeb815SJakub Kruzik 
72737eeb815SJakub Kruzik /*MC
728e26ad46dSJakub Kruzik    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
72937eeb815SJakub Kruzik 
730e26ad46dSJakub Kruzik    Options Database Keys:
731e26ad46dSJakub Kruzik +    -pc_deflation_init_only          <false> - if true computes only the special guess
732e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
73371f2caa7Sprj- .    -pc_deflation_reduction_factor <\-1>     - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
734e26ad46dSJakub Kruzik .    -pc_deflation_correction         <false> - if true apply coarse problem correction
735e26ad46dSJakub Kruzik .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
736e26ad46dSJakub Kruzik .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
737e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
73837eeb815SJakub Kruzik 
73937eeb815SJakub Kruzik    Notes:
7407e9012c5SJakub Kruzik     Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
741e26ad46dSJakub Kruzik     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
742e26ad46dSJakub Kruzik 
743e26ad46dSJakub Kruzik     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
744*f1580f4eSBarry Smith     If `PCDeflationSetInitOnly()` or -pc_deflation_init_only is set to `PETSC_TRUE` (InitDef scheme), the application of the
745e26ad46dSJakub Kruzik     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
7461fdb00f9SJakub Kruzik     application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
747e26ad46dSJakub Kruzik     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
748e26ad46dSJakub Kruzik     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
749e26ad46dSJakub Kruzik     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
750e26ad46dSJakub Kruzik 
751e26ad46dSJakub Kruzik     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
752*f1580f4eSBarry Smith     be controlled by `PCDeflationSetSpaceToCompute()` or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
753*f1580f4eSBarry Smith     User can set an arbitrary deflation space matrix with `PCDeflationSetSpace()`. If the deflation matrix
754*f1580f4eSBarry Smith     is a multiplicative `MATCOMPOSITE`, a multilevel deflation [3] is used. The first matrix in the composite is used as the
755*f1580f4eSBarry Smith     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by `KSPFCG` (if A is `MAT_SPD`) or `KSPFGMRES` preconditioned
756*f1580f4eSBarry Smith     by deflation with deflation matrix being the next matrix in the `MATCOMPOSITE`. This scheme repeats until the maximum
7571fdb00f9SJakub Kruzik     level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged
7581fdb00f9SJakub Kruzik     (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by
759*f1580f4eSBarry Smith     `PCDeflationSetLevels()` or by -pc_deflation_levels.
760e26ad46dSJakub Kruzik 
761*f1580f4eSBarry Smith     The coarse problem `KSP` can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1]
7626c93e71cSJakub Kruzik     from the second level onward. You can also use
763*f1580f4eSBarry Smith     `PCDeflationGetCoarseKSP()` to control it from code. The bottom level KSP defaults to
764*f1580f4eSBarry Smith     `KSPPREONLY` with `PCLU` direct solver (`MATSOLVERSUPERLU`/`MATSOLVERSUPERLU_DIST` if available) wrapped into `PCTELESCOPE`.
765*f1580f4eSBarry Smith     For convenience, the reduction factor can be set by `PCDeflationSetReductionFactor()`
766481b7641SJakub Kruzik     or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
767e26ad46dSJakub Kruzik 
768e63c2c6dSJakub Kruzik     The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for
769*f1580f4eSBarry Smith     coarse problem `KSP` apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use
770*f1580f4eSBarry Smith     `PCDeflationGetPC()` to control the additional preconditioner from code. It defaults to `PCNONE`.
771e26ad46dSJakub Kruzik 
772e26ad46dSJakub Kruzik     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
773*f1580f4eSBarry Smith     be set by pc_deflation_correction_factor or by `PCDeflationSetCorrectionFactor()`. The coarse problem can
774e26ad46dSJakub Kruzik     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
775e26ad46dSJakub Kruzik     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
776e26ad46dSJakub Kruzik     an isolated eigenvalue.
777e26ad46dSJakub Kruzik 
7781fdb00f9SJakub Kruzik     The options are automatically inherited from the previous deflation level.
7799f604af8SJakub Kruzik 
780*f1580f4eSBarry Smith     The preconditioner supports `KSPMonitorDynamicTolerance()`. This is useful for the multilevel scheme for which we also
7811fdb00f9SJakub Kruzik     recommend limiting the number of iterations for the coarse problems.
782e26ad46dSJakub Kruzik 
783e26ad46dSJakub Kruzik     See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
784e26ad46dSJakub Kruzik     Section 4 describes some possible choices for the deflation space.
785e26ad46dSJakub Kruzik 
786e26ad46dSJakub Kruzik      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
787e26ad46dSJakub Kruzik      Academy of Sciences and VSB - TU Ostrava.
788e26ad46dSJakub Kruzik 
789e26ad46dSJakub Kruzik      Developed from PERMON code used in [4] while on a research stay with
790e26ad46dSJakub Kruzik      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
791e26ad46dSJakub Kruzik 
792e26ad46dSJakub Kruzik    References:
793606c0280SSatish Balay +  * - A. Nicolaides. "Deflation of conjugate gradients with applications to boundary value problems", SIAM J. Numer. Anal. 24.2, 1987.
794606c0280SSatish Balay .  * - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
795606c0280SSatish Balay .  * - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
796606c0280SSatish Balay -  * - J. Kruzik "Implementation of the Deflated Variants of the Conjugate Gradient Method", Master's thesis, VSB-TUO, 2018 - http://dspace5.vsb.cz/bitstream/handle/10084/130303/KRU0097_USP_N2658_2612T078_2018.pdf
797e26ad46dSJakub Kruzik 
798e26ad46dSJakub Kruzik    Level: intermediate
79937eeb815SJakub Kruzik 
800db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
801db781477SPatrick Sanan           `PCDeflationSetInitOnly()`, `PCDeflationSetLevels()`, `PCDeflationSetReductionFactor()`,
802db781477SPatrick Sanan           `PCDeflationSetCorrectionFactor()`, `PCDeflationSetSpaceToCompute()`,
803db781477SPatrick Sanan           `PCDeflationSetSpace()`, `PCDeflationSpaceType`, `PCDeflationSetProjectionNullSpaceMat()`,
804db781477SPatrick Sanan           `PCDeflationSetCoarseMat()`, `PCDeflationGetCoarseKSP()`, `PCDeflationGetPC()`
80537eeb815SJakub Kruzik M*/
80637eeb815SJakub Kruzik 
8079371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc) {
80837eeb815SJakub Kruzik   PC_Deflation *def;
80937eeb815SJakub Kruzik 
81037eeb815SJakub Kruzik   PetscFunctionBegin;
8119566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &def));
81237eeb815SJakub Kruzik   pc->data = (void *)def;
81337eeb815SJakub Kruzik 
814e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
815e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
816fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
817e662bc50SJakub Kruzik   def->reductionfact = -1;
818e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
819e662bc50SJakub Kruzik   def->spacesize     = 1;
820e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
8216c93e71cSJakub Kruzik   def->lvl           = 0;
8226c93e71cSJakub Kruzik   def->maxlvl        = 0;
82370f9219eSJakub Kruzik   def->W             = NULL;
82470f9219eSJakub Kruzik   def->Wt            = NULL;
82537eeb815SJakub Kruzik 
82637eeb815SJakub Kruzik   pc->ops->apply          = PCApply_Deflation;
827b27c8092SJakub Kruzik   pc->ops->presolve       = PCPreSolve_Deflation;
82837eeb815SJakub Kruzik   pc->ops->setup          = PCSetUp_Deflation;
82937eeb815SJakub Kruzik   pc->ops->reset          = PCReset_Deflation;
83037eeb815SJakub Kruzik   pc->ops->destroy        = PCDestroy_Deflation;
83137eeb815SJakub Kruzik   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
8328513960bSJakub Kruzik   pc->ops->view           = PCView_Deflation;
83337eeb815SJakub Kruzik 
8349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", PCDeflationSetInitOnly_Deflation));
8359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", PCDeflationSetLevels_Deflation));
8369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", PCDeflationSetReductionFactor_Deflation));
8379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", PCDeflationSetCorrectionFactor_Deflation));
8389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", PCDeflationSetSpaceToCompute_Deflation));
8399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", PCDeflationSetSpace_Deflation));
8409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", PCDeflationSetProjectionNullSpaceMat_Deflation));
8419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", PCDeflationSetCoarseMat_Deflation));
8429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", PCDeflationGetCoarseKSP_Deflation));
8439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", PCDeflationGetPC_Deflation));
84437eeb815SJakub Kruzik   PetscFunctionReturn(0);
84537eeb815SJakub Kruzik }
846