xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision b2f12e249665d9a78b5d6d17280df9aa60e648ef)
15170378fSJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/  /* includes for fortran wrappers */
237eeb815SJakub Kruzik 
38513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = {
48513960bSJakub Kruzik   "haar",
58513960bSJakub Kruzik   "db2",
68513960bSJakub Kruzik   "db4",
78513960bSJakub Kruzik   "db8",
88513960bSJakub Kruzik   "db16",
98513960bSJakub Kruzik   "biorth22",
108513960bSJakub Kruzik   "meyer",
118513960bSJakub Kruzik   "aggregation",
128513960bSJakub Kruzik   "user",
130a78af22SJakub Kruzik   "PCDeflationSpaceType",
140a78af22SJakub Kruzik   "PC_DEFLATION_SPACE_",
158513960bSJakub Kruzik   0
168513960bSJakub Kruzik };
178513960bSJakub Kruzik 
18a122ebaeSJakub Kruzik static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc,PetscBool flg)
19a122ebaeSJakub Kruzik {
20a122ebaeSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
21a122ebaeSJakub Kruzik 
22a122ebaeSJakub Kruzik   PetscFunctionBegin;
23a122ebaeSJakub Kruzik   def->init = flg;
24a122ebaeSJakub Kruzik   PetscFunctionReturn(0);
25a122ebaeSJakub Kruzik }
26a122ebaeSJakub Kruzik 
27a122ebaeSJakub Kruzik /*@
28a122ebaeSJakub Kruzik    PCDeflationSetInitOnly - Do only initialization step.
29e26ad46dSJakub Kruzik     Sets initial guess to the solution on the deflation space but does not apply
30e26ad46dSJakub Kruzik     the deflation preconditioner. The additional preconditioner is still applied.
31a122ebaeSJakub Kruzik 
32e26ad46dSJakub Kruzik    Logically Collective
33a122ebaeSJakub Kruzik 
34a122ebaeSJakub Kruzik    Input Parameters:
35a122ebaeSJakub Kruzik +  pc  - the preconditioner context
36a122ebaeSJakub Kruzik -  flg - default PETSC_FALSE
37a122ebaeSJakub Kruzik 
38e26ad46dSJakub Kruzik    Options Database Keys:
39e26ad46dSJakub Kruzik .    -pc_deflation_init_only <false> - if true computes only the special guess
40e26ad46dSJakub Kruzik 
41a122ebaeSJakub Kruzik    Level: intermediate
42a122ebaeSJakub Kruzik 
43a122ebaeSJakub Kruzik .seealso: PCDEFLATION
44a122ebaeSJakub Kruzik @*/
45a122ebaeSJakub Kruzik PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg)
46a122ebaeSJakub Kruzik {
47a122ebaeSJakub Kruzik   PetscErrorCode ierr;
48a122ebaeSJakub Kruzik 
49a122ebaeSJakub Kruzik   PetscFunctionBegin;
50a122ebaeSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
51a122ebaeSJakub Kruzik   PetscValidLogicalCollectiveBool(pc,flg,2);
52a122ebaeSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
53a122ebaeSJakub Kruzik   PetscFunctionReturn(0);
54a122ebaeSJakub Kruzik }
55a122ebaeSJakub Kruzik 
56a122ebaeSJakub Kruzik 
5793b79a5aSJakub Kruzik static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc,PetscInt current,PetscInt max)
5898d6e3deSJakub Kruzik {
5998d6e3deSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
6098d6e3deSJakub Kruzik 
6198d6e3deSJakub Kruzik   PetscFunctionBegin;
626c93e71cSJakub Kruzik   if (current) def->lvl = current;
636c93e71cSJakub Kruzik   def->maxlvl = max;
6498d6e3deSJakub Kruzik   PetscFunctionReturn(0);
6598d6e3deSJakub Kruzik }
6698d6e3deSJakub Kruzik 
6798d6e3deSJakub Kruzik /*@
6893b79a5aSJakub Kruzik    PCDeflationSetLevels - Set the maximum level of deflation nesting.
6998d6e3deSJakub Kruzik 
70e26ad46dSJakub Kruzik    Logically Collective
7198d6e3deSJakub Kruzik 
7298d6e3deSJakub Kruzik    Input Parameters:
7398d6e3deSJakub Kruzik +  pc  - the preconditioner context
7498d6e3deSJakub Kruzik -  max - maximum deflation level
7598d6e3deSJakub Kruzik 
76e26ad46dSJakub Kruzik    Options Database Keys:
77e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation
78e26ad46dSJakub Kruzik 
7998d6e3deSJakub Kruzik    Level: intermediate
8098d6e3deSJakub Kruzik 
81e26ad46dSJakub Kruzik .seealso: PCDeflationSetSpaceToCompute(), PCDeflationSetSpace(), PCDEFLATION
8298d6e3deSJakub Kruzik @*/
8393b79a5aSJakub Kruzik PetscErrorCode PCDeflationSetLevels(PC pc,PetscInt max)
8498d6e3deSJakub Kruzik {
8598d6e3deSJakub Kruzik   PetscErrorCode ierr;
8698d6e3deSJakub Kruzik 
8798d6e3deSJakub Kruzik   PetscFunctionBegin;
8898d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8998d6e3deSJakub Kruzik   PetscValidLogicalCollectiveInt(pc,max,2);
9093b79a5aSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLevels_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
9198d6e3deSJakub Kruzik   PetscFunctionReturn(0);
9298d6e3deSJakub Kruzik }
9398d6e3deSJakub Kruzik 
94859a873cSJakub Kruzik static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red)
95859a873cSJakub Kruzik {
96859a873cSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
97859a873cSJakub Kruzik 
98859a873cSJakub Kruzik   PetscFunctionBegin;
99859a873cSJakub Kruzik   def->reductionfact = red;
100859a873cSJakub Kruzik   PetscFunctionReturn(0);
101859a873cSJakub Kruzik }
102859a873cSJakub Kruzik 
103859a873cSJakub Kruzik /*@
104e26ad46dSJakub Kruzik    PCDeflationSetReductionFactor - Set reduction factor for the bottom PCTELESCOPE coarse problem solver.
105859a873cSJakub Kruzik 
106e26ad46dSJakub Kruzik    Logically Collective
107859a873cSJakub Kruzik 
108859a873cSJakub Kruzik    Input Parameters:
109859a873cSJakub Kruzik +  pc  - the preconditioner context
110859a873cSJakub Kruzik -  red - reduction factor (or PETSC_DETERMINE)
111859a873cSJakub Kruzik 
112e26ad46dSJakub Kruzik    Options Database Keys:
113e26ad46dSJakub Kruzik .    -pc_deflation_reduction_factor <-1> - reduction factor on bottom level coarse problem for PCTELESCOPE
114e26ad46dSJakub Kruzik 
115e26ad46dSJakub Kruzik    Notes:
116e26ad46dSJakub Kruzik      Default is computed based on the size of the coarse problem.
117e26ad46dSJakub Kruzik 
118859a873cSJakub Kruzik    Level: intermediate
119859a873cSJakub Kruzik 
120e26ad46dSJakub Kruzik .seealso: PCTELESCOPE, PCDEFLATION
121859a873cSJakub Kruzik @*/
122859a873cSJakub Kruzik PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red)
123859a873cSJakub Kruzik {
124859a873cSJakub Kruzik   PetscErrorCode ierr;
125859a873cSJakub Kruzik 
126859a873cSJakub Kruzik   PetscFunctionBegin;
127859a873cSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
128859a873cSJakub Kruzik   PetscValidLogicalCollectiveInt(pc,red,2);
129859a873cSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr);
130859a873cSJakub Kruzik   PetscFunctionReturn(0);
131859a873cSJakub Kruzik }
132859a873cSJakub Kruzik 
1338a71cb68SJakub Kruzik static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact)
1348a71cb68SJakub Kruzik {
1358a71cb68SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
1368a71cb68SJakub Kruzik 
1378a71cb68SJakub Kruzik   PetscFunctionBegin;
1388a71cb68SJakub Kruzik   /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */
1398a71cb68SJakub Kruzik   def->correct = PETSC_TRUE;
1408a71cb68SJakub Kruzik   def->correctfact = fact;
141e26ad46dSJakub Kruzik   if (def->correct == 0.0) {
1428a71cb68SJakub Kruzik     def->correct = PETSC_FALSE;
1438a71cb68SJakub Kruzik   }
1448a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1458a71cb68SJakub Kruzik }
1468a71cb68SJakub Kruzik 
1478a71cb68SJakub Kruzik /*@
1488a71cb68SJakub Kruzik    PCDeflationSetCorrectionFactor - Set coarse problem correction factor.
149e26ad46dSJakub Kruzik     The Preconditioner becomes P*M^{-1} + fact*Q.
1508a71cb68SJakub Kruzik 
151e26ad46dSJakub Kruzik    Logically Collective
1528a71cb68SJakub Kruzik 
1538a71cb68SJakub Kruzik    Input Parameters:
1548a71cb68SJakub Kruzik +  pc   - the preconditioner context
155e26ad46dSJakub Kruzik -  fact - correction factor
156e26ad46dSJakub Kruzik 
157e26ad46dSJakub Kruzik    Options Database Keys:
158e26ad46dSJakub Kruzik +    -pc_deflation_correction        <false> - if true apply coarse problem correction
159e26ad46dSJakub Kruzik -    -pc_deflation_correction_factor <1.0>   - sets coarse problem correction factor
160e26ad46dSJakub Kruzik 
161e26ad46dSJakub Kruzik    Notes:
162e26ad46dSJakub Kruzik     Any non-zero fact enables the coarse problem correction.
1638a71cb68SJakub Kruzik 
1648a71cb68SJakub Kruzik    Level: intermediate
1658a71cb68SJakub Kruzik 
1668a71cb68SJakub Kruzik .seealso: PCDEFLATION
1678a71cb68SJakub Kruzik @*/
1688a71cb68SJakub Kruzik PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact)
1698a71cb68SJakub Kruzik {
1708a71cb68SJakub Kruzik   PetscErrorCode ierr;
1718a71cb68SJakub Kruzik 
1728a71cb68SJakub Kruzik   PetscFunctionBegin;
1738a71cb68SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1748a71cb68SJakub Kruzik   PetscValidLogicalCollectiveScalar(pc,fact,2);
1758a71cb68SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr);
1768a71cb68SJakub Kruzik   PetscFunctionReturn(0);
1778a71cb68SJakub Kruzik }
1788a71cb68SJakub Kruzik 
17939417d7eSJakub Kruzik static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size)
18039417d7eSJakub Kruzik {
18139417d7eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
18239417d7eSJakub Kruzik 
18339417d7eSJakub Kruzik   PetscFunctionBegin;
18439417d7eSJakub Kruzik   if (type) def->spacetype = type;
18539417d7eSJakub Kruzik   if (size > 0) def->spacesize = size;
18639417d7eSJakub Kruzik   PetscFunctionReturn(0);
18739417d7eSJakub Kruzik }
18839417d7eSJakub Kruzik 
18939417d7eSJakub Kruzik /*@
19039417d7eSJakub Kruzik    PCDeflationSetSpaceToCompute - Set deflation space type and size to compute.
19139417d7eSJakub Kruzik 
192e26ad46dSJakub Kruzik    Logically Collective
19339417d7eSJakub Kruzik 
19439417d7eSJakub Kruzik    Input Parameters:
19539417d7eSJakub Kruzik +  pc   - the preconditioner context
19639417d7eSJakub Kruzik .  type - deflation space type to compute (or PETSC_IGNORE)
19739417d7eSJakub Kruzik -  size - size of the space to compute (or PETSC_DEFAULT)
19839417d7eSJakub Kruzik 
199e26ad46dSJakub Kruzik    Options Database Keys:
200e26ad46dSJakub Kruzik +    -pc_deflation_compute_space      <haar> - compute PCDeflationSpaceType deflation space
201e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>    - size of the deflation space
202e26ad46dSJakub Kruzik 
20339417d7eSJakub Kruzik    Notes:
20439417d7eSJakub Kruzik     For wavelet-based deflation, size represents number of levels.
205e26ad46dSJakub Kruzik 
20639417d7eSJakub Kruzik     The deflation space is computed in PCSetUP().
20739417d7eSJakub Kruzik 
20839417d7eSJakub Kruzik    Level: intermediate
20939417d7eSJakub Kruzik 
21093b79a5aSJakub Kruzik .seealso: PCDeflationSetLevels(), PCDEFLATION
21139417d7eSJakub Kruzik @*/
21239417d7eSJakub Kruzik PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size)
21339417d7eSJakub Kruzik {
21439417d7eSJakub Kruzik   PetscErrorCode ierr;
21539417d7eSJakub Kruzik 
21639417d7eSJakub Kruzik   PetscFunctionBegin;
21739417d7eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
21839417d7eSJakub Kruzik   if (type) PetscValidLogicalCollectiveEnum(pc,type,2);
21939417d7eSJakub Kruzik   if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3);
22039417d7eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr);
22139417d7eSJakub Kruzik   PetscFunctionReturn(0);
22239417d7eSJakub Kruzik }
2238513960bSJakub Kruzik 
224e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
225e662bc50SJakub Kruzik {
226e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
227e662bc50SJakub Kruzik   PetscErrorCode ierr;
228e662bc50SJakub Kruzik 
229e662bc50SJakub Kruzik   PetscFunctionBegin;
23070f9219eSJakub Kruzik   /* possibly allows W' = Wt (which is valid but not tested) */
23170f9219eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
232e662bc50SJakub Kruzik   if (transpose) {
23370f9219eSJakub Kruzik     ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
234e662bc50SJakub Kruzik     def->Wt = W;
235e662bc50SJakub Kruzik   } else {
23670f9219eSJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
237e662bc50SJakub Kruzik     def->W = W;
238e662bc50SJakub Kruzik   }
239e662bc50SJakub Kruzik   PetscFunctionReturn(0);
240e662bc50SJakub Kruzik }
241e662bc50SJakub Kruzik 
242e662bc50SJakub Kruzik /*@
2437e9012c5SJakub Kruzik    PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose).
244e662bc50SJakub Kruzik 
245e26ad46dSJakub Kruzik    Logically Collective
246e662bc50SJakub Kruzik 
247e662bc50SJakub Kruzik    Input Parameters:
248e662bc50SJakub Kruzik +  pc        - the preconditioner context
249e662bc50SJakub Kruzik .  W         - deflation matrix
250e26ad46dSJakub Kruzik -  transpose - indicates that W is an explicit transpose of the deflation matrix
251e26ad46dSJakub Kruzik 
252e26ad46dSJakub Kruzik    Notes:
253e26ad46dSJakub Kruzik     Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel
254e26ad46dSJakub Kruzik     deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and
255e26ad46dSJakub Kruzik     the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with
256e26ad46dSJakub Kruzik     W1 as the deflation matrix. This repeats until the maximum level set by
25793b79a5aSJakub Kruzik     PCDeflationSetLevels() is reached or there are no more matrices available.
258e26ad46dSJakub Kruzik     If there are matrices left after reaching the maximum level,
259e26ad46dSJakub Kruzik     they are merged into a deflation matrix ...*W{n-1}*Wn.
260e662bc50SJakub Kruzik 
261e662bc50SJakub Kruzik    Level: intermediate
262e662bc50SJakub Kruzik 
26393b79a5aSJakub Kruzik .seealso: PCDeflationSetLevels(), PCDEFLATION
264e662bc50SJakub Kruzik @*/
265e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
266e662bc50SJakub Kruzik {
267e662bc50SJakub Kruzik   PetscErrorCode ierr;
268e662bc50SJakub Kruzik 
269e662bc50SJakub Kruzik   PetscFunctionBegin;
270e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
271e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
272e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
273e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
274e662bc50SJakub Kruzik   PetscFunctionReturn(0);
275e662bc50SJakub Kruzik }
276e662bc50SJakub Kruzik 
2773cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat)
2783cf3a049SJakub Kruzik {
2793cf3a049SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
2803cf3a049SJakub Kruzik   PetscErrorCode   ierr;
2813cf3a049SJakub Kruzik 
2823cf3a049SJakub Kruzik   PetscFunctionBegin;
28370f9219eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2843cf3a049SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
2853cf3a049SJakub Kruzik   def->WtA = mat;
2863cf3a049SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr);
2873cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2883cf3a049SJakub Kruzik }
2893cf3a049SJakub Kruzik 
2903cf3a049SJakub Kruzik /*@
291e26ad46dSJakub Kruzik    PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A).
2923cf3a049SJakub Kruzik 
293e26ad46dSJakub Kruzik    Collective
2943cf3a049SJakub Kruzik 
2953cf3a049SJakub Kruzik    Input Parameters:
2963cf3a049SJakub Kruzik +  pc  - preconditioner context
2973cf3a049SJakub Kruzik -  mat - projection null space matrix
2983cf3a049SJakub Kruzik 
2993cf3a049SJakub Kruzik    Level: developer
3003cf3a049SJakub Kruzik 
3013cf3a049SJakub Kruzik .seealso: PCDEFLATION
3023cf3a049SJakub Kruzik @*/
3033cf3a049SJakub Kruzik PetscErrorCode  PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat)
3043cf3a049SJakub Kruzik {
3053cf3a049SJakub Kruzik   PetscErrorCode ierr;
3063cf3a049SJakub Kruzik 
3073cf3a049SJakub Kruzik   PetscFunctionBegin;
3083cf3a049SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3093cf3a049SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
3103cf3a049SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
3113cf3a049SJakub Kruzik   PetscFunctionReturn(0);
3123cf3a049SJakub Kruzik }
3133cf3a049SJakub Kruzik 
314e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
315e924b002SJakub Kruzik {
316e924b002SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
317e924b002SJakub Kruzik   PetscErrorCode   ierr;
318e924b002SJakub Kruzik 
319e924b002SJakub Kruzik   PetscFunctionBegin;
32020cd032fSJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
321e924b002SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
322e924b002SJakub Kruzik   def->WtAW = mat;
323e924b002SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
324e924b002SJakub Kruzik   PetscFunctionReturn(0);
325e924b002SJakub Kruzik }
326e924b002SJakub Kruzik 
327e924b002SJakub Kruzik /*@
328e26ad46dSJakub Kruzik    PCDeflationSetCoarseMat - Set the coarse problem Mat.
329e924b002SJakub Kruzik 
330e26ad46dSJakub Kruzik    Collective
331e924b002SJakub Kruzik 
332e924b002SJakub Kruzik    Input Parameters:
333e924b002SJakub Kruzik +  pc  - preconditioner context
334e924b002SJakub Kruzik -  mat - coarse problem mat
335e924b002SJakub Kruzik 
336e924b002SJakub Kruzik    Level: developer
337e924b002SJakub Kruzik 
338e924b002SJakub Kruzik .seealso: PCDEFLATION
339e924b002SJakub Kruzik @*/
340e924b002SJakub Kruzik PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
341e924b002SJakub Kruzik {
342e924b002SJakub Kruzik   PetscErrorCode ierr;
343e924b002SJakub Kruzik 
344e924b002SJakub Kruzik   PetscFunctionBegin;
345e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
346e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
347e924b002SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
348e924b002SJakub Kruzik   PetscFunctionReturn(0);
349e924b002SJakub Kruzik }
350e924b002SJakub Kruzik 
35198d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
3525c0c31c5SJakub Kruzik {
3535c0c31c5SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
3545c0c31c5SJakub Kruzik 
3555c0c31c5SJakub Kruzik   PetscFunctionBegin;
35698d6e3deSJakub Kruzik   *ksp = def->WtAWinv;
3575c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
3585c0c31c5SJakub Kruzik }
3595c0c31c5SJakub Kruzik 
3605c0c31c5SJakub Kruzik /*@
3617e9012c5SJakub Kruzik    PCDeflationGetCoarseKSP - Returns the coarse problem KSP.
3625c0c31c5SJakub Kruzik 
36398d6e3deSJakub Kruzik    Not Collective
3645c0c31c5SJakub Kruzik 
3655c0c31c5SJakub Kruzik    Input Parameters:
36698d6e3deSJakub Kruzik .  pc - preconditioner context
3675c0c31c5SJakub Kruzik 
368e26ad46dSJakub Kruzik    Output Parameters:
36998d6e3deSJakub Kruzik .  ksp - coarse problem KSP context
3705c0c31c5SJakub Kruzik 
371e26ad46dSJakub Kruzik    Level: advanced
37298d6e3deSJakub Kruzik 
3738c4db21fSJakub Kruzik .seealso: PCDEFLATION
3745c0c31c5SJakub Kruzik @*/
37598d6e3deSJakub Kruzik PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
3765c0c31c5SJakub Kruzik {
3775c0c31c5SJakub Kruzik   PetscErrorCode ierr;
3785c0c31c5SJakub Kruzik 
3795c0c31c5SJakub Kruzik   PetscFunctionBegin;
38022b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
38198d6e3deSJakub Kruzik   PetscValidPointer(ksp,2);
38298d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
38398d6e3deSJakub Kruzik   PetscFunctionReturn(0);
38498d6e3deSJakub Kruzik }
38598d6e3deSJakub Kruzik 
386268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
387268c6673SJakub Kruzik {
388268c6673SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
389268c6673SJakub Kruzik 
390268c6673SJakub Kruzik   PetscFunctionBegin;
391268c6673SJakub Kruzik   *apc = def->pc;
392268c6673SJakub Kruzik   PetscFunctionReturn(0);
393268c6673SJakub Kruzik }
394268c6673SJakub Kruzik 
395268c6673SJakub Kruzik /*@
3967e9012c5SJakub Kruzik    PCDeflationGetPC - Returns the additional preconditioner M^{-1}.
397268c6673SJakub Kruzik 
398268c6673SJakub Kruzik    Not Collective
399268c6673SJakub Kruzik 
400268c6673SJakub Kruzik    Input Parameters:
401268c6673SJakub Kruzik .  pc  - the preconditioner context
402268c6673SJakub Kruzik 
403e26ad46dSJakub Kruzik    Output Parameters:
404268c6673SJakub Kruzik .  apc - additional preconditioner
405268c6673SJakub Kruzik 
406268c6673SJakub Kruzik    Level: advanced
407268c6673SJakub Kruzik 
4088c4db21fSJakub Kruzik .seealso: PCDEFLATION
409268c6673SJakub Kruzik @*/
410268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
411268c6673SJakub Kruzik {
412268c6673SJakub Kruzik   PetscErrorCode ierr;
413268c6673SJakub Kruzik 
414268c6673SJakub Kruzik   PetscFunctionBegin;
415268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
416268c6673SJakub Kruzik   PetscValidPointer(pc,2);
417268c6673SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
418268c6673SJakub Kruzik   PetscFunctionReturn(0);
419268c6673SJakub Kruzik }
420268c6673SJakub Kruzik 
42137eeb815SJakub Kruzik /*
422b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
423b27c8092SJakub Kruzik */
424b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
425b27c8092SJakub Kruzik {
426b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
427b27c8092SJakub Kruzik   Mat              A;
428b27c8092SJakub Kruzik   Vec              r,w1,w2;
429b27c8092SJakub Kruzik   PetscBool        nonzero;
430b27c8092SJakub Kruzik   PetscErrorCode   ierr;
43137eeb815SJakub Kruzik 
432b27c8092SJakub Kruzik   PetscFunctionBegin;
433b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
434b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
435b27c8092SJakub Kruzik   r  = def->work;
436b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
43737eeb815SJakub Kruzik 
438b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
439b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
440b27c8092SJakub Kruzik   if (nonzero) {
441b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);                          /*    r  <- b - Ax              */
442b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
443b27c8092SJakub Kruzik   } else {
444b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                            /*    r  <- b (x is 0)          */
445b27c8092SJakub Kruzik   }
446b27c8092SJakub Kruzik 
447a3177931SJakub Kruzik   if (def->Wt) {
448a3177931SJakub Kruzik     ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr);                   /*    w1 <- W'*r                */
449a3177931SJakub Kruzik   } else {
450a3177931SJakub Kruzik     ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
451a3177931SJakub Kruzik   }
452b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);              /*    w2 <- (W'*A*W)^{-1}*w1    */
453b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);                      /*    r  <- W*w2                */
454b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
455b27c8092SJakub Kruzik   PetscFunctionReturn(0);
456b27c8092SJakub Kruzik }
45737eeb815SJakub Kruzik 
458f8bf229cSJakub Kruzik /*
459f8bf229cSJakub Kruzik   if (def->correct) {
460ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
461f8bf229cSJakub Kruzik   } else {
462ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
463f8bf229cSJakub Kruzik   }
46437eeb815SJakub Kruzik */
465f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
466f8bf229cSJakub Kruzik {
467f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
468f8bf229cSJakub Kruzik   Mat              A;
469f8bf229cSJakub Kruzik   Vec              u,w1,w2;
470f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
471f8bf229cSJakub Kruzik 
472f8bf229cSJakub Kruzik   PetscFunctionBegin;
473f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
474f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
475f8bf229cSJakub Kruzik   u  = def->work;
476f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
477f8bf229cSJakub Kruzik 
478ae029463SJakub Kruzik   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                          /*    z <- M^{-1}*r             */
47922b0793eSJakub Kruzik   if (!def->init) {
480ae029463SJakub Kruzik     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                      /*    w1 <- W'*A*z              */
481f8bf229cSJakub Kruzik     if (def->correct) {
482ae029463SJakub Kruzik       if (def->Wt) {
483ae029463SJakub Kruzik         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);                   /*    w2 <- W'*r                */
484ae029463SJakub Kruzik       } else {
485a3177931SJakub Kruzik         ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr);  /*    w2 <- W'*r                */
486f8bf229cSJakub Kruzik       }
487ae029463SJakub Kruzik       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);      /*    w1 <- w1 - l*w2           */
488f8bf229cSJakub Kruzik     }
489f8bf229cSJakub Kruzik     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);                /*    w2 <- (W'*A*W)^{-1}*w1    */
490f8bf229cSJakub Kruzik     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                        /*    u  <- W*w2                */
49122b0793eSJakub Kruzik     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                           /*    z  <- z - u               */
49222b0793eSJakub Kruzik   }
493f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
494f8bf229cSJakub Kruzik }
495f8bf229cSJakub Kruzik 
49637eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
49737eeb815SJakub Kruzik {
498409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
499409a357bSJakub Kruzik   KSP              innerksp;
500409a357bSJakub Kruzik   PC               pcinner;
501409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
502409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
503409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
504409a357bSJakub Kruzik   MatCompositeType ctype;
505409a357bSJakub Kruzik   MPI_Comm         comm;
5066c93e71cSJakub Kruzik   char             prefix[128]="";
50737eeb815SJakub Kruzik   PetscErrorCode   ierr;
50837eeb815SJakub Kruzik 
50937eeb815SJakub Kruzik   PetscFunctionBegin;
510409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
511409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
512f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
5136c93e71cSJakub Kruzik   if (!def->lvl && !def->prefix) {
5146c93e71cSJakub Kruzik     ierr = PCGetOptionsPrefix(pc,&def->prefix);CHKERRQ(ierr);
5156c93e71cSJakub Kruzik   }
5166c93e71cSJakub Kruzik   if (def->lvl) {
517f44b4ef6SJakub Kruzik     ierr = PetscSNPrintf(prefix,sizeof(prefix),"%d_",(int)def->lvl);CHKERRQ(ierr);
5186c93e71cSJakub Kruzik   }
519f8bf229cSJakub Kruzik 
520f8bf229cSJakub Kruzik   /* compute a deflation space */
521409a357bSJakub Kruzik   if (def->W || def->Wt) {
522409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
523409a357bSJakub Kruzik   } else {
524e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
52537eeb815SJakub Kruzik   }
52637eeb815SJakub Kruzik 
527409a357bSJakub Kruzik   /* nested deflation */
528409a357bSJakub Kruzik   if (def->W) {
529409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
530409a357bSJakub Kruzik     if (match) {
531409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
532409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
53337eeb815SJakub Kruzik     }
534409a357bSJakub Kruzik   } else {
535a3177931SJakub Kruzik     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
536409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
537409a357bSJakub Kruzik     if (match) {
538409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
539409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
540409a357bSJakub Kruzik     }
541409a357bSJakub Kruzik     transp = PETSC_TRUE;
542409a357bSJakub Kruzik   }
543409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
544409a357bSJakub Kruzik     if (!transp) {
5456c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
54628da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
547409a357bSJakub Kruzik         for (i=0; i<size; i++) {
548409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
549409a357bSJakub Kruzik         }
550409a357bSJakub Kruzik         size -= 1;
551409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
552409a357bSJakub Kruzik         def->W = mats[size];
553409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
554409a357bSJakub Kruzik         if (size > 1) {
555409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
556409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
557409a357bSJakub Kruzik         } else {
558409a357bSJakub Kruzik           nextDef = mats[0];
559409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
560409a357bSJakub Kruzik         }
56128da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
562409a357bSJakub Kruzik       } else {
5631fdb00f9SJakub Kruzik         /* TODO test merge side performance */
564409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
565409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
566409a357bSJakub Kruzik       }
56728da0170SJakub Kruzik     } else {
5686c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
56928da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
57028da0170SJakub Kruzik         for (i=0; i<size; i++) {
57128da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
57228da0170SJakub Kruzik         }
57328da0170SJakub Kruzik         size -= 1;
57428da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
57528da0170SJakub Kruzik         def->Wt = mats[0];
57628da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
57728da0170SJakub Kruzik         if (size > 1) {
57828da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
57928da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
58028da0170SJakub Kruzik         } else {
58128da0170SJakub Kruzik           nextDef = mats[1];
58228da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
583409a357bSJakub Kruzik         }
584409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
58528da0170SJakub Kruzik       } else {
58628da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
58728da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
58828da0170SJakub Kruzik       }
58928da0170SJakub Kruzik     }
59028da0170SJakub Kruzik   }
59128da0170SJakub Kruzik 
59228da0170SJakub Kruzik   if (transp) {
59328da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
594a3177931SJakub Kruzik     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
595409a357bSJakub Kruzik   }
596409a357bSJakub Kruzik 
597ae029463SJakub Kruzik   /* assemble WtA */
598ae029463SJakub Kruzik   if (!def->WtA) {
599ae029463SJakub Kruzik     if (def->Wt) {
600ae029463SJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
601ae029463SJakub Kruzik     } else {
602a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
603a3177931SJakub Kruzik       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
6049e56ec8aSJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
605a3177931SJakub Kruzik #else
606ae029463SJakub Kruzik       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
607a3177931SJakub Kruzik #endif
608ae029463SJakub Kruzik     }
609ae029463SJakub Kruzik   }
610409a357bSJakub Kruzik   /* setup coarse problem */
611409a357bSJakub Kruzik   if (!def->WtAWinv) {
61228da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
613409a357bSJakub Kruzik     if (!def->WtAW) {
614ae029463SJakub Kruzik       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
615409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
616409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
617409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
618409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
619ae029463SJakub Kruzik       /* Check columns of W are not in kernel of A */
620409a357bSJakub Kruzik       PetscReal *norms;
621409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
622409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
623409a357bSJakub Kruzik       for (i=0; i<m; i++) {
624409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
625409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
626409a357bSJakub Kruzik         }
627409a357bSJakub Kruzik       }
628409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
629409a357bSJakub Kruzik #endif
630409a357bSJakub Kruzik     } else {
631409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
632409a357bSJakub Kruzik     }
6331fdb00f9SJakub Kruzik     /* TODO use MATINV ? */
634409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
635409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
636409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
637557a2f7dSJakub Kruzik     /* Setup KSP and PC */
638557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
639557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
640557a2f7dSJakub Kruzik       /* set default KSPtype */
641557a2f7dSJakub Kruzik       if (!def->ksptype) {
642557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
643557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
644557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
645557a2f7dSJakub Kruzik         }
646557a2f7dSJakub Kruzik       }
647ae029463SJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
648557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
649557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
65093b79a5aSJakub Kruzik       ierr = PCDeflationSetLevels_Deflation(pcinner,def->lvl+1,def->maxlvl);CHKERRQ(ierr);
651ae029463SJakub Kruzik       /* inherit options */
6526c93e71cSJakub Kruzik       if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix;
6539f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->init          = def->init;
654557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype       = def->ksptype;
655557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct       = def->correct;
6569f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correctfact   = def->correctfact;
6579f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact;
658557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
659557a2f7dSJakub Kruzik     } else { /* the last level */
660557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
661409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
6626c93e71cSJakub Kruzik       /* do not overwrite PCTELESCOPE */
6636c93e71cSJakub Kruzik       if (def->prefix) {
6646c93e71cSJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,def->prefix);CHKERRQ(ierr);
665ac66f006SJakub Kruzik       }
666e63c2c6dSJakub Kruzik       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"deflation_tel_");CHKERRQ(ierr);
667ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
668*b2f12e24SJakub Kruzik       ierr = PetscObjectTypeCompare((PetscObject)pcinner,PCTELESCOPE,&match);CHKERRQ(ierr);
669*b2f12e24SJakub Kruzik       if (!match) SETERRQ(comm,PETSC_ERR_SUP,"User can not owerwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead.");
670409a357bSJakub Kruzik       /* Reduction factor choice */
671409a357bSJakub Kruzik       red = def->reductionfact;
672409a357bSJakub Kruzik       if (red < 0) {
673409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
674409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
675409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
676409a357bSJakub Kruzik         if (match) red = commsize;
6776c93e71cSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr);
678409a357bSJakub Kruzik       }
679409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
680ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
681409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
682ac66f006SJakub Kruzik       if (innerksp) {
683409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
684ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
685481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU)
686481b7641SJakub Kruzik         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
687481b7641SJakub Kruzik         if (match) {
688ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
689481b7641SJakub Kruzik         }
690481b7641SJakub Kruzik #endif
691481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU_DIST)
692481b7641SJakub Kruzik         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
693481b7641SJakub Kruzik         if (match) {
694ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
695409a357bSJakub Kruzik         }
696481b7641SJakub Kruzik #endif
697409a357bSJakub Kruzik       }
698557a2f7dSJakub Kruzik     }
699557a2f7dSJakub Kruzik 
700557a2f7dSJakub Kruzik     if (innerksp) {
7016c93e71cSJakub Kruzik       if (def->prefix) {
7026c93e71cSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,def->prefix);CHKERRQ(ierr);
703e63c2c6dSJakub Kruzik         ierr = KSPAppendOptionsPrefix(innerksp,"deflation_");CHKERRQ(ierr);
7046c93e71cSJakub Kruzik       } else {
705e63c2c6dSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,"deflation_");CHKERRQ(ierr);
7066c93e71cSJakub Kruzik       }
7076c93e71cSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
708557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
709557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
710ac66f006SJakub Kruzik     }
711409a357bSJakub Kruzik   }
712557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
713557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
714409a357bSJakub Kruzik 
7151fdb00f9SJakub Kruzik   /* create preconditioner */
71622b0793eSJakub Kruzik   if (!def->pc) {
71722b0793eSJakub Kruzik     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
71822b0793eSJakub Kruzik     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
71922b0793eSJakub Kruzik     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
7206c93e71cSJakub Kruzik     if (def->prefix) {
7216c93e71cSJakub Kruzik       ierr = PCSetOptionsPrefix(def->pc,def->prefix);CHKERRQ(ierr);
72222b0793eSJakub Kruzik     }
723e63c2c6dSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"deflation_");CHKERRQ(ierr);
7246c93e71cSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
7256c93e71cSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"pc_");CHKERRQ(ierr);
72622b0793eSJakub Kruzik     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
72722b0793eSJakub Kruzik     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
72822b0793eSJakub Kruzik   }
72922b0793eSJakub Kruzik 
730f8bf229cSJakub Kruzik   /* create work vecs */
731f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
732f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
73337eeb815SJakub Kruzik   PetscFunctionReturn(0);
73437eeb815SJakub Kruzik }
735b30d4775SJakub Kruzik 
73637eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
73737eeb815SJakub Kruzik {
738b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
73937eeb815SJakub Kruzik   PetscErrorCode    ierr;
74037eeb815SJakub Kruzik 
74137eeb815SJakub Kruzik   PetscFunctionBegin;
742b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
743b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
744b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
745b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
746ae029463SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
747b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
748b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
74922b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
75037eeb815SJakub Kruzik   PetscFunctionReturn(0);
75137eeb815SJakub Kruzik }
75237eeb815SJakub Kruzik 
75337eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
75437eeb815SJakub Kruzik {
75537eeb815SJakub Kruzik   PetscErrorCode ierr;
75637eeb815SJakub Kruzik 
75737eeb815SJakub Kruzik   PetscFunctionBegin;
75837eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
75937eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
76037eeb815SJakub Kruzik   PetscFunctionReturn(0);
76137eeb815SJakub Kruzik }
76237eeb815SJakub Kruzik 
7638513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
7648513960bSJakub Kruzik {
7658513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
7661ac6250aSJakub Kruzik   PetscInt          its;
7678513960bSJakub Kruzik   PetscBool         iascii;
7688513960bSJakub Kruzik   PetscErrorCode    ierr;
7698513960bSJakub Kruzik 
7708513960bSJakub Kruzik   PetscFunctionBegin;
7718513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7728513960bSJakub Kruzik   if (iascii) {
7738513960bSJakub Kruzik     if (def->correct) {
7741ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n",
7751ac6250aSJakub Kruzik                                     (double)PetscRealPart(def->correctfact),
7761ac6250aSJakub Kruzik                                     (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr);
7778513960bSJakub Kruzik     }
7786c93e71cSJakub Kruzik     if (!def->lvl) {
7791ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
7808513960bSJakub Kruzik     }
7818513960bSJakub Kruzik 
7821ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr);
78322b0793eSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
78422b0793eSJakub Kruzik     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
78522b0793eSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
78622b0793eSJakub Kruzik 
7871ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr);
7888513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7891ac6250aSJakub Kruzik     ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr);
7901ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr);
7918513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
7928513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7938513960bSJakub Kruzik   }
7948513960bSJakub Kruzik   PetscFunctionReturn(0);
7958513960bSJakub Kruzik }
7968513960bSJakub Kruzik 
79737eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
79837eeb815SJakub Kruzik {
799880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
80037eeb815SJakub Kruzik   PetscErrorCode    ierr;
80137eeb815SJakub Kruzik 
80237eeb815SJakub Kruzik   PetscFunctionBegin;
80337eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
804a122ebaeSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr);
80593b79a5aSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_levels","Maximum of deflation levels","PCDeflationSetLevels",def->maxlvl,&def->maxlvl,NULL);CHKERRQ(ierr);
806859a873cSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
8078a71cb68SJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr);
8088a71cb68SJakub Kruzik   ierr = PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
809880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
810880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
811880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
81237eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
81337eeb815SJakub Kruzik   PetscFunctionReturn(0);
81437eeb815SJakub Kruzik }
81537eeb815SJakub Kruzik 
81637eeb815SJakub Kruzik /*MC
817e26ad46dSJakub Kruzik    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
81837eeb815SJakub Kruzik 
819e26ad46dSJakub Kruzik    Options Database Keys:
820e26ad46dSJakub Kruzik +    -pc_deflation_init_only          <false> - if true computes only the special guess
821e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
822e26ad46dSJakub Kruzik .    -pc_deflation_reduction_factor <-1>      - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
823e26ad46dSJakub Kruzik .    -pc_deflation_correction         <false> - if true apply coarse problem correction
824e26ad46dSJakub Kruzik .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
825e26ad46dSJakub Kruzik .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
826e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
82737eeb815SJakub Kruzik 
82837eeb815SJakub Kruzik    Notes:
8297e9012c5SJakub Kruzik     Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
830e26ad46dSJakub Kruzik     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
831e26ad46dSJakub Kruzik 
832e26ad46dSJakub Kruzik     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
833e26ad46dSJakub Kruzik     If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the
834e26ad46dSJakub Kruzik     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
8351fdb00f9SJakub Kruzik     application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
836e26ad46dSJakub Kruzik     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
837e26ad46dSJakub Kruzik     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
838e26ad46dSJakub Kruzik     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
839e26ad46dSJakub Kruzik 
840e26ad46dSJakub Kruzik     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
841e26ad46dSJakub Kruzik     be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
842e26ad46dSJakub Kruzik     User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix
843e26ad46dSJakub Kruzik     is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the
8441fdb00f9SJakub Kruzik     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned
845e26ad46dSJakub Kruzik     by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum
8461fdb00f9SJakub Kruzik     level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged
8471fdb00f9SJakub Kruzik     (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by
84893b79a5aSJakub Kruzik     PCDeflationSetLevels() or by -pc_deflation_levels.
849e26ad46dSJakub Kruzik 
850e63c2c6dSJakub Kruzik     The coarse problem KSP can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1]
8516c93e71cSJakub Kruzik     from the second level onward. You can also use
8528c4db21fSJakub Kruzik     PCDeflationGetCoarseKSP() to control it from code. The bottom level KSP defaults to
8531fdb00f9SJakub Kruzik     KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE.
854481b7641SJakub Kruzik     For convenience, the reduction factor can be set by PCDeflationSetReductionFactor()
855481b7641SJakub Kruzik     or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
856e26ad46dSJakub Kruzik 
857e63c2c6dSJakub Kruzik     The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for
858e63c2c6dSJakub Kruzik     coarse problem KSP apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use
8598c4db21fSJakub Kruzik     PCDeflationGetPC() to control the additional preconditioner from code. It defaults to PCNONE.
860e26ad46dSJakub Kruzik 
861e26ad46dSJakub Kruzik     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
862e26ad46dSJakub Kruzik     be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can
863e26ad46dSJakub Kruzik     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
864e26ad46dSJakub Kruzik     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
865e26ad46dSJakub Kruzik     an isolated eigenvalue.
866e26ad46dSJakub Kruzik 
8671fdb00f9SJakub Kruzik     The options are automatically inherited from the previous deflation level.
8689f604af8SJakub Kruzik 
869e26ad46dSJakub Kruzik     The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also
8701fdb00f9SJakub Kruzik     recommend limiting the number of iterations for the coarse problems.
871e26ad46dSJakub Kruzik 
872e26ad46dSJakub Kruzik     See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
873e26ad46dSJakub Kruzik     Section 4 describes some possible choices for the deflation space.
874e26ad46dSJakub Kruzik 
875e26ad46dSJakub Kruzik    Developer Notes:
876e26ad46dSJakub Kruzik      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
877e26ad46dSJakub Kruzik      Academy of Sciences and VSB - TU Ostrava.
878e26ad46dSJakub Kruzik 
879e26ad46dSJakub Kruzik      Developed from PERMON code used in [4] while on a research stay with
880e26ad46dSJakub Kruzik      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
881e26ad46dSJakub Kruzik 
882e26ad46dSJakub Kruzik    References:
883e26ad46dSJakub Kruzik +    [1] - A. Nicolaides. “Deflation of conjugate gradients with applications to boundary valueproblems”, SIAM J. Numer. Anal. 24.2, 1987.
884e26ad46dSJakub Kruzik .    [2] - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
885e26ad46dSJakub Kruzik .    [3] - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
8861fdb00f9SJakub Kruzik -    [4] - 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
887e26ad46dSJakub Kruzik 
888e26ad46dSJakub Kruzik    Level: intermediate
88937eeb815SJakub Kruzik 
89037eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
89193b79a5aSJakub Kruzik            PCDeflationSetInitOnly(), PCDeflationSetLevels(), PCDeflationSetReductionFactor(),
8921fdb00f9SJakub Kruzik            PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompute(),
893e26ad46dSJakub Kruzik            PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(),
8948c4db21fSJakub Kruzik            PCDeflationSetCoarseMat(), PCDeflationGetCoarseKSP(), PCDeflationGetPC()
89537eeb815SJakub Kruzik M*/
89637eeb815SJakub Kruzik 
89737eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
89837eeb815SJakub Kruzik {
89937eeb815SJakub Kruzik   PC_Deflation   *def;
90037eeb815SJakub Kruzik   PetscErrorCode ierr;
90137eeb815SJakub Kruzik 
90237eeb815SJakub Kruzik   PetscFunctionBegin;
90337eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
90437eeb815SJakub Kruzik   pc->data = (void*)def;
90537eeb815SJakub Kruzik 
906e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
907e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
908fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
909e662bc50SJakub Kruzik   def->reductionfact = -1;
910e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
911e662bc50SJakub Kruzik   def->spacesize     = 1;
912e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
9136c93e71cSJakub Kruzik   def->lvl           = 0;
9146c93e71cSJakub Kruzik   def->maxlvl        = 0;
91570f9219eSJakub Kruzik   def->W             = NULL;
91670f9219eSJakub Kruzik   def->Wt            = NULL;
91737eeb815SJakub Kruzik 
91837eeb815SJakub Kruzik   pc->ops->apply          = PCApply_Deflation;
919b27c8092SJakub Kruzik   pc->ops->presolve       = PCPreSolve_Deflation;
92037eeb815SJakub Kruzik   pc->ops->setup          = PCSetUp_Deflation;
92137eeb815SJakub Kruzik   pc->ops->reset          = PCReset_Deflation;
92237eeb815SJakub Kruzik   pc->ops->destroy        = PCDestroy_Deflation;
92337eeb815SJakub Kruzik   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
9248513960bSJakub Kruzik   pc->ops->view           = PCView_Deflation;
92537eeb815SJakub Kruzik 
926a122ebaeSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr);
92793b79a5aSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLevels_C",PCDeflationSetLevels_Deflation);CHKERRQ(ierr);
928859a873cSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr);
9298a71cb68SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr);
93039417d7eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr);
931e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
9323cf3a049SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr);
933e924b002SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
9344a99276eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
93598d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
93637eeb815SJakub Kruzik   PetscFunctionReturn(0);
93737eeb815SJakub Kruzik }
93837eeb815SJakub Kruzik 
939