xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
1dd01b7e5SBarry Smith #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/ /* header file 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 
5d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc, PetscBool flg)
6d71ae5a4SJacob Faibussowitsch {
7a122ebaeSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
8a122ebaeSJakub Kruzik 
9a122ebaeSJakub Kruzik   PetscFunctionBegin;
10a122ebaeSJakub Kruzik   def->init = flg;
113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12a122ebaeSJakub Kruzik }
13a122ebaeSJakub Kruzik 
14a122ebaeSJakub Kruzik /*@
15a122ebaeSJakub Kruzik   PCDeflationSetInitOnly - Do only initialization step.
16e26ad46dSJakub Kruzik   Sets initial guess to the solution on the deflation space but does not apply
17e26ad46dSJakub Kruzik   the deflation preconditioner. The additional preconditioner is still applied.
18a122ebaeSJakub Kruzik 
19e26ad46dSJakub Kruzik   Logically Collective
20a122ebaeSJakub Kruzik 
21a122ebaeSJakub Kruzik   Input Parameters:
22a122ebaeSJakub Kruzik + pc  - the preconditioner context
23f1580f4eSBarry Smith - flg - default `PETSC_FALSE`
24a122ebaeSJakub Kruzik 
25f1580f4eSBarry Smith   Options Database Key:
26e26ad46dSJakub Kruzik . -pc_deflation_init_only <false> - if true computes only the special guess
27e26ad46dSJakub Kruzik 
28a122ebaeSJakub Kruzik   Level: intermediate
29a122ebaeSJakub Kruzik 
30562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`
31a122ebaeSJakub Kruzik @*/
32d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetInitOnly(PC pc, PetscBool flg)
33d71ae5a4SJacob Faibussowitsch {
34a122ebaeSJakub Kruzik   PetscFunctionBegin;
35a122ebaeSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
36a122ebaeSJakub Kruzik   PetscValidLogicalCollectiveBool(pc, flg, 2);
37cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetInitOnly_C", (PC, PetscBool), (pc, flg));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39a122ebaeSJakub Kruzik }
40a122ebaeSJakub Kruzik 
41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc, PetscInt current, PetscInt max)
42d71ae5a4SJacob Faibussowitsch {
4398d6e3deSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
4498d6e3deSJakub Kruzik 
4598d6e3deSJakub Kruzik   PetscFunctionBegin;
466c93e71cSJakub Kruzik   if (current) def->lvl = current;
476c93e71cSJakub Kruzik   def->maxlvl = max;
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4998d6e3deSJakub Kruzik }
5098d6e3deSJakub Kruzik 
5198d6e3deSJakub Kruzik /*@
5293b79a5aSJakub Kruzik   PCDeflationSetLevels - Set the maximum level of deflation nesting.
5398d6e3deSJakub Kruzik 
54e26ad46dSJakub Kruzik   Logically Collective
5598d6e3deSJakub Kruzik 
5698d6e3deSJakub Kruzik   Input Parameters:
5798d6e3deSJakub Kruzik + pc  - the preconditioner context
5898d6e3deSJakub Kruzik - max - maximum deflation level
5998d6e3deSJakub Kruzik 
60f1580f4eSBarry Smith   Options Database Key:
61e26ad46dSJakub Kruzik . -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation
62e26ad46dSJakub Kruzik 
6398d6e3deSJakub Kruzik   Level: intermediate
6498d6e3deSJakub Kruzik 
65562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDeflationSetSpaceToCompute()`, `PCDeflationSetSpace()`, `PCDEFLATION`
6698d6e3deSJakub Kruzik @*/
67d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetLevels(PC pc, PetscInt max)
68d71ae5a4SJacob Faibussowitsch {
6998d6e3deSJakub Kruzik   PetscFunctionBegin;
7098d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7198d6e3deSJakub Kruzik   PetscValidLogicalCollectiveInt(pc, max, 2);
72cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetLevels_C", (PC, PetscInt, PetscInt), (pc, 0, max));
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7498d6e3deSJakub Kruzik }
7598d6e3deSJakub Kruzik 
76d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc, PetscInt red)
77d71ae5a4SJacob Faibussowitsch {
78859a873cSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
79859a873cSJakub Kruzik 
80859a873cSJakub Kruzik   PetscFunctionBegin;
81859a873cSJakub Kruzik   def->reductionfact = red;
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83859a873cSJakub Kruzik }
84859a873cSJakub Kruzik 
85859a873cSJakub Kruzik /*@
86f1580f4eSBarry Smith   PCDeflationSetReductionFactor - Set reduction factor for the `PCDEFLATION`
87859a873cSJakub Kruzik 
88e26ad46dSJakub Kruzik   Logically Collective
89859a873cSJakub Kruzik 
90859a873cSJakub Kruzik   Input Parameters:
91859a873cSJakub Kruzik + pc  - the preconditioner context
92f1580f4eSBarry Smith - red - reduction factor (or `PETSC_DETERMINE`)
93859a873cSJakub Kruzik 
94f1580f4eSBarry Smith   Options Database Key:
95f1580f4eSBarry Smith . -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for `PCDEFLATION`
96e26ad46dSJakub Kruzik 
97f1580f4eSBarry Smith   Note:
98e26ad46dSJakub Kruzik   Default is computed based on the size of the coarse problem.
99e26ad46dSJakub Kruzik 
100859a873cSJakub Kruzik   Level: intermediate
101859a873cSJakub Kruzik 
102562efe2eSBarry Smith .seealso: [](ch_ksp), `PCTELESCOPE`, `PCDEFLATION`, `PCDeflationSetLevels()`
103859a873cSJakub Kruzik @*/
104d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetReductionFactor(PC pc, PetscInt red)
105d71ae5a4SJacob Faibussowitsch {
106859a873cSJakub Kruzik   PetscFunctionBegin;
107859a873cSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108859a873cSJakub Kruzik   PetscValidLogicalCollectiveInt(pc, red, 2);
109cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetReductionFactor_C", (PC, PetscInt), (pc, red));
1103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111859a873cSJakub Kruzik }
112859a873cSJakub Kruzik 
113d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc, PetscScalar fact)
114d71ae5a4SJacob Faibussowitsch {
1158a71cb68SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
1168a71cb68SJakub Kruzik 
1178a71cb68SJakub Kruzik   PetscFunctionBegin;
1188a71cb68SJakub Kruzik   /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */
1198a71cb68SJakub Kruzik   def->correct     = PETSC_TRUE;
1208a71cb68SJakub Kruzik   def->correctfact = fact;
121d6b1b4b3SYANG Zongze   if (def->correctfact == 0.0) def->correct = PETSC_FALSE;
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1238a71cb68SJakub Kruzik }
1248a71cb68SJakub Kruzik 
1258a71cb68SJakub Kruzik /*@
1268a71cb68SJakub Kruzik   PCDeflationSetCorrectionFactor - Set coarse problem correction factor.
127f1580f4eSBarry Smith   The preconditioner becomes P*M^{-1} + fact*Q.
1288a71cb68SJakub Kruzik 
129e26ad46dSJakub Kruzik   Logically Collective
1308a71cb68SJakub Kruzik 
1318a71cb68SJakub Kruzik   Input Parameters:
1328a71cb68SJakub Kruzik + pc   - the preconditioner context
133e26ad46dSJakub Kruzik - fact - correction factor
134e26ad46dSJakub Kruzik 
135e26ad46dSJakub Kruzik   Options Database Keys:
136e26ad46dSJakub Kruzik + -pc_deflation_correction        <false> - if true apply coarse problem correction
137e26ad46dSJakub Kruzik - -pc_deflation_correction_factor <1.0>   - sets coarse problem correction factor
138e26ad46dSJakub Kruzik 
139f1580f4eSBarry Smith   Note:
140e26ad46dSJakub Kruzik   Any non-zero fact enables the coarse problem correction.
1418a71cb68SJakub Kruzik 
1428a71cb68SJakub Kruzik   Level: intermediate
1438a71cb68SJakub Kruzik 
144562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`, `PCDeflationSetLevels()`, `PCDeflationSetReductionFactor()`
1458a71cb68SJakub Kruzik @*/
146d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetCorrectionFactor(PC pc, PetscScalar fact)
147d71ae5a4SJacob Faibussowitsch {
1488a71cb68SJakub Kruzik   PetscFunctionBegin;
1498a71cb68SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1508a71cb68SJakub Kruzik   PetscValidLogicalCollectiveScalar(pc, fact, 2);
151cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetCorrectionFactor_C", (PC, PetscScalar), (pc, fact));
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1538a71cb68SJakub Kruzik }
1548a71cb68SJakub Kruzik 
155d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc, PCDeflationSpaceType type, PetscInt size)
156d71ae5a4SJacob Faibussowitsch {
15739417d7eSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
15839417d7eSJakub Kruzik 
15939417d7eSJakub Kruzik   PetscFunctionBegin;
16039417d7eSJakub Kruzik   if (type) def->spacetype = type;
16139417d7eSJakub Kruzik   if (size > 0) def->spacesize = size;
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16339417d7eSJakub Kruzik }
16439417d7eSJakub Kruzik 
16539417d7eSJakub Kruzik /*@
16639417d7eSJakub Kruzik   PCDeflationSetSpaceToCompute - Set deflation space type and size to compute.
16739417d7eSJakub Kruzik 
168e26ad46dSJakub Kruzik   Logically Collective
16939417d7eSJakub Kruzik 
17039417d7eSJakub Kruzik   Input Parameters:
17139417d7eSJakub Kruzik + pc   - the preconditioner context
172f1580f4eSBarry Smith . type - deflation space type to compute (or `PETSC_IGNORE`)
173f1580f4eSBarry Smith - size - size of the space to compute (or `PETSC_DEFAULT`)
17439417d7eSJakub Kruzik 
175e26ad46dSJakub Kruzik   Options Database Keys:
176f1580f4eSBarry Smith + -pc_deflation_compute_space      <haar> - compute `PCDeflationSpaceType` deflation space
177e26ad46dSJakub Kruzik - -pc_deflation_compute_space_size <1>    - size of the deflation space
178e26ad46dSJakub Kruzik 
17939417d7eSJakub Kruzik   Notes:
18039417d7eSJakub Kruzik   For wavelet-based deflation, size represents number of levels.
181e26ad46dSJakub Kruzik 
182f1580f4eSBarry Smith   The deflation space is computed in `PCSetUp()`.
18339417d7eSJakub Kruzik 
18439417d7eSJakub Kruzik   Level: intermediate
18539417d7eSJakub Kruzik 
186562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDeflationSetLevels()`, `PCDEFLATION`
18739417d7eSJakub Kruzik @*/
188d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetSpaceToCompute(PC pc, PCDeflationSpaceType type, PetscInt size)
189d71ae5a4SJacob Faibussowitsch {
19039417d7eSJakub Kruzik   PetscFunctionBegin;
19139417d7eSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19239417d7eSJakub Kruzik   if (type) PetscValidLogicalCollectiveEnum(pc, type, 2);
19339417d7eSJakub Kruzik   if (size > 0) PetscValidLogicalCollectiveInt(pc, size, 3);
194cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetSpaceToCompute_C", (PC, PCDeflationSpaceType, PetscInt), (pc, type, size));
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19639417d7eSJakub Kruzik }
1978513960bSJakub Kruzik 
198d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc, Mat W, PetscBool transpose)
199d71ae5a4SJacob Faibussowitsch {
200e662bc50SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
201e662bc50SJakub Kruzik 
202e662bc50SJakub Kruzik   PetscFunctionBegin;
20370f9219eSJakub Kruzik   /* possibly allows W' = Wt (which is valid but not tested) */
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)W));
205e662bc50SJakub Kruzik   if (transpose) {
2069566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->Wt));
207e662bc50SJakub Kruzik     def->Wt = W;
208e662bc50SJakub Kruzik   } else {
2099566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->W));
210e662bc50SJakub Kruzik     def->W = W;
211e662bc50SJakub Kruzik   }
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
213e662bc50SJakub Kruzik }
214e662bc50SJakub Kruzik 
215e662bc50SJakub Kruzik /*@
2167e9012c5SJakub Kruzik   PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose).
217e662bc50SJakub Kruzik 
218e26ad46dSJakub Kruzik   Logically Collective
219e662bc50SJakub Kruzik 
220e662bc50SJakub Kruzik   Input Parameters:
221e662bc50SJakub Kruzik + pc        - the preconditioner context
222e662bc50SJakub Kruzik . W         - deflation matrix
223e26ad46dSJakub Kruzik - transpose - indicates that W is an explicit transpose of the deflation matrix
224e26ad46dSJakub Kruzik 
2252fe279fdSBarry Smith   Level: intermediate
2262fe279fdSBarry Smith 
227e26ad46dSJakub Kruzik   Notes:
228f1580f4eSBarry Smith   Setting W as a multipliplicative `MATCOMPOSITE` enables use of the multilevel
229e26ad46dSJakub Kruzik   deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and
230e26ad46dSJakub Kruzik   the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with
231e26ad46dSJakub Kruzik   W1 as the deflation matrix. This repeats until the maximum level set by
23293b79a5aSJakub Kruzik   PCDeflationSetLevels() is reached or there are no more matrices available.
233e26ad46dSJakub Kruzik   If there are matrices left after reaching the maximum level,
234e26ad46dSJakub Kruzik   they are merged into a deflation matrix ...*W{n-1}*Wn.
235e662bc50SJakub Kruzik 
236562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDeflationSetLevels()`, `PCDEFLATION`, `PCDeflationSetProjectionNullSpaceMat()`
237e662bc50SJakub Kruzik @*/
238d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetSpace(PC pc, Mat W, PetscBool transpose)
239d71ae5a4SJacob Faibussowitsch {
240e662bc50SJakub Kruzik   PetscFunctionBegin;
241e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
242e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W, MAT_CLASSID, 2);
243e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc, transpose, 3);
244cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetSpace_C", (PC, Mat, PetscBool), (pc, W, transpose));
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246e662bc50SJakub Kruzik }
247e662bc50SJakub Kruzik 
248d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc, Mat mat)
249d71ae5a4SJacob Faibussowitsch {
2503cf3a049SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
2513cf3a049SJakub Kruzik 
2523cf3a049SJakub Kruzik   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtA));
2553cf3a049SJakub Kruzik   def->WtA = mat;
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2573cf3a049SJakub Kruzik }
2583cf3a049SJakub Kruzik 
2593cf3a049SJakub Kruzik /*@
260e26ad46dSJakub Kruzik   PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A).
2613cf3a049SJakub Kruzik 
262e26ad46dSJakub Kruzik   Collective
2633cf3a049SJakub Kruzik 
2643cf3a049SJakub Kruzik   Input Parameters:
2653cf3a049SJakub Kruzik + pc  - preconditioner context
2663cf3a049SJakub Kruzik - mat - projection null space matrix
2673cf3a049SJakub Kruzik 
2683cf3a049SJakub Kruzik   Level: developer
2693cf3a049SJakub Kruzik 
270562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`, `PCDeflationSetSpace()`
2713cf3a049SJakub Kruzik @*/
272d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC pc, Mat mat)
273d71ae5a4SJacob Faibussowitsch {
2743cf3a049SJakub Kruzik   PetscFunctionBegin;
2753cf3a049SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2763cf3a049SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 2);
277cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetProjectionNullSpaceMat_C", (PC, Mat), (pc, mat));
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2793cf3a049SJakub Kruzik }
2803cf3a049SJakub Kruzik 
281d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc, Mat mat)
282d71ae5a4SJacob Faibussowitsch {
283e924b002SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
284e924b002SJakub Kruzik 
285e924b002SJakub Kruzik   PetscFunctionBegin;
2869566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtAW));
288e924b002SJakub Kruzik   def->WtAW = mat;
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
290e924b002SJakub Kruzik }
291e924b002SJakub Kruzik 
292e924b002SJakub Kruzik /*@
293f1580f4eSBarry Smith   PCDeflationSetCoarseMat - Set the coarse problem `Mat`.
294e924b002SJakub Kruzik 
295e26ad46dSJakub Kruzik   Collective
296e924b002SJakub Kruzik 
297e924b002SJakub Kruzik   Input Parameters:
298e924b002SJakub Kruzik + pc  - preconditioner context
299e924b002SJakub Kruzik - mat - coarse problem mat
300e924b002SJakub Kruzik 
301e924b002SJakub Kruzik   Level: developer
302e924b002SJakub Kruzik 
303562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`, `PCDeflationGetCoarseKSP()`
304e924b002SJakub Kruzik @*/
305d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationSetCoarseMat(PC pc, Mat mat)
306d71ae5a4SJacob Faibussowitsch {
307e924b002SJakub Kruzik   PetscFunctionBegin;
308e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
309e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 2);
310cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationSetCoarseMat_C", (PC, Mat), (pc, mat));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312e924b002SJakub Kruzik }
313e924b002SJakub Kruzik 
314d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc, KSP *ksp)
315d71ae5a4SJacob Faibussowitsch {
3165c0c31c5SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
3175c0c31c5SJakub Kruzik 
3185c0c31c5SJakub Kruzik   PetscFunctionBegin;
31998d6e3deSJakub Kruzik   *ksp = def->WtAWinv;
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3215c0c31c5SJakub Kruzik }
3225c0c31c5SJakub Kruzik 
3235c0c31c5SJakub Kruzik /*@
324f1580f4eSBarry Smith   PCDeflationGetCoarseKSP - Returns the coarse problem `KSP`.
3255c0c31c5SJakub Kruzik 
32698d6e3deSJakub Kruzik   Not Collective
3275c0c31c5SJakub Kruzik 
328f1580f4eSBarry Smith   Input Parameter:
32998d6e3deSJakub Kruzik . pc - preconditioner context
3305c0c31c5SJakub Kruzik 
331f1580f4eSBarry Smith   Output Parameter:
332f1580f4eSBarry Smith . ksp - coarse problem `KSP` context
3335c0c31c5SJakub Kruzik 
334e26ad46dSJakub Kruzik   Level: advanced
33598d6e3deSJakub Kruzik 
336562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`, `PCDeflationSetCoarseMat()`
3375c0c31c5SJakub Kruzik @*/
338d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationGetCoarseKSP(PC pc, KSP *ksp)
339d71ae5a4SJacob Faibussowitsch {
3405c0c31c5SJakub Kruzik   PetscFunctionBegin;
34122b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3424f572ea9SToby Isaac   PetscAssertPointer(ksp, 2);
343cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationGetCoarseKSP_C", (PC, KSP *), (pc, ksp));
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34598d6e3deSJakub Kruzik }
34698d6e3deSJakub Kruzik 
347d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDeflationGetPC_Deflation(PC pc, PC *apc)
348d71ae5a4SJacob Faibussowitsch {
349268c6673SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
350268c6673SJakub Kruzik 
351268c6673SJakub Kruzik   PetscFunctionBegin;
352268c6673SJakub Kruzik   *apc = def->pc;
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354268c6673SJakub Kruzik }
355268c6673SJakub Kruzik 
356268c6673SJakub Kruzik /*@
3577e9012c5SJakub Kruzik   PCDeflationGetPC - Returns the additional preconditioner M^{-1}.
358268c6673SJakub Kruzik 
359268c6673SJakub Kruzik   Not Collective
360268c6673SJakub Kruzik 
361f1580f4eSBarry Smith   Input Parameter:
362268c6673SJakub Kruzik . pc - the preconditioner context
363268c6673SJakub Kruzik 
364f1580f4eSBarry Smith   Output Parameter:
365268c6673SJakub Kruzik . apc - additional preconditioner
366268c6673SJakub Kruzik 
367268c6673SJakub Kruzik   Level: advanced
368268c6673SJakub Kruzik 
369562efe2eSBarry Smith .seealso: [](ch_ksp), `PCDEFLATION`, `PCDeflationGetCoarseKSP()`
370268c6673SJakub Kruzik @*/
371d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDeflationGetPC(PC pc, PC *apc)
372d71ae5a4SJacob Faibussowitsch {
373268c6673SJakub Kruzik   PetscFunctionBegin;
374268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3754f572ea9SToby Isaac   PetscAssertPointer(pc, 1);
376cac4c232SBarry Smith   PetscTryMethod(pc, "PCDeflationGetPC_C", (PC, PC *), (pc, apc));
3773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
378268c6673SJakub Kruzik }
379268c6673SJakub Kruzik 
38037eeb815SJakub Kruzik /*
381b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
382b27c8092SJakub Kruzik */
383d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Deflation(PC pc, KSP ksp, Vec b, Vec x)
384d71ae5a4SJacob Faibussowitsch {
385b27c8092SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
386b27c8092SJakub Kruzik   Mat           A;
387b27c8092SJakub Kruzik   Vec           r, w1, w2;
388b27c8092SJakub Kruzik   PetscBool     nonzero;
38937eeb815SJakub Kruzik 
390b27c8092SJakub Kruzik   PetscFunctionBegin;
391b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
392b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
393b27c8092SJakub Kruzik   r  = def->work;
3949566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
39537eeb815SJakub Kruzik 
3969566063dSJacob Faibussowitsch   PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero));
3979566063dSJacob Faibussowitsch   PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
398b27c8092SJakub Kruzik   if (nonzero) {
3999566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, r)); /*    r  <- b - Ax              */
4009566063dSJacob Faibussowitsch     PetscCall(VecAYPX(r, -1.0, b));
401b27c8092SJakub Kruzik   } else {
4029566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, r)); /*    r  <- b (x is 0)          */
403b27c8092SJakub Kruzik   }
404b27c8092SJakub Kruzik 
405a3177931SJakub Kruzik   if (def->Wt) {
4069566063dSJacob Faibussowitsch     PetscCall(MatMult(def->Wt, r, w1)); /*    w1 <- W'*r                */
407a3177931SJakub Kruzik   } else {
4089566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(def->W, r, w1)); /*    w1 <- W'*r                */
409a3177931SJakub Kruzik   }
4109566063dSJacob Faibussowitsch   PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /*    w2 <- (W'*A*W)^{-1}*w1    */
4119566063dSJacob Faibussowitsch   PetscCall(MatMult(def->W, w2, r));         /*    r  <- W*w2                */
4129566063dSJacob Faibussowitsch   PetscCall(VecAYPX(x, 1.0, r));
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
414b27c8092SJakub Kruzik }
41537eeb815SJakub Kruzik 
416f8bf229cSJakub Kruzik /*
417f8bf229cSJakub Kruzik   if (def->correct) {
418ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
419f8bf229cSJakub Kruzik   } else {
420ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
421f8bf229cSJakub Kruzik   }
42237eeb815SJakub Kruzik */
423d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Deflation(PC pc, Vec r, Vec z)
424d71ae5a4SJacob Faibussowitsch {
425f8bf229cSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
426f8bf229cSJakub Kruzik   Mat           A;
427f8bf229cSJakub Kruzik   Vec           u, w1, w2;
428f8bf229cSJakub Kruzik 
429f8bf229cSJakub Kruzik   PetscFunctionBegin;
430f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
431f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
432f8bf229cSJakub Kruzik   u  = def->work;
4339566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
434f8bf229cSJakub Kruzik 
4359566063dSJacob Faibussowitsch   PetscCall(PCApply(def->pc, r, z)); /*    z <- M^{-1}*r             */
43622b0793eSJakub Kruzik   if (!def->init) {
4379566063dSJacob Faibussowitsch     PetscCall(MatMult(def->WtA, z, w1)); /*    w1 <- W'*A*z              */
438f8bf229cSJakub Kruzik     if (def->correct) {
439ae029463SJakub Kruzik       if (def->Wt) {
4409566063dSJacob Faibussowitsch         PetscCall(MatMult(def->Wt, r, w2)); /*    w2 <- W'*r                */
441ae029463SJakub Kruzik       } else {
4429566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(def->W, r, w2)); /*    w2 <- W'*r                */
443f8bf229cSJakub Kruzik       }
4449566063dSJacob Faibussowitsch       PetscCall(VecAXPY(w1, -1.0 * def->correctfact, w2)); /*    w1 <- w1 - l*w2           */
445f8bf229cSJakub Kruzik     }
4469566063dSJacob Faibussowitsch     PetscCall(KSPSolve(def->WtAWinv, w1, w2)); /*    w2 <- (W'*A*W)^{-1}*w1    */
4479566063dSJacob Faibussowitsch     PetscCall(MatMult(def->W, w2, u));         /*    u  <- W*w2                */
4489566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, -1.0, u));            /*    z  <- z - u               */
44922b0793eSJakub Kruzik   }
4503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
451f8bf229cSJakub Kruzik }
452f8bf229cSJakub Kruzik 
453d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Deflation(PC pc)
454d71ae5a4SJacob Faibussowitsch {
455409a357bSJakub Kruzik   PC_Deflation    *def = (PC_Deflation *)pc->data;
456d6b1b4b3SYANG Zongze   DM               dm;
457409a357bSJakub Kruzik   KSP              innerksp;
458409a357bSJakub Kruzik   PC               pcinner;
459409a357bSJakub Kruzik   Mat              Amat, nextDef = NULL, *mats;
4607b3faf33SJakub Kruzik   PetscInt         i, m, red, size;
4617b3faf33SJakub Kruzik   PetscMPIInt      commsize;
462b94d7dedSBarry Smith   PetscBool        match, flgspd, isset, transp = PETSC_FALSE;
463409a357bSJakub Kruzik   MatCompositeType ctype;
464409a357bSJakub Kruzik   MPI_Comm         comm;
4656c93e71cSJakub Kruzik   char             prefix[128] = "";
46637eeb815SJakub Kruzik 
46737eeb815SJakub Kruzik   PetscFunctionBegin;
4683ba16761SJacob Faibussowitsch   if (pc->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
4709566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &Amat));
47148a46eb9SPierre Jolivet   if (!def->lvl && !def->prefix) PetscCall(PCGetOptionsPrefix(pc, &def->prefix));
47248a46eb9SPierre Jolivet   if (def->lvl) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%d_", (int)def->lvl));
473f8bf229cSJakub Kruzik 
474f8bf229cSJakub Kruzik   /* compute a deflation space */
475409a357bSJakub Kruzik   if (def->W || def->Wt) {
476409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
477409a357bSJakub Kruzik   } else {
4789566063dSJacob Faibussowitsch     PetscCall(PCDeflationComputeSpace(pc));
47937eeb815SJakub Kruzik   }
48037eeb815SJakub Kruzik 
481409a357bSJakub Kruzik   /* nested deflation */
482409a357bSJakub Kruzik   if (def->W) {
4839566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)def->W, MATCOMPOSITE, &match));
484409a357bSJakub Kruzik     if (match) {
4859566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetType(def->W, &ctype));
4869566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetNumberMat(def->W, &size));
48737eeb815SJakub Kruzik     }
488409a357bSJakub Kruzik   } else {
4899566063dSJacob Faibussowitsch     PetscCall(MatCreateHermitianTranspose(def->Wt, &def->W));
4909566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)def->Wt, MATCOMPOSITE, &match));
491409a357bSJakub Kruzik     if (match) {
4929566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetType(def->Wt, &ctype));
4939566063dSJacob Faibussowitsch       PetscCall(MatCompositeGetNumberMat(def->Wt, &size));
494409a357bSJakub Kruzik     }
495409a357bSJakub Kruzik     transp = PETSC_TRUE;
496409a357bSJakub Kruzik   }
497409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
498409a357bSJakub Kruzik     if (!transp) {
4996c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
5009566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(size, &mats));
50148a46eb9SPierre Jolivet         for (i = 0; i < size; i++) PetscCall(MatCompositeGetMat(def->W, i, &mats[i]));
502409a357bSJakub Kruzik         size -= 1;
5039566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&def->W));
504409a357bSJakub Kruzik         def->W = mats[size];
5059566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)mats[size]));
506409a357bSJakub Kruzik         if (size > 1) {
5079566063dSJacob Faibussowitsch           PetscCall(MatCreateComposite(comm, size, mats, &nextDef));
5089566063dSJacob Faibussowitsch           PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE));
509409a357bSJakub Kruzik         } else {
510409a357bSJakub Kruzik           nextDef = mats[0];
5119566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)mats[0]));
512409a357bSJakub Kruzik         }
5139566063dSJacob Faibussowitsch         PetscCall(PetscFree(mats));
514409a357bSJakub Kruzik       } else {
5151fdb00f9SJakub Kruzik         /* TODO test merge side performance */
5169566063dSJacob Faibussowitsch         /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */
5179566063dSJacob Faibussowitsch         PetscCall(MatCompositeMerge(def->W));
518409a357bSJakub Kruzik       }
51928da0170SJakub Kruzik     } else {
5206c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
5219566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(size, &mats));
52248a46eb9SPierre Jolivet         for (i = 0; i < size; i++) PetscCall(MatCompositeGetMat(def->Wt, i, &mats[i]));
52328da0170SJakub Kruzik         size -= 1;
5249566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&def->Wt));
52528da0170SJakub Kruzik         def->Wt = mats[0];
5269566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)mats[0]));
52728da0170SJakub Kruzik         if (size > 1) {
5289566063dSJacob Faibussowitsch           PetscCall(MatCreateComposite(comm, size, &mats[1], &nextDef));
5299566063dSJacob Faibussowitsch           PetscCall(MatCompositeSetType(nextDef, MAT_COMPOSITE_MULTIPLICATIVE));
53028da0170SJakub Kruzik         } else {
53128da0170SJakub Kruzik           nextDef = mats[1];
5329566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)mats[1]));
533409a357bSJakub Kruzik         }
5349566063dSJacob Faibussowitsch         PetscCall(PetscFree(mats));
53528da0170SJakub Kruzik       } else {
5369566063dSJacob Faibussowitsch         /* PetscCall(MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT)); */
5379566063dSJacob Faibussowitsch         PetscCall(MatCompositeMerge(def->Wt));
53828da0170SJakub Kruzik       }
53928da0170SJakub Kruzik     }
54028da0170SJakub Kruzik   }
54128da0170SJakub Kruzik 
54228da0170SJakub Kruzik   if (transp) {
5439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&def->W));
5449566063dSJacob Faibussowitsch     PetscCall(MatHermitianTranspose(def->Wt, MAT_INITIAL_MATRIX, &def->W));
545409a357bSJakub Kruzik   }
546409a357bSJakub Kruzik 
547ae029463SJakub Kruzik   /* assemble WtA */
548ae029463SJakub Kruzik   if (!def->WtA) {
549ae029463SJakub Kruzik     if (def->Wt) {
5509566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
551ae029463SJakub Kruzik     } else {
552a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
5539566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(def->W, MAT_INITIAL_MATRIX, &def->Wt));
5549566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->Wt, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
555a3177931SJakub Kruzik #else
5569566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(def->W, Amat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtA));
557a3177931SJakub Kruzik #endif
558ae029463SJakub Kruzik     }
559ae029463SJakub Kruzik   }
560409a357bSJakub Kruzik   /* setup coarse problem */
561409a357bSJakub Kruzik   if (!def->WtAWinv) {
5629566063dSJacob Faibussowitsch     PetscCall(MatGetSize(def->W, NULL, &m));
563409a357bSJakub Kruzik     if (!def->WtAW) {
5649566063dSJacob Faibussowitsch       PetscCall(MatMatMult(def->WtA, def->W, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &def->WtAW));
565409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
566b94d7dedSBarry Smith       PetscCall(MatIsSPDKnown(Amat, &isset, &flgspd));
567b94d7dedSBarry Smith       if (isset) PetscCall(MatSetOption(def->WtAW, MAT_SPD, flgspd));
56876bd3646SJed Brown       if (PetscDefined(USE_DEBUG)) {
569ae029463SJakub Kruzik         /* Check columns of W are not in kernel of A */
570409a357bSJakub Kruzik         PetscReal *norms;
571b94d7dedSBarry Smith 
5729566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m, &norms));
5739566063dSJacob Faibussowitsch         PetscCall(MatGetColumnNorms(def->WtAW, NORM_INFINITY, norms));
574ad540459SPierre 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);
5759566063dSJacob Faibussowitsch         PetscCall(PetscFree(norms));
576377006c5SJakub Kruzik       }
577b94d7dedSBarry Smith     } else PetscCall(MatIsSPDKnown(def->WtAW, &isset, &flgspd));
578b94d7dedSBarry Smith 
5791fdb00f9SJakub Kruzik     /* TODO use MATINV ? */
5809566063dSJacob Faibussowitsch     PetscCall(KSPCreate(comm, &def->WtAWinv));
5813821be0aSBarry Smith     PetscCall(KSPSetNestLevel(def->WtAWinv, pc->kspnestlevel));
5829566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(def->WtAWinv, def->WtAW, def->WtAW));
5839566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(def->WtAWinv, &pcinner));
584557a2f7dSJakub Kruzik     /* Setup KSP and PC */
585557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
586557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
587557a2f7dSJakub Kruzik       /* set default KSPtype */
588557a2f7dSJakub Kruzik       if (!def->ksptype) {
589557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
590b94d7dedSBarry Smith         if (isset && flgspd) { /* SPD system */
591557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
592557a2f7dSJakub Kruzik         }
593557a2f7dSJakub Kruzik       }
5949566063dSJacob Faibussowitsch       PetscCall(KSPSetType(innerksp, def->ksptype)); /* TODO iherit from KSP + tolerances */
5959566063dSJacob Faibussowitsch       PetscCall(PCSetType(pcinner, PCDEFLATION));    /* TODO create coarse preconditinoner M_c = WtMW ? */
5969566063dSJacob Faibussowitsch       PetscCall(PCDeflationSetSpace(pcinner, nextDef, transp));
5979566063dSJacob Faibussowitsch       PetscCall(PCDeflationSetLevels_Deflation(pcinner, def->lvl + 1, def->maxlvl));
598ae029463SJakub Kruzik       /* inherit options */
599*f4f49eeaSPierre Jolivet       if (def->prefix) ((PC_Deflation *)pcinner->data)->prefix = def->prefix;
600*f4f49eeaSPierre Jolivet       ((PC_Deflation *)pcinner->data)->init          = def->init;
601*f4f49eeaSPierre Jolivet       ((PC_Deflation *)pcinner->data)->ksptype       = def->ksptype;
602*f4f49eeaSPierre Jolivet       ((PC_Deflation *)pcinner->data)->correct       = def->correct;
603*f4f49eeaSPierre Jolivet       ((PC_Deflation *)pcinner->data)->correctfact   = def->correctfact;
604*f4f49eeaSPierre Jolivet       ((PC_Deflation *)pcinner->data)->reductionfact = def->reductionfact;
6059566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&nextDef));
606557a2f7dSJakub Kruzik     } else { /* the last level */
6079566063dSJacob Faibussowitsch       PetscCall(KSPSetType(def->WtAWinv, KSPPREONLY));
6089566063dSJacob Faibussowitsch       PetscCall(PCSetType(pcinner, PCTELESCOPE));
6096c93e71cSJakub Kruzik       /* do not overwrite PCTELESCOPE */
6101baa6e33SBarry Smith       if (def->prefix) PetscCall(KSPSetOptionsPrefix(def->WtAWinv, def->prefix));
6119566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(def->WtAWinv, "deflation_tel_"));
6129566063dSJacob Faibussowitsch       PetscCall(PCSetFromOptions(pcinner));
6139566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pcinner, PCTELESCOPE, &match));
614da81f932SPierre Jolivet       PetscCheck(match, comm, PETSC_ERR_SUP, "User can not overwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead.");
615409a357bSJakub Kruzik       /* Reduction factor choice */
616409a357bSJakub Kruzik       red = def->reductionfact;
617409a357bSJakub Kruzik       if (red < 0) {
6189566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm, &commsize));
619faa75363SBarry Smith         red = PetscCeilInt(commsize, PetscCeilInt(m, commsize));
620*f4f49eeaSPierre Jolivet         PetscCall(PetscObjectTypeCompareAny((PetscObject)def->WtAW, &match, MATSEQDENSE, MATMPIDENSE, MATDENSE, ""));
621409a357bSJakub Kruzik         if (match) red = commsize;
62263a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "Auto choosing reduction factor %" PetscInt_FMT "\n", red));
623409a357bSJakub Kruzik       }
6249566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSetReductionFactor(pcinner, red));
6259566063dSJacob Faibussowitsch       PetscCall(PCSetUp(pcinner));
6269566063dSJacob Faibussowitsch       PetscCall(PCTelescopeGetKSP(pcinner, &innerksp));
627ac66f006SJakub Kruzik       if (innerksp) {
6289566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(innerksp, &pcinner));
6299566063dSJacob Faibussowitsch         PetscCall(PCSetType(pcinner, PCLU));
630481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU)
6319566063dSJacob Faibussowitsch         PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU, MAT_FACTOR_LU, &match));
6321baa6e33SBarry Smith         if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU));
633481b7641SJakub Kruzik #endif
634481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU_DIST)
6359566063dSJacob Faibussowitsch         PetscCall(MatGetFactorAvailable(def->WtAW, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &match));
6361baa6e33SBarry Smith         if (match) PetscCall(PCFactorSetMatSolverType(pcinner, MATSOLVERSUPERLU_DIST));
637481b7641SJakub Kruzik #endif
638409a357bSJakub Kruzik       }
639557a2f7dSJakub Kruzik     }
640557a2f7dSJakub Kruzik 
641557a2f7dSJakub Kruzik     if (innerksp) {
6426c93e71cSJakub Kruzik       if (def->prefix) {
6439566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(innerksp, def->prefix));
6449566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(innerksp, "deflation_"));
6456c93e71cSJakub Kruzik       } else {
6469566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(innerksp, "deflation_"));
6476c93e71cSJakub Kruzik       }
6489566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(innerksp, prefix));
6499566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(innerksp));
6509566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(innerksp));
651ac66f006SJakub Kruzik     }
652409a357bSJakub Kruzik   }
6539566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(def->WtAWinv));
6549566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(def->WtAWinv));
655409a357bSJakub Kruzik 
6561fdb00f9SJakub Kruzik   /* create preconditioner */
65722b0793eSJakub Kruzik   if (!def->pc) {
6589566063dSJacob Faibussowitsch     PetscCall(PCCreate(comm, &def->pc));
6599566063dSJacob Faibussowitsch     PetscCall(PCSetOperators(def->pc, Amat, Amat));
6609566063dSJacob Faibussowitsch     PetscCall(PCSetType(def->pc, PCNONE));
6611baa6e33SBarry Smith     if (def->prefix) PetscCall(PCSetOptionsPrefix(def->pc, def->prefix));
6629566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, "deflation_"));
6639566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, prefix));
6649566063dSJacob Faibussowitsch     PetscCall(PCAppendOptionsPrefix(def->pc, "pc_"));
665d6b1b4b3SYANG Zongze     PetscCall(PCGetDM(pc, &dm));
666d6b1b4b3SYANG Zongze     PetscCall(PCSetDM(def->pc, dm));
6679566063dSJacob Faibussowitsch     PetscCall(PCSetFromOptions(def->pc));
6689566063dSJacob Faibussowitsch     PetscCall(PCSetUp(def->pc));
66922b0793eSJakub Kruzik   }
67022b0793eSJakub Kruzik 
671f8bf229cSJakub Kruzik   /* create work vecs */
6729566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Amat, NULL, &def->work));
6739566063dSJacob Faibussowitsch   PetscCall(KSPCreateVecs(def->WtAWinv, 2, &def->workcoarse, 0, NULL));
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67537eeb815SJakub Kruzik }
676b30d4775SJakub Kruzik 
677d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Deflation(PC pc)
678d71ae5a4SJacob Faibussowitsch {
679b30d4775SJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
68037eeb815SJakub Kruzik 
68137eeb815SJakub Kruzik   PetscFunctionBegin;
6829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&def->work));
6839566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2, &def->workcoarse));
6849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->W));
6859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->Wt));
6869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtA));
6879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&def->WtAW));
6889566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&def->WtAWinv));
6899566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&def->pc));
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69137eeb815SJakub Kruzik }
69237eeb815SJakub Kruzik 
693d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Deflation(PC pc)
694d71ae5a4SJacob Faibussowitsch {
69537eeb815SJakub Kruzik   PetscFunctionBegin;
6969566063dSJacob Faibussowitsch   PetscCall(PCReset_Deflation(pc));
6972e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", NULL));
6982e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", NULL));
6992e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", NULL));
7002e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", NULL));
7012e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", NULL));
7022e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", NULL));
7032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", NULL));
7042e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", NULL));
7052e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", NULL));
7062e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", NULL));
7079566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
7083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70937eeb815SJakub Kruzik }
71037eeb815SJakub Kruzik 
711d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Deflation(PC pc, PetscViewer viewer)
712d71ae5a4SJacob Faibussowitsch {
7138513960bSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
7141ac6250aSJakub Kruzik   PetscInt      its;
7158513960bSJakub Kruzik   PetscBool     iascii;
7168513960bSJakub Kruzik 
7178513960bSJakub Kruzik   PetscFunctionBegin;
7189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7198513960bSJakub Kruzik   if (iascii) {
72048a46eb9SPierre Jolivet     if (def->correct) PetscCall(PetscViewerASCIIPrintf(viewer, "using CP correction, factor = %g+%gi\n", (double)PetscRealPart(def->correctfact), (double)PetscImaginaryPart(def->correctfact)));
72148a46eb9SPierre Jolivet     if (!def->lvl) PetscCall(PetscViewerASCIIPrintf(viewer, "deflation space type: %s\n", PCDeflationSpaceTypes[def->spacetype]));
7228513960bSJakub Kruzik 
7239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "--- Additional PC:\n"));
7249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
7259566063dSJacob Faibussowitsch     PetscCall(PCView(def->pc, viewer));
7269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
72722b0793eSJakub Kruzik 
7289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "--- Coarse problem solver:\n"));
7299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
7309566063dSJacob Faibussowitsch     PetscCall(KSPGetTotalIterations(def->WtAWinv, &its));
73163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations: %" PetscInt_FMT "\n", its));
7329566063dSJacob Faibussowitsch     PetscCall(KSPView(def->WtAWinv, viewer));
7339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
7348513960bSJakub Kruzik   }
7353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7368513960bSJakub Kruzik }
7378513960bSJakub Kruzik 
738d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Deflation(PC pc, PetscOptionItems *PetscOptionsObject)
739d71ae5a4SJacob Faibussowitsch {
740880d05ceSJakub Kruzik   PC_Deflation *def = (PC_Deflation *)pc->data;
74137eeb815SJakub Kruzik 
74237eeb815SJakub Kruzik   PetscFunctionBegin;
743d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Deflation options");
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_init_only", "Use only initialization step - Initdef", "PCDeflationSetInitOnly", def->init, &def->init, NULL));
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_levels", "Maximum of deflation levels", "PCDeflationSetLevels", def->maxlvl, &def->maxlvl, NULL));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_reduction_factor", "Reduction factor for coarse problem solution using PCTELESCOPE", "PCDeflationSetReductionFactor", def->reductionfact, &def->reductionfact, NULL));
7479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_correction", "Add coarse problem correction Q to P", "PCDeflationSetCorrectionFactor", def->correct, &def->correct, NULL));
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsScalar("-pc_deflation_correction_factor", "Set multiple of Q to use as coarse problem correction", "PCDeflationSetCorrectionFactor", def->correctfact, &def->correctfact, NULL));
7499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_deflation_compute_space", "Compute deflation space", "PCDeflationSetSpace", PCDeflationSpaceTypes, (PetscEnum)def->spacetype, (PetscEnum *)&def->spacetype, NULL));
7509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_deflation_compute_space_size", "Set size of the deflation space to compute", "PCDeflationSetSpace", def->spacesize, &def->spacesize, NULL));
7519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_deflation_space_extend", "Extend deflation space instead of truncating (wavelets)", "PCDeflation", def->extendsp, &def->extendsp, NULL));
752d0609cedSBarry Smith   PetscOptionsHeadEnd();
7533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75437eeb815SJakub Kruzik }
75537eeb815SJakub Kruzik 
75637eeb815SJakub Kruzik /*MC
757e26ad46dSJakub Kruzik    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
75837eeb815SJakub Kruzik 
759e26ad46dSJakub Kruzik    Options Database Keys:
760e26ad46dSJakub Kruzik +    -pc_deflation_init_only          <false> - if true computes only the special guess
761e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
76271f2caa7Sprj- .    -pc_deflation_reduction_factor <\-1>     - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
763e26ad46dSJakub Kruzik .    -pc_deflation_correction         <false> - if true apply coarse problem correction
764e26ad46dSJakub Kruzik .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
765e26ad46dSJakub Kruzik .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
766e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
76737eeb815SJakub Kruzik 
76837eeb815SJakub Kruzik    Notes:
7697e9012c5SJakub Kruzik     Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
770e26ad46dSJakub Kruzik     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
771e26ad46dSJakub Kruzik 
772e26ad46dSJakub Kruzik     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
773f1580f4eSBarry Smith     If `PCDeflationSetInitOnly()` or -pc_deflation_init_only is set to `PETSC_TRUE` (InitDef scheme), the application of the
774e26ad46dSJakub Kruzik     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
7751fdb00f9SJakub Kruzik     application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
776e26ad46dSJakub Kruzik     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
777e26ad46dSJakub Kruzik     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
778e26ad46dSJakub Kruzik     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
779e26ad46dSJakub Kruzik 
780e26ad46dSJakub Kruzik     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
781f1580f4eSBarry Smith     be controlled by `PCDeflationSetSpaceToCompute()` or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
782f1580f4eSBarry Smith     User can set an arbitrary deflation space matrix with `PCDeflationSetSpace()`. If the deflation matrix
783f1580f4eSBarry Smith     is a multiplicative `MATCOMPOSITE`, a multilevel deflation [3] is used. The first matrix in the composite is used as the
784f1580f4eSBarry Smith     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by `KSPFCG` (if A is `MAT_SPD`) or `KSPFGMRES` preconditioned
785f1580f4eSBarry Smith     by deflation with deflation matrix being the next matrix in the `MATCOMPOSITE`. This scheme repeats until the maximum
7861fdb00f9SJakub Kruzik     level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged
7871fdb00f9SJakub Kruzik     (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by
788f1580f4eSBarry Smith     `PCDeflationSetLevels()` or by -pc_deflation_levels.
789e26ad46dSJakub Kruzik 
790f1580f4eSBarry Smith     The coarse problem `KSP` can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1]
7916c93e71cSJakub Kruzik     from the second level onward. You can also use
792f1580f4eSBarry Smith     `PCDeflationGetCoarseKSP()` to control it from code. The bottom level KSP defaults to
793f1580f4eSBarry Smith     `KSPPREONLY` with `PCLU` direct solver (`MATSOLVERSUPERLU`/`MATSOLVERSUPERLU_DIST` if available) wrapped into `PCTELESCOPE`.
794f1580f4eSBarry Smith     For convenience, the reduction factor can be set by `PCDeflationSetReductionFactor()`
795481b7641SJakub Kruzik     or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
796e26ad46dSJakub Kruzik 
797e63c2c6dSJakub Kruzik     The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for
798f1580f4eSBarry Smith     coarse problem `KSP` apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use
799f1580f4eSBarry Smith     `PCDeflationGetPC()` to control the additional preconditioner from code. It defaults to `PCNONE`.
800e26ad46dSJakub Kruzik 
801e26ad46dSJakub Kruzik     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
802f1580f4eSBarry Smith     be set by pc_deflation_correction_factor or by `PCDeflationSetCorrectionFactor()`. The coarse problem can
803e26ad46dSJakub Kruzik     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
804e26ad46dSJakub Kruzik     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
805e26ad46dSJakub Kruzik     an isolated eigenvalue.
806e26ad46dSJakub Kruzik 
8071fdb00f9SJakub Kruzik     The options are automatically inherited from the previous deflation level.
8089f604af8SJakub Kruzik 
809f1580f4eSBarry Smith     The preconditioner supports `KSPMonitorDynamicTolerance()`. This is useful for the multilevel scheme for which we also
8101fdb00f9SJakub Kruzik     recommend limiting the number of iterations for the coarse problems.
811e26ad46dSJakub Kruzik 
812aaa8cc7dSPierre Jolivet     See section 3 of [4] for additional references and description of the algorithm when used for conjugate gradients.
813e26ad46dSJakub Kruzik     Section 4 describes some possible choices for the deflation space.
814e26ad46dSJakub Kruzik 
815e26ad46dSJakub Kruzik      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
816e26ad46dSJakub Kruzik      Academy of Sciences and VSB - TU Ostrava.
817e26ad46dSJakub Kruzik 
818e26ad46dSJakub Kruzik      Developed from PERMON code used in [4] while on a research stay with
819e26ad46dSJakub Kruzik      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
820e26ad46dSJakub Kruzik 
821e26ad46dSJakub Kruzik    References:
822606c0280SSatish Balay +  * - A. Nicolaides. "Deflation of conjugate gradients with applications to boundary value problems", SIAM J. Numer. Anal. 24.2, 1987.
823606c0280SSatish Balay .  * - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
824606c0280SSatish Balay .  * - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
825606c0280SSatish 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
826e26ad46dSJakub Kruzik 
827e26ad46dSJakub Kruzik    Level: intermediate
82837eeb815SJakub Kruzik 
829562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
830db781477SPatrick Sanan           `PCDeflationSetInitOnly()`, `PCDeflationSetLevels()`, `PCDeflationSetReductionFactor()`,
831db781477SPatrick Sanan           `PCDeflationSetCorrectionFactor()`, `PCDeflationSetSpaceToCompute()`,
832db781477SPatrick Sanan           `PCDeflationSetSpace()`, `PCDeflationSpaceType`, `PCDeflationSetProjectionNullSpaceMat()`,
833db781477SPatrick Sanan           `PCDeflationSetCoarseMat()`, `PCDeflationGetCoarseKSP()`, `PCDeflationGetPC()`
83437eeb815SJakub Kruzik M*/
83537eeb815SJakub Kruzik 
836d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
837d71ae5a4SJacob Faibussowitsch {
83837eeb815SJakub Kruzik   PC_Deflation *def;
83937eeb815SJakub Kruzik 
84037eeb815SJakub Kruzik   PetscFunctionBegin;
8414dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&def));
84237eeb815SJakub Kruzik   pc->data = (void *)def;
84337eeb815SJakub Kruzik 
844e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
845e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
846fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
847e662bc50SJakub Kruzik   def->reductionfact = -1;
848e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
849e662bc50SJakub Kruzik   def->spacesize     = 1;
850e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
8516c93e71cSJakub Kruzik   def->lvl           = 0;
8526c93e71cSJakub Kruzik   def->maxlvl        = 0;
85370f9219eSJakub Kruzik   def->W             = NULL;
85470f9219eSJakub Kruzik   def->Wt            = NULL;
85537eeb815SJakub Kruzik 
85637eeb815SJakub Kruzik   pc->ops->apply          = PCApply_Deflation;
857b27c8092SJakub Kruzik   pc->ops->presolve       = PCPreSolve_Deflation;
85837eeb815SJakub Kruzik   pc->ops->setup          = PCSetUp_Deflation;
85937eeb815SJakub Kruzik   pc->ops->reset          = PCReset_Deflation;
86037eeb815SJakub Kruzik   pc->ops->destroy        = PCDestroy_Deflation;
86137eeb815SJakub Kruzik   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
8628513960bSJakub Kruzik   pc->ops->view           = PCView_Deflation;
86337eeb815SJakub Kruzik 
8649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetInitOnly_C", PCDeflationSetInitOnly_Deflation));
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetLevels_C", PCDeflationSetLevels_Deflation));
8669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetReductionFactor_C", PCDeflationSetReductionFactor_Deflation));
8679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCorrectionFactor_C", PCDeflationSetCorrectionFactor_Deflation));
8689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpaceToCompute_C", PCDeflationSetSpaceToCompute_Deflation));
8699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetSpace_C", PCDeflationSetSpace_Deflation));
8709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetProjectionNullSpaceMat_C", PCDeflationSetProjectionNullSpaceMat_Deflation));
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationSetCoarseMat_C", PCDeflationSetCoarseMat_Deflation));
8729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetCoarseKSP_C", PCDeflationGetCoarseKSP_Deflation));
8739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCDeflationGetPC_C", PCDeflationGetPC_Deflation));
8743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87537eeb815SJakub Kruzik }
876