xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision a317793185f69e008bb18a3bdf6e2d19aed1bae2)
137eeb815SJakub Kruzik 
237eeb815SJakub Kruzik /*  --------------------------------------------------------------------
337eeb815SJakub Kruzik 
437eeb815SJakub Kruzik      This file implements a Deflation preconditioner in PETSc as part of PC.
537eeb815SJakub Kruzik      You can use this as a starting point for implementing your own
637eeb815SJakub Kruzik      preconditioner that is not provided with PETSc. (You might also consider
737eeb815SJakub Kruzik      just using PCSHELL)
837eeb815SJakub Kruzik 
937eeb815SJakub Kruzik      The following basic routines are required for each preconditioner.
1037eeb815SJakub Kruzik           PCCreate_XXX()          - Creates a preconditioner context
1137eeb815SJakub Kruzik           PCSetFromOptions_XXX()  - Sets runtime options
1237eeb815SJakub Kruzik           PCApply_XXX()           - Applies the preconditioner
1337eeb815SJakub Kruzik           PCDestroy_XXX()         - Destroys the preconditioner context
1437eeb815SJakub Kruzik      where the suffix "_XXX" denotes a particular implementation, in
1537eeb815SJakub Kruzik      this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation).
1637eeb815SJakub Kruzik      These routines are actually called via the common user interface
1737eeb815SJakub Kruzik      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
1837eeb815SJakub Kruzik      so the application code interface remains identical for all
1937eeb815SJakub Kruzik      preconditioners.
2037eeb815SJakub Kruzik 
2137eeb815SJakub Kruzik      Another key routine is:
2237eeb815SJakub Kruzik           PCSetUp_XXX()           - Prepares for the use of a preconditioner
2337eeb815SJakub Kruzik      by setting data structures and options.   The interface routine PCSetUp()
2437eeb815SJakub Kruzik      is not usually called directly by the user, but instead is called by
2537eeb815SJakub Kruzik      PCApply() if necessary.
2637eeb815SJakub Kruzik 
2737eeb815SJakub Kruzik      Additional basic routines are:
2837eeb815SJakub Kruzik           PCView_XXX()            - Prints details of runtime options that
2937eeb815SJakub Kruzik                                     have actually been used.
3037eeb815SJakub Kruzik      These are called by application codes via the interface routines
3137eeb815SJakub Kruzik      PCView().
3237eeb815SJakub Kruzik 
3337eeb815SJakub Kruzik      The various types of solvers (preconditioners, Krylov subspace methods,
3437eeb815SJakub Kruzik      nonlinear solvers, timesteppers) are all organized similarly, so the
3537eeb815SJakub Kruzik      above description applies to these categories also.  One exception is
3637eeb815SJakub Kruzik      that the analogues of PCApply() for these components are KSPSolve(),
3737eeb815SJakub Kruzik      SNESSolve(), and TSSolve().
3837eeb815SJakub Kruzik 
3937eeb815SJakub Kruzik      Additional optional functionality unique to preconditioners is left and
4037eeb815SJakub Kruzik      right symmetric preconditioner application via PCApplySymmetricLeft()
4137eeb815SJakub Kruzik      and PCApplySymmetricRight().  The Deflation implementation is
4237eeb815SJakub Kruzik      PCApplySymmetricLeftOrRight_Deflation().
4337eeb815SJakub Kruzik 
4437eeb815SJakub Kruzik     -------------------------------------------------------------------- */
4537eeb815SJakub Kruzik 
4637eeb815SJakub Kruzik /*
4737eeb815SJakub Kruzik    Include files needed for the Deflation preconditioner:
4837eeb815SJakub Kruzik      pcimpl.h - private include file intended for use by all preconditioners
4937eeb815SJakub Kruzik */
5037eeb815SJakub Kruzik 
51292e2e67SJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
5237eeb815SJakub Kruzik 
538513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = {
548513960bSJakub Kruzik   "haar",
558513960bSJakub Kruzik   "jacket-haar",
568513960bSJakub Kruzik   "db2",
578513960bSJakub Kruzik   "db4",
588513960bSJakub Kruzik   "db8",
598513960bSJakub Kruzik   "db16",
608513960bSJakub Kruzik   "biorth22",
618513960bSJakub Kruzik   "meyer",
628513960bSJakub Kruzik   "aggregation",
638513960bSJakub Kruzik   "slepc",
648513960bSJakub Kruzik   "slepc-cheap",
658513960bSJakub Kruzik   "user",
668513960bSJakub Kruzik   "DdefSpaceType",
678513960bSJakub Kruzik   "Ddef_SPACE_",
688513960bSJakub Kruzik   0
698513960bSJakub Kruzik };
708513960bSJakub Kruzik 
718513960bSJakub Kruzik 
72e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
73e662bc50SJakub Kruzik {
74e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
75e662bc50SJakub Kruzik   PetscErrorCode ierr;
76e662bc50SJakub Kruzik 
77e662bc50SJakub Kruzik   PetscFunctionBegin;
78e662bc50SJakub Kruzik   if (transpose) {
79e662bc50SJakub Kruzik     def->Wt = W;
80e662bc50SJakub Kruzik     def->W = NULL;
81e662bc50SJakub Kruzik   } else {
82e662bc50SJakub Kruzik     def->W = W;
83e662bc50SJakub Kruzik   }
84e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
85e662bc50SJakub Kruzik   PetscFunctionReturn(0);
86e662bc50SJakub Kruzik }
87e662bc50SJakub Kruzik 
88e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */
89e662bc50SJakub Kruzik /*@
90e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
91e662bc50SJakub Kruzik 
92e662bc50SJakub Kruzik    Logically Collective on PC
93e662bc50SJakub Kruzik 
94e662bc50SJakub Kruzik    Input Parameters:
95e662bc50SJakub Kruzik +  pc - the preconditioner context
96e662bc50SJakub Kruzik .  W  - deflation matrix
97e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
98e662bc50SJakub Kruzik 
99e662bc50SJakub Kruzik    Level: intermediate
100e662bc50SJakub Kruzik 
101e662bc50SJakub Kruzik .seealso: PCDEFLATION
102e662bc50SJakub Kruzik @*/
103e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
104e662bc50SJakub Kruzik {
105e662bc50SJakub Kruzik   PetscErrorCode ierr;
106e662bc50SJakub Kruzik 
107e662bc50SJakub Kruzik   PetscFunctionBegin;
108e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
109e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
110e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
111e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
112e662bc50SJakub Kruzik   PetscFunctionReturn(0);
113e662bc50SJakub Kruzik }
114e662bc50SJakub Kruzik 
115e924b002SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
116e924b002SJakub Kruzik {
117e924b002SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
118e924b002SJakub Kruzik   PetscErrorCode   ierr;
119e924b002SJakub Kruzik 
120e924b002SJakub Kruzik   PetscFunctionBegin;
121e924b002SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
122e924b002SJakub Kruzik   def->WtAW = mat;
123e924b002SJakub Kruzik   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
124e924b002SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
125e924b002SJakub Kruzik   PetscFunctionReturn(0);
126e924b002SJakub Kruzik }
127e924b002SJakub Kruzik 
128e924b002SJakub Kruzik /*@
129e924b002SJakub Kruzik    PCDeflationSetCoarseMat - Set coarse problem Mat.
130e924b002SJakub Kruzik 
131e924b002SJakub Kruzik    Not Collective
132e924b002SJakub Kruzik 
133e924b002SJakub Kruzik    Input Parameters:
134e924b002SJakub Kruzik +  pc - preconditioner context
135e924b002SJakub Kruzik -  mat - coarse problem mat
136e924b002SJakub Kruzik 
137e924b002SJakub Kruzik    Level: developer
138e924b002SJakub Kruzik 
139e924b002SJakub Kruzik .seealso: PCDEFLATION
140e924b002SJakub Kruzik @*/
141e924b002SJakub Kruzik PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
142e924b002SJakub Kruzik {
143e924b002SJakub Kruzik   PetscErrorCode ierr;
144e924b002SJakub Kruzik 
145e924b002SJakub Kruzik   PetscFunctionBegin;
146e924b002SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
147e924b002SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
148e924b002SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
149e924b002SJakub Kruzik   PetscFunctionReturn(0);
150e924b002SJakub Kruzik }
151e924b002SJakub Kruzik 
1525c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
1535c0c31c5SJakub Kruzik {
1545c0c31c5SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
1555c0c31c5SJakub Kruzik 
1565c0c31c5SJakub Kruzik   PetscFunctionBegin;
15781e2347eSJakub Kruzik   if (current) def->nestedlvl = current;
1585c0c31c5SJakub Kruzik   def->maxnestedlvl = max;
1595c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1605c0c31c5SJakub Kruzik }
1615c0c31c5SJakub Kruzik 
1625c0c31c5SJakub Kruzik /*@
1635c0c31c5SJakub Kruzik    PCDeflationSetMaxLvl - Set maximum level of deflation.
1645c0c31c5SJakub Kruzik 
1655c0c31c5SJakub Kruzik    Logically Collective on PC
1665c0c31c5SJakub Kruzik 
1675c0c31c5SJakub Kruzik    Input Parameters:
1685c0c31c5SJakub Kruzik +  pc  - the preconditioner context
16922b0793eSJakub Kruzik -  max - maximum deflation level
1705c0c31c5SJakub Kruzik 
1715c0c31c5SJakub Kruzik    Level: intermediate
1725c0c31c5SJakub Kruzik 
1735c0c31c5SJakub Kruzik .seealso: PCDEFLATION
1745c0c31c5SJakub Kruzik @*/
1755c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
1765c0c31c5SJakub Kruzik {
1775c0c31c5SJakub Kruzik   PetscErrorCode ierr;
1785c0c31c5SJakub Kruzik 
1795c0c31c5SJakub Kruzik   PetscFunctionBegin;
18022b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1815c0c31c5SJakub Kruzik   PetscValidLogicalCollectiveInt(pc,max,2);
1825c0c31c5SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
1835c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1845c0c31c5SJakub Kruzik }
185e662bc50SJakub Kruzik 
186268c6673SJakub Kruzik static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
187268c6673SJakub Kruzik {
188268c6673SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
189268c6673SJakub Kruzik 
190268c6673SJakub Kruzik   PetscFunctionBegin;
191268c6673SJakub Kruzik   *apc = def->pc;
192268c6673SJakub Kruzik   PetscFunctionReturn(0);
193268c6673SJakub Kruzik }
194268c6673SJakub Kruzik 
195268c6673SJakub Kruzik /*@
196268c6673SJakub Kruzik    PCDeflationGetPC - Returns a pointer to additional preconditioner.
197268c6673SJakub Kruzik 
198268c6673SJakub Kruzik    Not Collective
199268c6673SJakub Kruzik 
200268c6673SJakub Kruzik    Input Parameters:
201268c6673SJakub Kruzik .  pc  - the preconditioner context
202268c6673SJakub Kruzik 
203268c6673SJakub Kruzik    Output Parameter:
204268c6673SJakub Kruzik .  apc - additional preconditioner
205268c6673SJakub Kruzik 
206268c6673SJakub Kruzik    Level: advanced
207268c6673SJakub Kruzik 
208268c6673SJakub Kruzik .seealso: PCDeflationSetPC(), PCDEFLATION
209268c6673SJakub Kruzik @*/
210268c6673SJakub Kruzik PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
211268c6673SJakub Kruzik {
212268c6673SJakub Kruzik   PetscErrorCode ierr;
213268c6673SJakub Kruzik 
214268c6673SJakub Kruzik   PetscFunctionBegin;
215268c6673SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
216268c6673SJakub Kruzik   PetscValidPointer(pc,2);
217268c6673SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
218268c6673SJakub Kruzik   PetscFunctionReturn(0);
219268c6673SJakub Kruzik }
220268c6673SJakub Kruzik 
22122b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
22222b0793eSJakub Kruzik {
22322b0793eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
22422b0793eSJakub Kruzik   PetscErrorCode ierr;
22522b0793eSJakub Kruzik 
22622b0793eSJakub Kruzik   PetscFunctionBegin;
22722b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
22822b0793eSJakub Kruzik   def->pc = apc;
22922b0793eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
23022b0793eSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
23122b0793eSJakub Kruzik   PetscFunctionReturn(0);
23222b0793eSJakub Kruzik }
23322b0793eSJakub Kruzik 
23422b0793eSJakub Kruzik /*@
23522b0793eSJakub Kruzik    PCDeflationSetPC - Set additional preconditioner.
23622b0793eSJakub Kruzik 
237268c6673SJakub Kruzik    Collective on PC
23822b0793eSJakub Kruzik 
23922b0793eSJakub Kruzik    Input Parameters:
24022b0793eSJakub Kruzik +  pc  - the preconditioner context
24122b0793eSJakub Kruzik -  apc - additional preconditioner
24222b0793eSJakub Kruzik 
243268c6673SJakub Kruzik    Level: developer
24422b0793eSJakub Kruzik 
245268c6673SJakub Kruzik .seealso: PCDeflationGetPC(), PCDEFLATION
24622b0793eSJakub Kruzik @*/
24722b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
24822b0793eSJakub Kruzik {
24922b0793eSJakub Kruzik   PetscErrorCode ierr;
25022b0793eSJakub Kruzik 
25122b0793eSJakub Kruzik   PetscFunctionBegin;
25222b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
25322b0793eSJakub Kruzik   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
25422b0793eSJakub Kruzik   PetscCheckSameComm(pc,1,apc,2);
25522b0793eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
25622b0793eSJakub Kruzik   PetscFunctionReturn(0);
25722b0793eSJakub Kruzik }
25822b0793eSJakub Kruzik 
2594a99276eSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
2604a99276eSJakub Kruzik {
2614a99276eSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
2624a99276eSJakub Kruzik 
2634a99276eSJakub Kruzik   PetscFunctionBegin;
2644a99276eSJakub Kruzik   *ksp = def->WtAWinv;
2654a99276eSJakub Kruzik   PetscFunctionReturn(0);
2664a99276eSJakub Kruzik }
2674a99276eSJakub Kruzik 
2684a99276eSJakub Kruzik /*@
2694a99276eSJakub Kruzik    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
2704a99276eSJakub Kruzik 
2714a99276eSJakub Kruzik    Not Collective
2724a99276eSJakub Kruzik 
2734a99276eSJakub Kruzik    Input Parameters:
2744a99276eSJakub Kruzik .  pc - preconditioner context
2754a99276eSJakub Kruzik 
2764a99276eSJakub Kruzik    Output Parameter:
2774a99276eSJakub Kruzik .  ksp - coarse problem KSP context
2784a99276eSJakub Kruzik 
2794a99276eSJakub Kruzik    Level: developer
2804a99276eSJakub Kruzik 
281f3eaa4f0SJakub Kruzik .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION
2824a99276eSJakub Kruzik @*/
2834a99276eSJakub Kruzik PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
2844a99276eSJakub Kruzik {
2854a99276eSJakub Kruzik   PetscErrorCode ierr;
2864a99276eSJakub Kruzik 
2874a99276eSJakub Kruzik   PetscFunctionBegin;
2884a99276eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2894a99276eSJakub Kruzik   PetscValidPointer(ksp,2);
2904a99276eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
2914a99276eSJakub Kruzik   PetscFunctionReturn(0);
2924a99276eSJakub Kruzik }
2934a99276eSJakub Kruzik 
294f3eaa4f0SJakub Kruzik static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp)
295f3eaa4f0SJakub Kruzik {
296f3eaa4f0SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
297f3eaa4f0SJakub Kruzik   PetscErrorCode   ierr;
298f3eaa4f0SJakub Kruzik 
299f3eaa4f0SJakub Kruzik   PetscFunctionBegin;
300f3eaa4f0SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
301f3eaa4f0SJakub Kruzik   def->WtAWinv = ksp;
302f3eaa4f0SJakub Kruzik   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
303f3eaa4f0SJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr);
304f3eaa4f0SJakub Kruzik   PetscFunctionReturn(0);
305f3eaa4f0SJakub Kruzik }
306f3eaa4f0SJakub Kruzik 
307f3eaa4f0SJakub Kruzik /*@
308f3eaa4f0SJakub Kruzik    PCDeflationSetCoarseKSP - Set coarse problem KSP.
309f3eaa4f0SJakub Kruzik 
310f3eaa4f0SJakub Kruzik    Collective on PC
311f3eaa4f0SJakub Kruzik 
312f3eaa4f0SJakub Kruzik    Input Parameters:
313f3eaa4f0SJakub Kruzik +  pc - preconditioner context
314f3eaa4f0SJakub Kruzik -  ksp - coarse problem KSP context
315f3eaa4f0SJakub Kruzik 
316f3eaa4f0SJakub Kruzik    Level: developer
317f3eaa4f0SJakub Kruzik 
318f3eaa4f0SJakub Kruzik .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION
319f3eaa4f0SJakub Kruzik @*/
320f3eaa4f0SJakub Kruzik PetscErrorCode  PCDeflationSetCoarseKSP(PC pc,KSP ksp)
321f3eaa4f0SJakub Kruzik {
322f3eaa4f0SJakub Kruzik   PetscErrorCode ierr;
323f3eaa4f0SJakub Kruzik 
324f3eaa4f0SJakub Kruzik   PetscFunctionBegin;
325f3eaa4f0SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
326f3eaa4f0SJakub Kruzik   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
327f3eaa4f0SJakub Kruzik   PetscCheckSameComm(pc,1,ksp,2);
328f3eaa4f0SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr);
329f3eaa4f0SJakub Kruzik   PetscFunctionReturn(0);
330f3eaa4f0SJakub Kruzik }
331f3eaa4f0SJakub Kruzik 
33237eeb815SJakub Kruzik /*
333b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
334b27c8092SJakub Kruzik */
335b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
336b27c8092SJakub Kruzik {
337b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
338b27c8092SJakub Kruzik   Mat              A;
339b27c8092SJakub Kruzik   Vec              r,w1,w2;
340b27c8092SJakub Kruzik   PetscBool        nonzero;
341b27c8092SJakub Kruzik   PetscErrorCode   ierr;
34237eeb815SJakub Kruzik 
343b27c8092SJakub Kruzik   PetscFunctionBegin;
344b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
345b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
346b27c8092SJakub Kruzik   r  = def->work;
347b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
34837eeb815SJakub Kruzik 
349b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
350b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
351b27c8092SJakub Kruzik   if (nonzero) {
352b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);                          /*    r  <- b - Ax              */
353b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
354b27c8092SJakub Kruzik   } else {
355b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                            /*    r  <- b (x is 0)          */
356b27c8092SJakub Kruzik   }
357b27c8092SJakub Kruzik 
358*a3177931SJakub Kruzik   if (def->Wt) {
359*a3177931SJakub Kruzik     ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr);                   /*    w1 <- W'*r                */
360*a3177931SJakub Kruzik   } else {
361*a3177931SJakub Kruzik     ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
362*a3177931SJakub Kruzik   }
363b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);              /*    w2 <- (W'*A*W)^{-1}*w1    */
364b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);                      /*    r  <- W*w2                */
365b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
366b27c8092SJakub Kruzik   PetscFunctionReturn(0);
367b27c8092SJakub Kruzik }
36837eeb815SJakub Kruzik 
369f8bf229cSJakub Kruzik /*
370f8bf229cSJakub Kruzik   if (def->correct) {
371ae029463SJakub Kruzik     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
372f8bf229cSJakub Kruzik   } else {
373ae029463SJakub Kruzik     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
374f8bf229cSJakub Kruzik   }
37537eeb815SJakub Kruzik */
376f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
377f8bf229cSJakub Kruzik {
378f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
379f8bf229cSJakub Kruzik   Mat              A;
380f8bf229cSJakub Kruzik   Vec              u,w1,w2;
381f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
382f8bf229cSJakub Kruzik 
383f8bf229cSJakub Kruzik   PetscFunctionBegin;
384f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
385f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
386f8bf229cSJakub Kruzik   u  = def->work;
387f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
388f8bf229cSJakub Kruzik 
389ae029463SJakub Kruzik   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                      /*    z <- M^{-1}*r             */
39022b0793eSJakub Kruzik   if (!def->init) {
391ae029463SJakub Kruzik     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                  /*    w1 <- W'*A*z              */
392f8bf229cSJakub Kruzik     if (def->correct) {
393ae029463SJakub Kruzik       if (def->Wt) {
394ae029463SJakub Kruzik         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);               /*    w2 <- W'*r                */
395ae029463SJakub Kruzik       } else {
396*a3177931SJakub Kruzik         ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr);       /*    w2 <- W'*r                */
397f8bf229cSJakub Kruzik       }
398ae029463SJakub Kruzik       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);  /*    w1 <- w1 - l*w2           */
399f8bf229cSJakub Kruzik     }
400f8bf229cSJakub Kruzik     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);            /*    w2 <- (W'*A*W)^{-1}*w1    */
401f8bf229cSJakub Kruzik     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                    /*    u  <- W*w2                */
40222b0793eSJakub Kruzik     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                       /*    z  <- z - u               */
40322b0793eSJakub Kruzik   }
404f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
405f8bf229cSJakub Kruzik }
406f8bf229cSJakub Kruzik 
40737eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
40837eeb815SJakub Kruzik {
409409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
410409a357bSJakub Kruzik   KSP              innerksp;
411409a357bSJakub Kruzik   PC               pcinner;
412409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
413409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
414409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
415409a357bSJakub Kruzik   MatCompositeType ctype;
416409a357bSJakub Kruzik   MPI_Comm         comm;
417409a357bSJakub Kruzik   const char       *prefix;
41837eeb815SJakub Kruzik   PetscErrorCode   ierr;
41937eeb815SJakub Kruzik 
42037eeb815SJakub Kruzik   PetscFunctionBegin;
421409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
422409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
423f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
424f8bf229cSJakub Kruzik 
425f8bf229cSJakub Kruzik   /* compute a deflation space */
426409a357bSJakub Kruzik   if (def->W || def->Wt) {
427409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
428409a357bSJakub Kruzik   } else {
429e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
43037eeb815SJakub Kruzik   }
43137eeb815SJakub Kruzik 
432409a357bSJakub Kruzik   /* nested deflation */
433409a357bSJakub Kruzik   if (def->W) {
434409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
435409a357bSJakub Kruzik     if (match) {
436409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
437409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
43837eeb815SJakub Kruzik     }
439409a357bSJakub Kruzik   } else {
440*a3177931SJakub Kruzik     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
441409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
442409a357bSJakub Kruzik     if (match) {
443409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
444409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
445409a357bSJakub Kruzik     }
446409a357bSJakub Kruzik     transp = PETSC_TRUE;
447409a357bSJakub Kruzik   }
448409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
449409a357bSJakub Kruzik     if (!transp) {
45028da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
45128da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
452409a357bSJakub Kruzik         for (i=0; i<size; i++) {
453409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
454409a357bSJakub Kruzik         }
455409a357bSJakub Kruzik         size -= 1;
456409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
457409a357bSJakub Kruzik         def->W = mats[size];
458409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
459409a357bSJakub Kruzik         if (size > 1) {
460409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
461409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
462409a357bSJakub Kruzik         } else {
463409a357bSJakub Kruzik           nextDef = mats[0];
464409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
465409a357bSJakub Kruzik         }
46628da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
467409a357bSJakub Kruzik       } else {
468409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
469409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
470409a357bSJakub Kruzik       }
47128da0170SJakub Kruzik     } else {
47228da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
47328da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
47428da0170SJakub Kruzik         for (i=0; i<size; i++) {
47528da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
47628da0170SJakub Kruzik         }
47728da0170SJakub Kruzik         size -= 1;
47828da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
47928da0170SJakub Kruzik         def->Wt = mats[0];
48028da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
48128da0170SJakub Kruzik         if (size > 1) {
48228da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
48328da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
48428da0170SJakub Kruzik         } else {
48528da0170SJakub Kruzik           nextDef = mats[1];
48628da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
487409a357bSJakub Kruzik         }
488409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
48928da0170SJakub Kruzik       } else {
49028da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
49128da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
49228da0170SJakub Kruzik       }
49328da0170SJakub Kruzik     }
49428da0170SJakub Kruzik   }
49528da0170SJakub Kruzik 
49628da0170SJakub Kruzik   if (transp) {
49728da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
498*a3177931SJakub Kruzik     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
499409a357bSJakub Kruzik   }
500409a357bSJakub Kruzik 
50122b0793eSJakub Kruzik 
50222b0793eSJakub Kruzik   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
50322b0793eSJakub Kruzik 
504ae029463SJakub Kruzik   /* assemble WtA */
505ae029463SJakub Kruzik   if (!def->WtA) {
506ae029463SJakub Kruzik     if (def->Wt) {
507ae029463SJakub Kruzik       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
508ae029463SJakub Kruzik     } else {
509*a3177931SJakub Kruzik #if defined(PETSC_USE_COMPLEX)
510*a3177931SJakub Kruzik       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
511*a3177931SJakub Kruzik       ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
512*a3177931SJakub Kruzik #else
513ae029463SJakub Kruzik       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
514*a3177931SJakub Kruzik #endif
515ae029463SJakub Kruzik     }
516ae029463SJakub Kruzik   }
517409a357bSJakub Kruzik   /* setup coarse problem */
518409a357bSJakub Kruzik   if (!def->WtAWinv) {
51928da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
520409a357bSJakub Kruzik     if (!def->WtAW) {
521ae029463SJakub Kruzik       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
522409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
523409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
524409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
525409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
526ae029463SJakub Kruzik       /* Check columns of W are not in kernel of A */
527409a357bSJakub Kruzik       PetscReal *norms;
528409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
529409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
530409a357bSJakub Kruzik       for (i=0; i<m; i++) {
531409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
532409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
533409a357bSJakub Kruzik         }
534409a357bSJakub Kruzik       }
535409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
536409a357bSJakub Kruzik #endif
537409a357bSJakub Kruzik     } else {
538409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
539409a357bSJakub Kruzik     }
540409a357bSJakub Kruzik     /* TODO use MATINV */
541409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
542409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
543409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
544557a2f7dSJakub Kruzik     /* Setup KSP and PC */
545557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
546557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
547557a2f7dSJakub Kruzik       /* set default KSPtype */
548557a2f7dSJakub Kruzik       if (!def->ksptype) {
549557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
550557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
551557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
552557a2f7dSJakub Kruzik         }
553557a2f7dSJakub Kruzik       }
554ae029463SJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
555557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
556557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
557557a2f7dSJakub Kruzik       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
558ae029463SJakub Kruzik       /* inherit options */
559557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
560557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
561557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
562557a2f7dSJakub Kruzik     } else { /* the last level */
563557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
564409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
565ac66f006SJakub Kruzik       /* ugly hack to not have overwritten PCTELESCOPE */
566ac66f006SJakub Kruzik       if (prefix) {
567ac66f006SJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
568ac66f006SJakub Kruzik       }
56922b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
570ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
571409a357bSJakub Kruzik       /* Reduction factor choice */
572409a357bSJakub Kruzik       red = def->reductionfact;
573409a357bSJakub Kruzik       if (red < 0) {
574409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
575409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
576409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
577409a357bSJakub Kruzik         if (match) red = commsize;
578409a357bSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
579409a357bSJakub Kruzik       }
580409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
581ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
582409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
583ac66f006SJakub Kruzik       if (innerksp) {
584409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
585409a357bSJakub Kruzik         /* TODO Cholesky if flgspd? */
586ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
587409a357bSJakub Kruzik         //TODO remove explicit matSolverPackage
588409a357bSJakub Kruzik         if (commsize == red) {
589ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
590409a357bSJakub Kruzik         } else {
591ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
592409a357bSJakub Kruzik         }
593409a357bSJakub Kruzik       }
594557a2f7dSJakub Kruzik     }
595557a2f7dSJakub Kruzik 
596557a2f7dSJakub Kruzik     if (innerksp) {
597409a357bSJakub Kruzik       /* TODO use def_[lvl]_ if lvl > 0? */
598409a357bSJakub Kruzik       if (prefix) {
599409a357bSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
600409a357bSJakub Kruzik       }
60122b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
602557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
603557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
604ac66f006SJakub Kruzik     }
605409a357bSJakub Kruzik   }
606557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
607557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
608409a357bSJakub Kruzik 
60922b0793eSJakub Kruzik   if (!def->pc) {
61022b0793eSJakub Kruzik     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
61122b0793eSJakub Kruzik     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
61222b0793eSJakub Kruzik     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
61322b0793eSJakub Kruzik     if (prefix) {
61422b0793eSJakub Kruzik       ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
61522b0793eSJakub Kruzik     }
61622b0793eSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr);
61722b0793eSJakub Kruzik     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
61822b0793eSJakub Kruzik     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
61922b0793eSJakub Kruzik   }
62022b0793eSJakub Kruzik 
621f8bf229cSJakub Kruzik   /* create work vecs */
622f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
623f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
62437eeb815SJakub Kruzik   PetscFunctionReturn(0);
62537eeb815SJakub Kruzik }
626b30d4775SJakub Kruzik 
62737eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
62837eeb815SJakub Kruzik {
629b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
63037eeb815SJakub Kruzik   PetscErrorCode    ierr;
63137eeb815SJakub Kruzik 
63237eeb815SJakub Kruzik   PetscFunctionBegin;
633b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
634b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
635b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
636b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
637ae029463SJakub Kruzik   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
638b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
639b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
64022b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
64137eeb815SJakub Kruzik   PetscFunctionReturn(0);
64237eeb815SJakub Kruzik }
64337eeb815SJakub Kruzik 
64437eeb815SJakub Kruzik /*
64537eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
64637eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
64737eeb815SJakub Kruzik 
64837eeb815SJakub Kruzik    Input Parameter:
64937eeb815SJakub Kruzik .  pc - the preconditioner context
65037eeb815SJakub Kruzik 
65137eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
65237eeb815SJakub Kruzik */
65337eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
65437eeb815SJakub Kruzik {
65537eeb815SJakub Kruzik   PetscErrorCode ierr;
65637eeb815SJakub Kruzik 
65737eeb815SJakub Kruzik   PetscFunctionBegin;
65837eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
65937eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
66037eeb815SJakub Kruzik   PetscFunctionReturn(0);
66137eeb815SJakub Kruzik }
66237eeb815SJakub Kruzik 
6638513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
6648513960bSJakub Kruzik {
6658513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
6668513960bSJakub Kruzik   PetscBool         iascii;
6678513960bSJakub Kruzik   PetscErrorCode    ierr;
6688513960bSJakub Kruzik 
6698513960bSJakub Kruzik   PetscFunctionBegin;
6708513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6718513960bSJakub Kruzik   if (iascii) {
6728513960bSJakub Kruzik     //if (cg->adaptiveconv) {
6738513960bSJakub Kruzik     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
6748513960bSJakub Kruzik     //}
6758513960bSJakub Kruzik     if (def->correct) {
6768513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
6778513960bSJakub Kruzik     }
6788513960bSJakub Kruzik     if (!def->nestedlvl) {
6798513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
6808513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
6818513960bSJakub Kruzik     }
6828513960bSJakub Kruzik 
68322b0793eSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr);
68422b0793eSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
68522b0793eSJakub Kruzik     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
68622b0793eSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
68722b0793eSJakub Kruzik 
6888513960bSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
6898513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6908513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
6918513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
6928513960bSJakub Kruzik   }
6938513960bSJakub Kruzik   PetscFunctionReturn(0);
6948513960bSJakub Kruzik }
6958513960bSJakub Kruzik 
69637eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
69737eeb815SJakub Kruzik {
698880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
69937eeb815SJakub Kruzik   PetscErrorCode    ierr;
70037eeb815SJakub Kruzik 
70137eeb815SJakub Kruzik   PetscFunctionBegin;
70237eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
703880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
704880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
705880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
706880d05ceSJakub Kruzik //TODO add set function and fix manpages
707880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr);
708fcb31d99SJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr);
709fcb31d99SJakub Kruzik   ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
710880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
711880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
71237eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
71337eeb815SJakub Kruzik   PetscFunctionReturn(0);
71437eeb815SJakub Kruzik }
71537eeb815SJakub Kruzik 
71637eeb815SJakub Kruzik /*MC
717e361f8b5SJakub Kruzik      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
718e361f8b5SJakub Kruzik      or to a predefined value
71937eeb815SJakub Kruzik 
72037eeb815SJakub Kruzik    Options Database Key:
721e361f8b5SJakub Kruzik +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
72237eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
72337eeb815SJakub Kruzik 
72437eeb815SJakub Kruzik    Level: beginner
72537eeb815SJakub Kruzik 
72637eeb815SJakub Kruzik   Notes:
727e361f8b5SJakub Kruzik     todo
72837eeb815SJakub Kruzik 
72937eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
730e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
73137eeb815SJakub Kruzik M*/
73237eeb815SJakub Kruzik 
73337eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
73437eeb815SJakub Kruzik {
73537eeb815SJakub Kruzik   PC_Deflation   *def;
73637eeb815SJakub Kruzik   PetscErrorCode ierr;
73737eeb815SJakub Kruzik 
73837eeb815SJakub Kruzik   PetscFunctionBegin;
73937eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
74037eeb815SJakub Kruzik   pc->data = (void*)def;
74137eeb815SJakub Kruzik 
742e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
743e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
744fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
745e662bc50SJakub Kruzik   def->reductionfact = -1;
746e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
747e662bc50SJakub Kruzik   def->spacesize     = 1;
748e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
749e662bc50SJakub Kruzik   def->nestedlvl     = 0;
750e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
75137eeb815SJakub Kruzik 
75237eeb815SJakub Kruzik   /*
75337eeb815SJakub Kruzik       Set the pointers for the functions that are provided above.
75437eeb815SJakub Kruzik       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
75537eeb815SJakub Kruzik       are called, they will automatically call these functions.  Note we
75637eeb815SJakub Kruzik       choose not to provide a couple of these functions since they are
75737eeb815SJakub Kruzik       not needed.
75837eeb815SJakub Kruzik   */
75937eeb815SJakub Kruzik   pc->ops->apply               = PCApply_Deflation;
76037eeb815SJakub Kruzik   pc->ops->applytranspose      = PCApply_Deflation;
761b27c8092SJakub Kruzik   pc->ops->presolve            = PCPreSolve_Deflation;
762b27c8092SJakub Kruzik   pc->ops->postsolve           = 0;
76337eeb815SJakub Kruzik   pc->ops->setup               = PCSetUp_Deflation;
76437eeb815SJakub Kruzik   pc->ops->reset               = PCReset_Deflation;
76537eeb815SJakub Kruzik   pc->ops->destroy             = PCDestroy_Deflation;
76637eeb815SJakub Kruzik   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
7678513960bSJakub Kruzik   pc->ops->view                = PCView_Deflation;
76837eeb815SJakub Kruzik   pc->ops->applyrichardson     = 0;
76937eeb815SJakub Kruzik   pc->ops->applysymmetricleft  = 0;
77037eeb815SJakub Kruzik   pc->ops->applysymmetricright = 0;
77137eeb815SJakub Kruzik 
772e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
7735c0c31c5SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
774268c6673SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
77522b0793eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
776e924b002SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
7774a99276eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
778f3eaa4f0SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
77937eeb815SJakub Kruzik   PetscFunctionReturn(0);
78037eeb815SJakub Kruzik }
78137eeb815SJakub Kruzik 
782