xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 1fdb00f988be8a1800bb0e9f9be4e7bc64bf28d3)
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 
5798d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
5898d6e3deSJakub Kruzik {
5998d6e3deSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
6098d6e3deSJakub Kruzik 
6198d6e3deSJakub Kruzik   PetscFunctionBegin;
626c93e71cSJakub Kruzik   if (current) def->lvl = current;
636c93e71cSJakub Kruzik   def->maxlvl = max;
6498d6e3deSJakub Kruzik   PetscFunctionReturn(0);
6598d6e3deSJakub Kruzik }
6698d6e3deSJakub Kruzik 
6798d6e3deSJakub Kruzik /*@
68e26ad46dSJakub Kruzik    PCDeflationSetMaxLvl - 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 @*/
8398d6e3deSJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(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);
9098d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_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 
210e26ad46dSJakub Kruzik .seealso: PCDeflationSetMaxLvl(), 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;
230e662bc50SJakub Kruzik   if (transpose) {
231e662bc50SJakub Kruzik     def->Wt = W;
232e662bc50SJakub Kruzik     def->W = NULL;
233e662bc50SJakub Kruzik   } else {
234e662bc50SJakub Kruzik     def->W = W;
235e662bc50SJakub Kruzik   }
236e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
237e662bc50SJakub Kruzik   PetscFunctionReturn(0);
238e662bc50SJakub Kruzik }
239e662bc50SJakub Kruzik 
240e662bc50SJakub Kruzik /*@
241e26ad46dSJakub Kruzik    PCDeflationSetSpace - Set the deflation space matrix (or its (hermitian) transpose).
242e662bc50SJakub Kruzik 
243e26ad46dSJakub Kruzik    Logically Collective
244e662bc50SJakub Kruzik 
245e662bc50SJakub Kruzik    Input Parameters:
246e662bc50SJakub Kruzik +  pc        - the preconditioner context
247e662bc50SJakub Kruzik .  W         - deflation matrix
248e26ad46dSJakub Kruzik -  transpose - indicates that W is an explicit transpose of the deflation matrix
249e26ad46dSJakub Kruzik 
250e26ad46dSJakub Kruzik    Notes:
251e26ad46dSJakub Kruzik     Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel
252e26ad46dSJakub Kruzik     deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and
253e26ad46dSJakub Kruzik     the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with
254e26ad46dSJakub Kruzik     W1 as the deflation matrix. This repeats until the maximum level set by
255e26ad46dSJakub Kruzik     PCDeflationSetMaxLvl is reached or there are no more matrices available.
256e26ad46dSJakub Kruzik     If there are matrices left after reaching the maximum level,
257e26ad46dSJakub Kruzik     they are merged into a deflation matrix ...*W{n-1}*Wn.
258e662bc50SJakub Kruzik 
259e662bc50SJakub Kruzik    Level: intermediate
260e662bc50SJakub Kruzik 
261e26ad46dSJakub Kruzik .seealso: PCDeflationSetMaxLvl(), PCDEFLATION
262e662bc50SJakub Kruzik @*/
263e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
264e662bc50SJakub Kruzik {
265e662bc50SJakub Kruzik   PetscErrorCode ierr;
266e662bc50SJakub Kruzik 
267e662bc50SJakub Kruzik   PetscFunctionBegin;
268e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
269e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
270e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
271e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
272e662bc50SJakub Kruzik   PetscFunctionReturn(0);
273e662bc50SJakub Kruzik }
274e662bc50SJakub Kruzik 
2753cf3a049SJakub Kruzik static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat)
2763cf3a049SJakub Kruzik {
2773cf3a049SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
2783cf3a049SJakub Kruzik   PetscErrorCode   ierr;
2793cf3a049SJakub Kruzik 
2803cf3a049SJakub Kruzik   PetscFunctionBegin;
2813cf3a049SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
2823cf3a049SJakub Kruzik   def->WtA = mat;
2833cf3a049SJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2843cf3a049SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr);
2853cf3a049SJakub Kruzik   PetscFunctionReturn(0);
2863cf3a049SJakub Kruzik }
2873cf3a049SJakub Kruzik 
2883cf3a049SJakub Kruzik /*@
289e26ad46dSJakub Kruzik    PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A).
2903cf3a049SJakub Kruzik 
291e26ad46dSJakub Kruzik    Collective
2923cf3a049SJakub Kruzik 
2933cf3a049SJakub Kruzik    Input Parameters:
2943cf3a049SJakub Kruzik +  pc  - preconditioner context
2953cf3a049SJakub Kruzik -  mat - projection null space matrix
2963cf3a049SJakub Kruzik 
2973cf3a049SJakub Kruzik    Level: developer
2983cf3a049SJakub Kruzik 
2993cf3a049SJakub Kruzik .seealso: PCDEFLATION
3003cf3a049SJakub Kruzik @*/
3013cf3a049SJakub Kruzik PetscErrorCode  PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat)
3023cf3a049SJakub Kruzik {
3033cf3a049SJakub Kruzik   PetscErrorCode ierr;
3043cf3a049SJakub Kruzik 
3053cf3a049SJakub Kruzik   PetscFunctionBegin;
3063cf3a049SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3073cf3a049SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
3083cf3a049SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
3093cf3a049SJakub Kruzik   PetscFunctionReturn(0);
3103cf3a049SJakub Kruzik }
3113cf3a049SJakub Kruzik 
312e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
313e924b002SJakub Kruzik {
314e924b002SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
315e924b002SJakub Kruzik   PetscErrorCode   ierr;
316e924b002SJakub Kruzik 
317e924b002SJakub Kruzik   PetscFunctionBegin;
318e924b002SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
319e924b002SJakub Kruzik   def->WtAW = mat;
320e924b002SJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
321e924b002SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
322e924b002SJakub Kruzik   PetscFunctionReturn(0);
323e924b002SJakub Kruzik }
324e924b002SJakub Kruzik 
325e924b002SJakub Kruzik /*@
326e26ad46dSJakub Kruzik    PCDeflationSetCoarseMat - Set the coarse problem Mat.
327e924b002SJakub Kruzik 
328e26ad46dSJakub Kruzik    Collective
329e924b002SJakub Kruzik 
330e924b002SJakub Kruzik    Input Parameters:
331e924b002SJakub Kruzik +  pc  - preconditioner context
332e924b002SJakub Kruzik -  mat - coarse problem mat
333e924b002SJakub Kruzik 
334e924b002SJakub Kruzik    Level: developer
335e924b002SJakub Kruzik 
336e924b002SJakub Kruzik .seealso: PCDEFLATION
337e924b002SJakub Kruzik @*/
338e924b002SJakub Kruzik PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
339e924b002SJakub Kruzik {
340e924b002SJakub Kruzik   PetscErrorCode ierr;
341e924b002SJakub Kruzik 
342e924b002SJakub Kruzik   PetscFunctionBegin;
343e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
344e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
345e924b002SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
346e924b002SJakub Kruzik   PetscFunctionReturn(0);
347e924b002SJakub Kruzik }
348e924b002SJakub Kruzik 
34998d6e3deSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
3505c0c31c5SJakub Kruzik {
3515c0c31c5SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
3525c0c31c5SJakub Kruzik 
3535c0c31c5SJakub Kruzik   PetscFunctionBegin;
35498d6e3deSJakub Kruzik   *ksp = def->WtAWinv;
3555c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
3565c0c31c5SJakub Kruzik }
3575c0c31c5SJakub Kruzik 
3585c0c31c5SJakub Kruzik /*@
35998d6e3deSJakub Kruzik    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
3605c0c31c5SJakub Kruzik 
36198d6e3deSJakub Kruzik    Not Collective
3625c0c31c5SJakub Kruzik 
3635c0c31c5SJakub Kruzik    Input Parameters:
36498d6e3deSJakub Kruzik .  pc - preconditioner context
3655c0c31c5SJakub Kruzik 
366e26ad46dSJakub Kruzik    Output Parameters:
36798d6e3deSJakub Kruzik .  ksp - coarse problem KSP context
3685c0c31c5SJakub Kruzik 
369e26ad46dSJakub Kruzik    Level: advanced
37098d6e3deSJakub Kruzik 
37198d6e3deSJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION
3725c0c31c5SJakub Kruzik @*/
37398d6e3deSJakub Kruzik PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
3745c0c31c5SJakub Kruzik {
3755c0c31c5SJakub Kruzik   PetscErrorCode ierr;
3765c0c31c5SJakub Kruzik 
3775c0c31c5SJakub Kruzik   PetscFunctionBegin;
37822b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
37998d6e3deSJakub Kruzik   PetscValidPointer(ksp,2);
38098d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
38198d6e3deSJakub Kruzik   PetscFunctionReturn(0);
38298d6e3deSJakub Kruzik }
38398d6e3deSJakub Kruzik 
38498d6e3deSJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp)
38598d6e3deSJakub Kruzik {
38698d6e3deSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
38798d6e3deSJakub Kruzik   PetscErrorCode   ierr;
38898d6e3deSJakub Kruzik 
38998d6e3deSJakub Kruzik   PetscFunctionBegin;
39098d6e3deSJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
39198d6e3deSJakub Kruzik   def->WtAWinv = ksp;
39298d6e3deSJakub Kruzik   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
39398d6e3deSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr);
39498d6e3deSJakub Kruzik   PetscFunctionReturn(0);
39598d6e3deSJakub Kruzik }
39698d6e3deSJakub Kruzik 
39798d6e3deSJakub Kruzik /*@
398e26ad46dSJakub Kruzik    PCDeflationSetCoarseKSP - Set the coarse problem KSP.
39998d6e3deSJakub Kruzik 
400e26ad46dSJakub Kruzik    Collective
40198d6e3deSJakub Kruzik 
40298d6e3deSJakub Kruzik    Input Parameters:
40398d6e3deSJakub Kruzik +  pc  - preconditioner context
40498d6e3deSJakub Kruzik -  ksp - coarse problem KSP context
40598d6e3deSJakub Kruzik 
40698d6e3deSJakub Kruzik    Level: developer
40798d6e3deSJakub Kruzik 
40898d6e3deSJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION
40998d6e3deSJakub Kruzik @*/
41098d6e3deSJakub Kruzik PetscErrorCode  PCDeflationSetCoarseKSP(PC pc,KSP ksp)
41198d6e3deSJakub Kruzik {
41298d6e3deSJakub Kruzik   PetscErrorCode ierr;
41398d6e3deSJakub Kruzik 
41498d6e3deSJakub Kruzik   PetscFunctionBegin;
41598d6e3deSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
41698d6e3deSJakub Kruzik   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
41798d6e3deSJakub Kruzik   PetscCheckSameComm(pc,1,ksp,2);
41898d6e3deSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr);
4195c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
4205c0c31c5SJakub Kruzik }
421e662bc50SJakub Kruzik 
422268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
423268c6673SJakub Kruzik {
424268c6673SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
425268c6673SJakub Kruzik 
426268c6673SJakub Kruzik   PetscFunctionBegin;
427268c6673SJakub Kruzik   *apc = def->pc;
428268c6673SJakub Kruzik   PetscFunctionReturn(0);
429268c6673SJakub Kruzik }
430268c6673SJakub Kruzik 
431268c6673SJakub Kruzik /*@
432e26ad46dSJakub Kruzik    PCDeflationGetPC - Returns a pointer to the additional preconditioner.
433268c6673SJakub Kruzik 
434268c6673SJakub Kruzik    Not Collective
435268c6673SJakub Kruzik 
436268c6673SJakub Kruzik    Input Parameters:
437268c6673SJakub Kruzik .  pc  - the preconditioner context
438268c6673SJakub Kruzik 
439e26ad46dSJakub Kruzik    Output Parameters:
440268c6673SJakub Kruzik .  apc - additional preconditioner
441268c6673SJakub Kruzik 
442268c6673SJakub Kruzik    Level: advanced
443268c6673SJakub Kruzik 
444268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION
445268c6673SJakub Kruzik @*/
446268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
447268c6673SJakub Kruzik {
448268c6673SJakub Kruzik   PetscErrorCode ierr;
449268c6673SJakub Kruzik 
450268c6673SJakub Kruzik   PetscFunctionBegin;
451268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
452268c6673SJakub Kruzik   PetscValidPointer(pc,2);
453268c6673SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
454268c6673SJakub Kruzik   PetscFunctionReturn(0);
455268c6673SJakub Kruzik }
456268c6673SJakub Kruzik 
45722b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
45822b0793eSJakub Kruzik {
45922b0793eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
46022b0793eSJakub Kruzik   PetscErrorCode ierr;
46122b0793eSJakub Kruzik 
46222b0793eSJakub Kruzik   PetscFunctionBegin;
46322b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
46422b0793eSJakub Kruzik   def->pc = apc;
46522b0793eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
46622b0793eSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
46722b0793eSJakub Kruzik   PetscFunctionReturn(0);
46822b0793eSJakub Kruzik }
46922b0793eSJakub Kruzik 
47022b0793eSJakub Kruzik /*@
471e26ad46dSJakub Kruzik    PCDeflationSetPC - Set the additional preconditioner.
47222b0793eSJakub Kruzik 
473e26ad46dSJakub Kruzik    Collective
47422b0793eSJakub Kruzik 
47522b0793eSJakub Kruzik    Input Parameters:
47622b0793eSJakub Kruzik +  pc  - the preconditioner context
47722b0793eSJakub Kruzik -  apc - additional preconditioner
47822b0793eSJakub Kruzik 
479268c6673SJakub Kruzik    Level: developer
48022b0793eSJakub Kruzik 
481268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION
48222b0793eSJakub Kruzik @*/
48322b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
48422b0793eSJakub Kruzik {
48522b0793eSJakub Kruzik   PetscErrorCode ierr;
48622b0793eSJakub Kruzik 
48722b0793eSJakub Kruzik   PetscFunctionBegin;
48822b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
48922b0793eSJakub Kruzik   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
49022b0793eSJakub Kruzik   PetscCheckSameComm(pc,1,apc,2);
49122b0793eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
49222b0793eSJakub Kruzik   PetscFunctionReturn(0);
49322b0793eSJakub Kruzik }
49422b0793eSJakub Kruzik 
49537eeb815SJakub Kruzik /*
496b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
497b27c8092SJakub Kruzik */
498b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
499b27c8092SJakub Kruzik {
500b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
501b27c8092SJakub Kruzik   Mat              A;
502b27c8092SJakub Kruzik   Vec              r,w1,w2;
503b27c8092SJakub Kruzik   PetscBool        nonzero;
504b27c8092SJakub Kruzik   PetscErrorCode   ierr;
50537eeb815SJakub Kruzik 
506b27c8092SJakub Kruzik   PetscFunctionBegin;
507b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
508b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
509b27c8092SJakub Kruzik   r  = def->work;
510b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
51137eeb815SJakub Kruzik 
512b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
513b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
514b27c8092SJakub Kruzik   if (nonzero) {
515b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);                          /*    r  <- b - Ax              */
516b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
517b27c8092SJakub Kruzik   } else {
518b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                            /*    r  <- b (x is 0)          */
519b27c8092SJakub Kruzik   }
520b27c8092SJakub Kruzik 
521a3177931SJakub Kruzik   if (def->Wt) {
522a3177931SJakub Kruzik     ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr);                   /*    w1 <- W'*r                */
523a3177931SJakub Kruzik   } else {
524a3177931SJakub Kruzik     ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
525a3177931SJakub Kruzik   }
526b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);              /*    w2 <- (W'*A*W)^{-1}*w1    */
527b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);                      /*    r  <- W*w2                */
528b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
529b27c8092SJakub Kruzik   PetscFunctionReturn(0);
530b27c8092SJakub Kruzik }
53137eeb815SJakub Kruzik 
532f8bf229cSJakub Kruzik /*
533f8bf229cSJakub Kruzik   if (def->correct) {
534ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
535f8bf229cSJakub Kruzik   } else {
536ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
537f8bf229cSJakub Kruzik   }
53837eeb815SJakub Kruzik */
539f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
540f8bf229cSJakub Kruzik {
541f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
542f8bf229cSJakub Kruzik   Mat              A;
543f8bf229cSJakub Kruzik   Vec              u,w1,w2;
544f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
545f8bf229cSJakub Kruzik 
546f8bf229cSJakub Kruzik   PetscFunctionBegin;
547f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
548f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
549f8bf229cSJakub Kruzik   u  = def->work;
550f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
551f8bf229cSJakub Kruzik 
552ae029463SJakub Kruzik   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                          /*    z <- M^{-1}*r             */
55322b0793eSJakub Kruzik   if (!def->init) {
554ae029463SJakub Kruzik     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                      /*    w1 <- W'*A*z              */
555f8bf229cSJakub Kruzik     if (def->correct) {
556ae029463SJakub Kruzik       if (def->Wt) {
557ae029463SJakub Kruzik         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);                   /*    w2 <- W'*r                */
558ae029463SJakub Kruzik       } else {
559a3177931SJakub Kruzik         ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr);  /*    w2 <- W'*r                */
560f8bf229cSJakub Kruzik       }
561ae029463SJakub Kruzik       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);      /*    w1 <- w1 - l*w2           */
562f8bf229cSJakub Kruzik     }
563f8bf229cSJakub Kruzik     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);                /*    w2 <- (W'*A*W)^{-1}*w1    */
564f8bf229cSJakub Kruzik     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                        /*    u  <- W*w2                */
56522b0793eSJakub Kruzik     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                           /*    z  <- z - u               */
56622b0793eSJakub Kruzik   }
567f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
568f8bf229cSJakub Kruzik }
569f8bf229cSJakub Kruzik 
57037eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
57137eeb815SJakub Kruzik {
572409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
573409a357bSJakub Kruzik   KSP              innerksp;
574409a357bSJakub Kruzik   PC               pcinner;
575409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
576409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
577409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
578409a357bSJakub Kruzik   MatCompositeType ctype;
579409a357bSJakub Kruzik   MPI_Comm         comm;
5806c93e71cSJakub Kruzik   char             prefix[128]="";
58137eeb815SJakub Kruzik   PetscErrorCode   ierr;
58237eeb815SJakub Kruzik 
58337eeb815SJakub Kruzik   PetscFunctionBegin;
584409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
585409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
586f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
5876c93e71cSJakub Kruzik   if (!def->lvl && !def->prefix) {
5886c93e71cSJakub Kruzik     ierr = PCGetOptionsPrefix(pc,&def->prefix);CHKERRQ(ierr);
5896c93e71cSJakub Kruzik   }
5906c93e71cSJakub Kruzik   if (def->lvl) {
5916c93e71cSJakub Kruzik     sprintf(prefix,"%d_",(int)def->lvl);
5926c93e71cSJakub Kruzik   }
593f8bf229cSJakub Kruzik 
594f8bf229cSJakub Kruzik   /* compute a deflation space */
595409a357bSJakub Kruzik   if (def->W || def->Wt) {
596409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
597409a357bSJakub Kruzik   } else {
598e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
59937eeb815SJakub Kruzik   }
60037eeb815SJakub Kruzik 
601409a357bSJakub Kruzik   /* nested deflation */
602409a357bSJakub Kruzik   if (def->W) {
603409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
604409a357bSJakub Kruzik     if (match) {
605409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
606409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
60737eeb815SJakub Kruzik     }
608409a357bSJakub Kruzik   } else {
609a3177931SJakub Kruzik     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
610409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
611409a357bSJakub Kruzik     if (match) {
612409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
613409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
614409a357bSJakub Kruzik     }
615409a357bSJakub Kruzik     transp = PETSC_TRUE;
616409a357bSJakub Kruzik   }
617409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
618409a357bSJakub Kruzik     if (!transp) {
6196c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
62028da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
621409a357bSJakub Kruzik         for (i=0; i<size; i++) {
622409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
623409a357bSJakub Kruzik         }
624409a357bSJakub Kruzik         size -= 1;
625409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
626409a357bSJakub Kruzik         def->W = mats[size];
627409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
628409a357bSJakub Kruzik         if (size > 1) {
629409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
630409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
631409a357bSJakub Kruzik         } else {
632409a357bSJakub Kruzik           nextDef = mats[0];
633409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
634409a357bSJakub Kruzik         }
63528da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
636409a357bSJakub Kruzik       } else {
637*1fdb00f9SJakub Kruzik         /* TODO test merge side performance */
638409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
639409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
640409a357bSJakub Kruzik       }
64128da0170SJakub Kruzik     } else {
6426c93e71cSJakub Kruzik       if (def->lvl < def->maxlvl) {
64328da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
64428da0170SJakub Kruzik         for (i=0; i<size; i++) {
64528da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
64628da0170SJakub Kruzik         }
64728da0170SJakub Kruzik         size -= 1;
64828da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
64928da0170SJakub Kruzik         def->Wt = mats[0];
65028da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
65128da0170SJakub Kruzik         if (size > 1) {
65228da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
65328da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
65428da0170SJakub Kruzik         } else {
65528da0170SJakub Kruzik           nextDef = mats[1];
65628da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
657409a357bSJakub Kruzik         }
658409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
65928da0170SJakub Kruzik       } else {
66028da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
66128da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
66228da0170SJakub Kruzik       }
66328da0170SJakub Kruzik     }
66428da0170SJakub Kruzik   }
66528da0170SJakub Kruzik 
66628da0170SJakub Kruzik   if (transp) {
66728da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
668a3177931SJakub Kruzik     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
669409a357bSJakub Kruzik   }
670409a357bSJakub Kruzik 
671ae029463SJakub Kruzik   /* assemble WtA */
672ae029463SJakub Kruzik   if (!def->WtA) {
673ae029463SJakub Kruzik     if (def->Wt) {
674ae029463SJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
675ae029463SJakub Kruzik     } else {
676a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
677a3177931SJakub Kruzik       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
6789e56ec8aSJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
679a3177931SJakub Kruzik #else
680ae029463SJakub Kruzik       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
681a3177931SJakub Kruzik #endif
682ae029463SJakub Kruzik     }
683ae029463SJakub Kruzik   }
684409a357bSJakub Kruzik   /* setup coarse problem */
685409a357bSJakub Kruzik   if (!def->WtAWinv) {
68628da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
687409a357bSJakub Kruzik     if (!def->WtAW) {
688ae029463SJakub Kruzik       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
689409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
690409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
691409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
692409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
693ae029463SJakub Kruzik       /* Check columns of W are not in kernel of A */
694409a357bSJakub Kruzik       PetscReal *norms;
695409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
696409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
697409a357bSJakub Kruzik       for (i=0; i<m; i++) {
698409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
699409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
700409a357bSJakub Kruzik         }
701409a357bSJakub Kruzik       }
702409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
703409a357bSJakub Kruzik #endif
704409a357bSJakub Kruzik     } else {
705409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
706409a357bSJakub Kruzik     }
707*1fdb00f9SJakub Kruzik     /* TODO use MATINV ? */
708409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
709409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
710409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
711557a2f7dSJakub Kruzik     /* Setup KSP and PC */
712557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
713557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
714557a2f7dSJakub Kruzik       /* set default KSPtype */
715557a2f7dSJakub Kruzik       if (!def->ksptype) {
716557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
717557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
718557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
719557a2f7dSJakub Kruzik         }
720557a2f7dSJakub Kruzik       }
721ae029463SJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
722557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
723557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
7246c93e71cSJakub Kruzik       ierr = PCDeflationSetLvl_Deflation(pcinner,def->lvl+1,def->maxlvl);CHKERRQ(ierr);
725ae029463SJakub Kruzik       /* inherit options */
7266c93e71cSJakub Kruzik       if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix;
7279f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->init          = def->init;
728557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype       = def->ksptype;
729557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct       = def->correct;
7309f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correctfact   = def->correctfact;
7319f604af8SJakub Kruzik       ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact;
732557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
733557a2f7dSJakub Kruzik     } else { /* the last level */
734557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
735409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
7366c93e71cSJakub Kruzik       /* do not overwrite PCTELESCOPE */
7376c93e71cSJakub Kruzik       if (def->prefix) {
7386c93e71cSJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,def->prefix);CHKERRQ(ierr);
739ac66f006SJakub Kruzik       }
7406c93e71cSJakub Kruzik       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"def_tel_");CHKERRQ(ierr);
741ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
742409a357bSJakub Kruzik       /* Reduction factor choice */
743409a357bSJakub Kruzik       red = def->reductionfact;
744409a357bSJakub Kruzik       if (red < 0) {
745409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
746409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
747409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
748409a357bSJakub Kruzik         if (match) red = commsize;
7496c93e71cSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr);
750409a357bSJakub Kruzik       }
751409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
752ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
753409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
754ac66f006SJakub Kruzik       if (innerksp) {
755409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
756ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
757481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU)
758481b7641SJakub Kruzik         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
759481b7641SJakub Kruzik         if (match) {
760ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
761481b7641SJakub Kruzik         }
762481b7641SJakub Kruzik #endif
763481b7641SJakub Kruzik #if defined(PETSC_HAVE_SUPERLU_DIST)
764481b7641SJakub Kruzik         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
765481b7641SJakub Kruzik         if (match) {
766ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
767409a357bSJakub Kruzik         }
768481b7641SJakub Kruzik #endif
769409a357bSJakub Kruzik       }
770557a2f7dSJakub Kruzik     }
771557a2f7dSJakub Kruzik 
772557a2f7dSJakub Kruzik     if (innerksp) {
7736c93e71cSJakub Kruzik       if (def->prefix) {
7746c93e71cSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,def->prefix);CHKERRQ(ierr);
77522b0793eSJakub Kruzik         ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
7766c93e71cSJakub Kruzik       } else {
7776c93e71cSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
7786c93e71cSJakub Kruzik       }
7796c93e71cSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
780557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
781557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
782ac66f006SJakub Kruzik     }
783409a357bSJakub Kruzik   }
784557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
785557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
786409a357bSJakub Kruzik 
787*1fdb00f9SJakub Kruzik   /* create preconditioner */
78822b0793eSJakub Kruzik   if (!def->pc) {
78922b0793eSJakub Kruzik     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
79022b0793eSJakub Kruzik     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
79122b0793eSJakub Kruzik     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
7926c93e71cSJakub Kruzik     if (def->prefix) {
7936c93e71cSJakub Kruzik       ierr = PCSetOptionsPrefix(def->pc,def->prefix);CHKERRQ(ierr);
79422b0793eSJakub Kruzik     }
7956c93e71cSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"def_");CHKERRQ(ierr);
7966c93e71cSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
7976c93e71cSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"pc_");CHKERRQ(ierr);
79822b0793eSJakub Kruzik     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
79922b0793eSJakub Kruzik     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
80022b0793eSJakub Kruzik   }
80122b0793eSJakub Kruzik 
802f8bf229cSJakub Kruzik   /* create work vecs */
803f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
804f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
80537eeb815SJakub Kruzik   PetscFunctionReturn(0);
80637eeb815SJakub Kruzik }
807b30d4775SJakub Kruzik 
80837eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
80937eeb815SJakub Kruzik {
810b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
81137eeb815SJakub Kruzik   PetscErrorCode    ierr;
81237eeb815SJakub Kruzik 
81337eeb815SJakub Kruzik   PetscFunctionBegin;
814b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
815b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
816b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
817b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
818ae029463SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
819b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
820b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
82122b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
82237eeb815SJakub Kruzik   PetscFunctionReturn(0);
82337eeb815SJakub Kruzik }
82437eeb815SJakub Kruzik 
82537eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
82637eeb815SJakub Kruzik {
82737eeb815SJakub Kruzik   PetscErrorCode ierr;
82837eeb815SJakub Kruzik 
82937eeb815SJakub Kruzik   PetscFunctionBegin;
83037eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
83137eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
83237eeb815SJakub Kruzik   PetscFunctionReturn(0);
83337eeb815SJakub Kruzik }
83437eeb815SJakub Kruzik 
8358513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
8368513960bSJakub Kruzik {
8378513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
8381ac6250aSJakub Kruzik   PetscInt          its;
8398513960bSJakub Kruzik   PetscBool         iascii;
8408513960bSJakub Kruzik   PetscErrorCode    ierr;
8418513960bSJakub Kruzik 
8428513960bSJakub Kruzik   PetscFunctionBegin;
8438513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8448513960bSJakub Kruzik   if (iascii) {
8458513960bSJakub Kruzik     if (def->correct) {
8461ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n",
8471ac6250aSJakub Kruzik                                     (double)PetscRealPart(def->correctfact),
8481ac6250aSJakub Kruzik                                     (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr);
8498513960bSJakub Kruzik     }
8506c93e71cSJakub Kruzik     if (!def->lvl) {
8511ac6250aSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
8528513960bSJakub Kruzik     }
8538513960bSJakub Kruzik 
8541ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr);
85522b0793eSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
85622b0793eSJakub Kruzik     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
85722b0793eSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
85822b0793eSJakub Kruzik 
8591ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr);
8608513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
8611ac6250aSJakub Kruzik     ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr);
8621ac6250aSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr);
8638513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
8648513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
8658513960bSJakub Kruzik   }
8668513960bSJakub Kruzik   PetscFunctionReturn(0);
8678513960bSJakub Kruzik }
8688513960bSJakub Kruzik 
86937eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
87037eeb815SJakub Kruzik {
871880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
87237eeb815SJakub Kruzik   PetscErrorCode    ierr;
87337eeb815SJakub Kruzik 
87437eeb815SJakub Kruzik   PetscFunctionBegin;
87537eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
876a122ebaeSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr);
8776c93e71cSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxlvl,&def->maxlvl,NULL);CHKERRQ(ierr);
878859a873cSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
8798a71cb68SJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr);
8808a71cb68SJakub 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);
881880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
882880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
883880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
88437eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
88537eeb815SJakub Kruzik   PetscFunctionReturn(0);
88637eeb815SJakub Kruzik }
88737eeb815SJakub Kruzik 
88837eeb815SJakub Kruzik /*MC
889e26ad46dSJakub Kruzik    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
89037eeb815SJakub Kruzik 
891e26ad46dSJakub Kruzik    Options Database Keys:
892e26ad46dSJakub Kruzik +    -pc_deflation_init_only          <false> - if true computes only the special guess
893e26ad46dSJakub Kruzik .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
894e26ad46dSJakub Kruzik .    -pc_deflation_reduction_factor <-1>      - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
895e26ad46dSJakub Kruzik .    -pc_deflation_correction         <false> - if true apply coarse problem correction
896e26ad46dSJakub Kruzik .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
897e26ad46dSJakub Kruzik .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
898e26ad46dSJakub Kruzik -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
89937eeb815SJakub Kruzik 
90037eeb815SJakub Kruzik    Notes:
901e26ad46dSJakub Kruzik     Given a (complex - transpose is always hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
902e26ad46dSJakub Kruzik     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
903e26ad46dSJakub Kruzik 
904e26ad46dSJakub Kruzik     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
905e26ad46dSJakub Kruzik     If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the
906e26ad46dSJakub Kruzik     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
907*1fdb00f9SJakub Kruzik     application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
908e26ad46dSJakub Kruzik     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
909e26ad46dSJakub Kruzik     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
910e26ad46dSJakub Kruzik     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
911e26ad46dSJakub Kruzik 
912e26ad46dSJakub Kruzik     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
913e26ad46dSJakub Kruzik     be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
914e26ad46dSJakub Kruzik     User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix
915e26ad46dSJakub Kruzik     is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the
916*1fdb00f9SJakub Kruzik     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned
917e26ad46dSJakub Kruzik     by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum
918*1fdb00f9SJakub Kruzik     level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged
919*1fdb00f9SJakub Kruzik     (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by
920e26ad46dSJakub Kruzik     PCDeflationSetMaxLvl() or by -pc_deflation_max_lvl.
921e26ad46dSJakub Kruzik 
9226c93e71cSJakub Kruzik     The coarse problem KSP can be controlled from the command line with prefix -def_ for the first level and -def_[lvl-1]
9236c93e71cSJakub Kruzik     from the second level onward. You can also use
924e26ad46dSJakub Kruzik     PCDeflationGetCoarseKSP() or PCDeflationSetCoarseKSP() to control it from code. The bottom level KSP defaults to
925*1fdb00f9SJakub Kruzik     KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE.
926481b7641SJakub Kruzik     For convenience, the reduction factor can be set by PCDeflationSetReductionFactor()
927481b7641SJakub Kruzik     or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
928e26ad46dSJakub Kruzik 
9296c93e71cSJakub Kruzik     The additional preconditioner can be controlled from command line with prefix -def_[lvl]_pc (same rules used for
9306c93e71cSJakub Kruzik     coarse problem KSP apply for [lvl]_ part of prefix), e.g., -def_1_pc_pc_type bjacobi. You can also use
9316c93e71cSJakub Kruzik     PCDeflationGetPC() or PCDeflationSetPC() to control the additional preconditioner from code. It defaults to PCNONE.
932e26ad46dSJakub Kruzik 
933e26ad46dSJakub Kruzik     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
934e26ad46dSJakub Kruzik     be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can
935e26ad46dSJakub Kruzik     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
936e26ad46dSJakub Kruzik     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
937e26ad46dSJakub Kruzik     an isolated eigenvalue.
938e26ad46dSJakub Kruzik 
939*1fdb00f9SJakub Kruzik     The options are automatically inherited from the previous deflation level.
9409f604af8SJakub Kruzik 
941e26ad46dSJakub Kruzik     The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also
942*1fdb00f9SJakub Kruzik     recommend limiting the number of iterations for the coarse problems.
943e26ad46dSJakub Kruzik 
944e26ad46dSJakub Kruzik     See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
945e26ad46dSJakub Kruzik     Section 4 describes some possible choices for the deflation space.
946e26ad46dSJakub Kruzik 
947e26ad46dSJakub Kruzik    Developer Notes:
948e26ad46dSJakub Kruzik      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
949e26ad46dSJakub Kruzik      Academy of Sciences and VSB - TU Ostrava.
950e26ad46dSJakub Kruzik 
951e26ad46dSJakub Kruzik      Developed from PERMON code used in [4] while on a research stay with
952e26ad46dSJakub Kruzik      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
953e26ad46dSJakub Kruzik 
954e26ad46dSJakub Kruzik    References:
955e26ad46dSJakub Kruzik +    [1] - A. Nicolaides. “Deflation of conjugate gradients with applications to boundary valueproblems”, SIAM J. Numer. Anal. 24.2, 1987.
956e26ad46dSJakub Kruzik .    [2] - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
957e26ad46dSJakub 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.
958*1fdb00f9SJakub 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
959e26ad46dSJakub Kruzik 
960e26ad46dSJakub Kruzik    Level: intermediate
96137eeb815SJakub Kruzik 
96237eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
963e26ad46dSJakub Kruzik            PCDeflationSetInitOnly(), PCDeflationSetMaxLvl(), PCDeflationSetReductionFactor(),
964*1fdb00f9SJakub Kruzik            PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompute(),
965e26ad46dSJakub Kruzik            PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(),
966e26ad46dSJakub Kruzik            PCDeflationSetCoarseMat(), PCDeflationSetCoarseKSP(), PCDeflationGetCoarseKSP(),
967e26ad46dSJakub Kruzik            PCDeflationGetPC(), PCDeflationSetPC()
96837eeb815SJakub Kruzik M*/
96937eeb815SJakub Kruzik 
97037eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
97137eeb815SJakub Kruzik {
97237eeb815SJakub Kruzik   PC_Deflation   *def;
97337eeb815SJakub Kruzik   PetscErrorCode ierr;
97437eeb815SJakub Kruzik 
97537eeb815SJakub Kruzik   PetscFunctionBegin;
97637eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
97737eeb815SJakub Kruzik   pc->data = (void*)def;
97837eeb815SJakub Kruzik 
979e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
980e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
981fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
982e662bc50SJakub Kruzik   def->reductionfact = -1;
983e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
984e662bc50SJakub Kruzik   def->spacesize     = 1;
985e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
9866c93e71cSJakub Kruzik   def->lvl           = 0;
9876c93e71cSJakub Kruzik   def->maxlvl        = 0;
98837eeb815SJakub Kruzik 
98937eeb815SJakub Kruzik   pc->ops->apply          = PCApply_Deflation;
990b27c8092SJakub Kruzik   pc->ops->presolve       = PCPreSolve_Deflation;
99137eeb815SJakub Kruzik   pc->ops->setup          = PCSetUp_Deflation;
99237eeb815SJakub Kruzik   pc->ops->reset          = PCReset_Deflation;
99337eeb815SJakub Kruzik   pc->ops->destroy        = PCDestroy_Deflation;
99437eeb815SJakub Kruzik   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
9958513960bSJakub Kruzik   pc->ops->view           = PCView_Deflation;
99637eeb815SJakub Kruzik 
997a122ebaeSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr);
99898d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
999859a873cSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr);
10008a71cb68SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr);
100139417d7eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr);
1002e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
10033cf3a049SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr);
1004e924b002SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
10054a99276eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
1006f3eaa4f0SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
100798d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
100898d6e3deSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
100937eeb815SJakub Kruzik   PetscFunctionReturn(0);
101037eeb815SJakub Kruzik }
101137eeb815SJakub Kruzik 
1012