xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 409a357bae74e46200d26aa0a2f21a6fc06036bf)
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 
5137eeb815SJakub Kruzik #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
5237eeb815SJakub Kruzik 
5337eeb815SJakub Kruzik const char *const PCDeflationTypes[]    = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0};
5437eeb815SJakub Kruzik 
5537eeb815SJakub Kruzik /*
5637eeb815SJakub Kruzik    Private context (data structure) for the deflation preconditioner.
5737eeb815SJakub Kruzik */
5837eeb815SJakub Kruzik typedef struct {
5937eeb815SJakub Kruzik   PetscBool init;            /* do only init step - error correction of direction is omitted */
6037eeb815SJakub Kruzik   PetscBool pre;             /* start with x0 being the solution in the deflation space */
6137eeb815SJakub Kruzik   PetscBool correct;         /* add CP (Qr) correction to descent direction */
6237eeb815SJakub Kruzik   PetscBool truenorm;
6337eeb815SJakub Kruzik   PetscBool adaptiveconv;
6437eeb815SJakub Kruzik   PetscReal adaptiveconst;
6537eeb815SJakub Kruzik   PetscInt  reductionfact;
6637eeb815SJakub Kruzik   Mat       W,Wt,AW,WtAW;    /* deflation space, coarse problem mats */
6737eeb815SJakub Kruzik   KSP       WtAWinv;         /* deflation coarse problem */
68*409a357bSJakub Kruzik   KSPType   ksptype;
6937eeb815SJakub Kruzik   Vec       *work;
7037eeb815SJakub Kruzik 
71*409a357bSJakub Kruzik   PCDeflationSpaceType spacetype;
7237eeb815SJakub Kruzik   PetscInt             spacesize;
7337eeb815SJakub Kruzik   PetscInt             nestedlvl;
7437eeb815SJakub Kruzik   PetscInt             maxnestedlvl;
7537eeb815SJakub Kruzik   PetscBool            extendsp;
7637eeb815SJakub Kruzik } PC_Deflation;
7737eeb815SJakub Kruzik 
7837eeb815SJakub Kruzik static PetscErrorCode  PCDeflationSetType_Deflation(PC pc,PCDeflationType type)
7937eeb815SJakub Kruzik {
8037eeb815SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
8137eeb815SJakub Kruzik 
8237eeb815SJakub Kruzik   PetscFunctionBegin;
8337eeb815SJakub Kruzik   def->init = PETSC_FALSE;
8437eeb815SJakub Kruzik   def->pre = PETSC_FALSE;
8537eeb815SJakub Kruzik   if (type == PC_DEFLATION_INIT) {
8637eeb815SJakub Kruzik     def->init = PETSC_TRUE;
8737eeb815SJakub Kruzik     def->pre  = PETSC_TRUE;
8837eeb815SJakub Kruzik   } else if (type == PC_DEFLATION_PRE) {
8937eeb815SJakub Kruzik     def->pre  = PETSC_TRUE;
9037eeb815SJakub Kruzik   }
9137eeb815SJakub Kruzik   PetscFunctionReturn(0);
9237eeb815SJakub Kruzik }
9337eeb815SJakub Kruzik 
94bb9edbfaSJakub Kruzik /*@
95bb9edbfaSJakub Kruzik    PCDeflationSetType - Causes the deflation preconditioner to use only a special
96bb9edbfaSJakub Kruzik     initial gues or pre/post solve solution update
97bb9edbfaSJakub Kruzik 
98bb9edbfaSJakub Kruzik    Logically Collective on PC
99bb9edbfaSJakub Kruzik 
100bb9edbfaSJakub Kruzik    Input Parameters:
101bb9edbfaSJakub Kruzik +  pc - the preconditioner context
102bb9edbfaSJakub Kruzik -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
103bb9edbfaSJakub Kruzik 
104bb9edbfaSJakub Kruzik    Options Database Key:
105bb9edbfaSJakub Kruzik .  -pc_deflation_type <pre,init,post>
106bb9edbfaSJakub Kruzik 
107bb9edbfaSJakub Kruzik    Level: intermediate
108bb9edbfaSJakub Kruzik 
109bb9edbfaSJakub Kruzik    Concepts: Deflation preconditioner
110bb9edbfaSJakub Kruzik 
111bb9edbfaSJakub Kruzik .seealso: PCDeflationGetType()
112bb9edbfaSJakub Kruzik @*/
113bb9edbfaSJakub Kruzik PetscErrorCode  PCDeflationSetType(PC pc,PCDeflationType type)
114bb9edbfaSJakub Kruzik {
115bb9edbfaSJakub Kruzik   PetscErrorCode ierr;
116bb9edbfaSJakub Kruzik 
117bb9edbfaSJakub Kruzik   PetscFunctionBegin;
118bb9edbfaSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119bb9edbfaSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr);
120bb9edbfaSJakub Kruzik   PetscFunctionReturn(0);
121bb9edbfaSJakub Kruzik }
122bb9edbfaSJakub Kruzik 
12337eeb815SJakub Kruzik static PetscErrorCode  PCDeflationGetType_Deflation(PC pc,PCDeflationType *type)
12437eeb815SJakub Kruzik {
12537eeb815SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
12637eeb815SJakub Kruzik 
12737eeb815SJakub Kruzik   PetscFunctionBegin;
12837eeb815SJakub Kruzik   if (def->init) {
12937eeb815SJakub Kruzik     *type = PC_DEFLATION_INIT;
13037eeb815SJakub Kruzik   } else if (def->pre) {
13137eeb815SJakub Kruzik     *type = PC_DEFLATION_PRE;
13237eeb815SJakub Kruzik   } else {
13337eeb815SJakub Kruzik     *type = PC_DEFLATION_POST;
13437eeb815SJakub Kruzik   }
13537eeb815SJakub Kruzik   PetscFunctionReturn(0);
13637eeb815SJakub Kruzik }
13737eeb815SJakub Kruzik 
138bb9edbfaSJakub Kruzik /*@
139bb9edbfaSJakub Kruzik    PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner
140bb9edbfaSJakub Kruzik 
141bb9edbfaSJakub Kruzik    Not Collective on PC
142bb9edbfaSJakub Kruzik 
143bb9edbfaSJakub Kruzik    Input Parameter:
144bb9edbfaSJakub Kruzik .  pc - the preconditioner context
145bb9edbfaSJakub Kruzik 
146bb9edbfaSJakub Kruzik    Output Parameter:
147bb9edbfaSJakub Kruzik -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
148bb9edbfaSJakub Kruzik 
149bb9edbfaSJakub Kruzik    Level: intermediate
150bb9edbfaSJakub Kruzik 
151bb9edbfaSJakub Kruzik    Concepts: Deflation preconditioner
152bb9edbfaSJakub Kruzik 
153bb9edbfaSJakub Kruzik .seealso: PCDeflationSetType()
154bb9edbfaSJakub Kruzik @*/
155bb9edbfaSJakub Kruzik PetscErrorCode  PCDeflationGetType(PC pc,PCDeflationType *type)
156bb9edbfaSJakub Kruzik {
157bb9edbfaSJakub Kruzik   PetscErrorCode ierr;
158bb9edbfaSJakub Kruzik 
159bb9edbfaSJakub Kruzik   PetscFunctionBegin;
160bb9edbfaSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
161bb9edbfaSJakub Kruzik   ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr);
162bb9edbfaSJakub Kruzik   PetscFunctionReturn(0);
163bb9edbfaSJakub Kruzik }
164bb9edbfaSJakub Kruzik 
165e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
166e662bc50SJakub Kruzik {
167e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
168e662bc50SJakub Kruzik   PetscErrorCode ierr;
169e662bc50SJakub Kruzik 
170e662bc50SJakub Kruzik   PetscFunctionBegin;
171e662bc50SJakub Kruzik   if (transpose) {
172e662bc50SJakub Kruzik     def->Wt = W;
173e662bc50SJakub Kruzik     def->W = NULL;
174e662bc50SJakub Kruzik   } else {
175e662bc50SJakub Kruzik     def->W = W;
176e662bc50SJakub Kruzik   }
177e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
178e662bc50SJakub Kruzik   PetscFunctionReturn(0);
179e662bc50SJakub Kruzik }
180e662bc50SJakub Kruzik 
181e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */
182e662bc50SJakub Kruzik /*@
183e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
184e662bc50SJakub Kruzik 
185e662bc50SJakub Kruzik    Logically Collective on PC
186e662bc50SJakub Kruzik 
187e662bc50SJakub Kruzik    Input Parameters:
188e662bc50SJakub Kruzik +  pc - the preconditioner context
189e662bc50SJakub Kruzik .  W  - deflation matrix
190e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
191e662bc50SJakub Kruzik 
192e662bc50SJakub Kruzik    Level: intermediate
193e662bc50SJakub Kruzik 
194e662bc50SJakub Kruzik .seealso: PCDEFLATION
195e662bc50SJakub Kruzik @*/
196e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
197e662bc50SJakub Kruzik {
198e662bc50SJakub Kruzik   PetscErrorCode ierr;
199e662bc50SJakub Kruzik 
200e662bc50SJakub Kruzik   PetscFunctionBegin;
201e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
202e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
203e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
204e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
205e662bc50SJakub Kruzik   PetscFunctionReturn(0);
206e662bc50SJakub Kruzik }
207e662bc50SJakub Kruzik 
208e662bc50SJakub Kruzik 
20937eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
21037eeb815SJakub Kruzik /*
21137eeb815SJakub Kruzik    PCSetUp_Deflation - Prepares for the use of the Deflation preconditioner
21237eeb815SJakub Kruzik                     by setting data structures and options.
21337eeb815SJakub Kruzik 
21437eeb815SJakub Kruzik    Input Parameter:
21537eeb815SJakub Kruzik .  pc - the preconditioner context
21637eeb815SJakub Kruzik 
21737eeb815SJakub Kruzik    Application Interface Routine: PCSetUp()
21837eeb815SJakub Kruzik 
21937eeb815SJakub Kruzik    Notes:
22037eeb815SJakub Kruzik    The interface routine PCSetUp() is not usually called directly by
22137eeb815SJakub Kruzik    the user, but instead is called by PCApply() if necessary.
22237eeb815SJakub Kruzik */
22337eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
22437eeb815SJakub Kruzik {
225*409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
226*409a357bSJakub Kruzik   KSP              innerksp;
227*409a357bSJakub Kruzik   PC               pcinner;
228*409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
229*409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
230*409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
231*409a357bSJakub Kruzik   MatCompositeType ctype;
232*409a357bSJakub Kruzik   MPI_Comm         comm;
233*409a357bSJakub Kruzik   const char       *prefix;
23437eeb815SJakub Kruzik   PetscErrorCode   ierr;
23537eeb815SJakub Kruzik 
23637eeb815SJakub Kruzik   PetscFunctionBegin;
237*409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
238*409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
239*409a357bSJakub Kruzik   if (def->W || def->Wt) {
240*409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
241*409a357bSJakub Kruzik   } else {
242*409a357bSJakub Kruzik     //ierr = KSPDCGComputeDeflationSpace(ksp);CHKERRQ(ierr);
24337eeb815SJakub Kruzik   }
24437eeb815SJakub Kruzik 
245*409a357bSJakub Kruzik   /* nested deflation */
246*409a357bSJakub Kruzik   if (def->W) {
247*409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
248*409a357bSJakub Kruzik     if (match) {
249*409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
250*409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
25137eeb815SJakub Kruzik     }
252*409a357bSJakub Kruzik   } else {
253*409a357bSJakub Kruzik     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
254*409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
255*409a357bSJakub Kruzik     if (match) {
256*409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
257*409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
258*409a357bSJakub Kruzik     }
259*409a357bSJakub Kruzik     transp = PETSC_TRUE;
260*409a357bSJakub Kruzik   }
261*409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
262*409a357bSJakub Kruzik     ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
263*409a357bSJakub Kruzik     if (!transp) {
264*409a357bSJakub Kruzik       for (i=0; i<size; i++) {
265*409a357bSJakub Kruzik         ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
266*409a357bSJakub Kruzik         //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr);
267*409a357bSJakub Kruzik       }
268*409a357bSJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
269*409a357bSJakub Kruzik         size -= 1;
270*409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
271*409a357bSJakub Kruzik         def->W = mats[size];
272*409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
273*409a357bSJakub Kruzik         if (size > 1) {
274*409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
275*409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
276*409a357bSJakub Kruzik         } else {
277*409a357bSJakub Kruzik           nextDef = mats[0];
278*409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
279*409a357bSJakub Kruzik         }
280*409a357bSJakub Kruzik       } else {
281*409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
282*409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
283*409a357bSJakub Kruzik       }
284*409a357bSJakub Kruzik     }
285*409a357bSJakub Kruzik     ierr = PetscFree(mats);CHKERRQ(ierr);
286*409a357bSJakub Kruzik   }
287*409a357bSJakub Kruzik 
288*409a357bSJakub Kruzik   /* setup coarse problem */
289*409a357bSJakub Kruzik   if (!def->WtAWinv) {
290*409a357bSJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */
291*409a357bSJakub Kruzik     if (!def->WtAW) {
292*409a357bSJakub Kruzik       ierr = PCGetOperators(pc,&Amat,NULL);CHKERRQ(ierr); /* using Amat! */
293*409a357bSJakub Kruzik       /* TODO add implicit product version ? */
294*409a357bSJakub Kruzik       if (!def->AW) {
295*409a357bSJakub Kruzik         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
296*409a357bSJakub Kruzik       } else {
297*409a357bSJakub Kruzik         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
298*409a357bSJakub Kruzik       }
299*409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
300*409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
301*409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
302*409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
303*409a357bSJakub Kruzik       /* Check WtAW is not sigular */
304*409a357bSJakub Kruzik       PetscReal *norms;
305*409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
306*409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
307*409a357bSJakub Kruzik       for (i=0; i<m; i++) {
308*409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
309*409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
310*409a357bSJakub Kruzik         }
311*409a357bSJakub Kruzik       }
312*409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
313*409a357bSJakub Kruzik #endif
314*409a357bSJakub Kruzik     } else {
315*409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
316*409a357bSJakub Kruzik     }
317*409a357bSJakub Kruzik     /* TODO use MATINV */
318*409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
319*409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
320*409a357bSJakub Kruzik     ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
321*409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
322*409a357bSJakub Kruzik     ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
323*409a357bSJakub Kruzik     /* Reduction factor choice */
324*409a357bSJakub Kruzik     red = def->reductionfact;
325*409a357bSJakub Kruzik     if (red < 0) {
326*409a357bSJakub Kruzik       ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
327*409a357bSJakub Kruzik       red  = ceil((float)commsize/ceil((float)m/commsize));
328*409a357bSJakub Kruzik       ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
329*409a357bSJakub Kruzik       if (match) red = commsize;
330*409a357bSJakub Kruzik       ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
331*409a357bSJakub Kruzik     }
332*409a357bSJakub Kruzik     ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
333*409a357bSJakub Kruzik     ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
334*409a357bSJakub Kruzik     ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
335*409a357bSJakub Kruzik     /* Setup KSP and PC */
336*409a357bSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
337*409a357bSJakub Kruzik       /* set default KSPtype */
338*409a357bSJakub Kruzik       if (!def->ksptype) {
339*409a357bSJakub Kruzik         def->ksptype = KSPFGMRES;
340*409a357bSJakub Kruzik         if (flgspd) { /* SPD system */
341*409a357bSJakub Kruzik           def->ksptype = KSPFCG;
342*409a357bSJakub Kruzik         }
343*409a357bSJakub Kruzik       }
344*409a357bSJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
345*409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
346*409a357bSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
347*409a357bSJakub Kruzik       ierr = PCDeflationSetNestLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
348*409a357bSJakub Kruzik       /* inherit options TODO if not set */
349*409a357bSJakub Kruzik       ((PC_Deflation*)(pcinner))->ksptype = def->ksptype;
350*409a357bSJakub Kruzik       ((PC_Deflation*)(pcinner))->correct = def->correct;
351*409a357bSJakub Kruzik       ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv;
352*409a357bSJakub Kruzik       ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst;
353*409a357bSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
354*409a357bSJakub Kruzik     } else { /* the last level */
355*409a357bSJakub Kruzik       ierr = KSPSetType(innerksp,KSPPREONLY);CHKERRQ(ierr);
356*409a357bSJakub Kruzik       /* TODO Cholesky if flgspd? */
357*409a357bSJakub Kruzik       ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
358*409a357bSJakub Kruzik       //TODO remove explicit matSolverPackage
359*409a357bSJakub Kruzik       if (commsize == red) {
360*409a357bSJakub Kruzik         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);CHKERRQ(ierr);
361*409a357bSJakub Kruzik       } else {
362*409a357bSJakub Kruzik         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
363*409a357bSJakub Kruzik       }
364*409a357bSJakub Kruzik     }
365*409a357bSJakub Kruzik     /* TODO use def_[lvl]_ if lvl > 0? */
366*409a357bSJakub Kruzik     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
367*409a357bSJakub Kruzik     if (prefix) {
368*409a357bSJakub Kruzik       ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
369*409a357bSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
370*409a357bSJakub Kruzik     } else {
371*409a357bSJakub Kruzik       ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
372*409a357bSJakub Kruzik     }
373*409a357bSJakub Kruzik     /* TODO: check if WtAWinv is KSP and move following from this if */
374*409a357bSJakub Kruzik     ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
375*409a357bSJakub Kruzik     //if (def->adaptiveconv) {
376*409a357bSJakub Kruzik     //  PetscReal *rnorm;
377*409a357bSJakub Kruzik     //  PetscNew(&rnorm);
378*409a357bSJakub Kruzik     //  ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr);
379*409a357bSJakub Kruzik     //}
380*409a357bSJakub Kruzik     ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
381*409a357bSJakub Kruzik   }
382*409a357bSJakub Kruzik 
383*409a357bSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->work,0,NULL);CHKERRQ(ierr);
38437eeb815SJakub Kruzik   PetscFunctionReturn(0);
38537eeb815SJakub Kruzik }
38637eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
38737eeb815SJakub Kruzik /*
38837eeb815SJakub Kruzik    PCApply_Deflation - Applies the Deflation preconditioner to a vector.
38937eeb815SJakub Kruzik 
39037eeb815SJakub Kruzik    Input Parameters:
39137eeb815SJakub Kruzik .  pc - the preconditioner context
39237eeb815SJakub Kruzik .  x - input vector
39337eeb815SJakub Kruzik 
39437eeb815SJakub Kruzik    Output Parameter:
39537eeb815SJakub Kruzik .  y - output vector
39637eeb815SJakub Kruzik 
39737eeb815SJakub Kruzik    Application Interface Routine: PCApply()
39837eeb815SJakub Kruzik  */
39937eeb815SJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec x,Vec y)
40037eeb815SJakub Kruzik {
40137eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
40237eeb815SJakub Kruzik   PetscErrorCode ierr;
40337eeb815SJakub Kruzik 
40437eeb815SJakub Kruzik   PetscFunctionBegin;
40537eeb815SJakub Kruzik   PetscFunctionReturn(0);
40637eeb815SJakub Kruzik }
40737eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
40837eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
40937eeb815SJakub Kruzik {
41037eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
41137eeb815SJakub Kruzik   PetscErrorCode ierr;
41237eeb815SJakub Kruzik 
41337eeb815SJakub Kruzik   PetscFunctionBegin;
41437eeb815SJakub Kruzik   PetscFunctionReturn(0);
41537eeb815SJakub Kruzik }
41637eeb815SJakub Kruzik 
41737eeb815SJakub Kruzik /*
41837eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
41937eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
42037eeb815SJakub Kruzik 
42137eeb815SJakub Kruzik    Input Parameter:
42237eeb815SJakub Kruzik .  pc - the preconditioner context
42337eeb815SJakub Kruzik 
42437eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
42537eeb815SJakub Kruzik */
42637eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
42737eeb815SJakub Kruzik {
42837eeb815SJakub Kruzik   PetscErrorCode ierr;
42937eeb815SJakub Kruzik 
43037eeb815SJakub Kruzik   PetscFunctionBegin;
43137eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
43237eeb815SJakub Kruzik 
43337eeb815SJakub Kruzik   /*
43437eeb815SJakub Kruzik       Free the private data structure that was hanging off the PC
43537eeb815SJakub Kruzik   */
43637eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
43737eeb815SJakub Kruzik   PetscFunctionReturn(0);
43837eeb815SJakub Kruzik }
43937eeb815SJakub Kruzik 
44037eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
44137eeb815SJakub Kruzik {
44237eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
44337eeb815SJakub Kruzik   PetscErrorCode ierr;
44437eeb815SJakub Kruzik   PetscBool      flg;
44537eeb815SJakub Kruzik   PCDeflationType   deflt,type;
44637eeb815SJakub Kruzik 
44737eeb815SJakub Kruzik   PetscFunctionBegin;
44837eeb815SJakub Kruzik   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
44937eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
45037eeb815SJakub Kruzik   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
45137eeb815SJakub Kruzik   if (flg) {
45237eeb815SJakub Kruzik     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
45337eeb815SJakub Kruzik   }
45437eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
45537eeb815SJakub Kruzik   PetscFunctionReturn(0);
45637eeb815SJakub Kruzik }
45737eeb815SJakub Kruzik 
45837eeb815SJakub Kruzik /*MC
459e361f8b5SJakub Kruzik      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
460e361f8b5SJakub Kruzik      or to a predefined value
46137eeb815SJakub Kruzik 
46237eeb815SJakub Kruzik    Options Database Key:
463e361f8b5SJakub Kruzik +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
46437eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
46537eeb815SJakub Kruzik 
46637eeb815SJakub Kruzik    Level: beginner
46737eeb815SJakub Kruzik 
46837eeb815SJakub Kruzik   Notes:
469e361f8b5SJakub Kruzik     todo
47037eeb815SJakub Kruzik 
47137eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
472e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
47337eeb815SJakub Kruzik M*/
47437eeb815SJakub Kruzik 
47537eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
47637eeb815SJakub Kruzik {
47737eeb815SJakub Kruzik   PC_Deflation   *def;
47837eeb815SJakub Kruzik   PetscErrorCode ierr;
47937eeb815SJakub Kruzik 
48037eeb815SJakub Kruzik   PetscFunctionBegin;
48137eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
48237eeb815SJakub Kruzik   pc->data = (void*)def;
48337eeb815SJakub Kruzik 
484e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
485e662bc50SJakub Kruzik   def->pre           = PETSC_TRUE;
486e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
487e662bc50SJakub Kruzik   def->truenorm      = PETSC_TRUE;
488e662bc50SJakub Kruzik   def->reductionfact = -1;
489e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
490e662bc50SJakub Kruzik   def->spacesize     = 1;
491e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
492e662bc50SJakub Kruzik   def->nestedlvl     = 0;
493e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
494e662bc50SJakub Kruzik   def->adaptiveconv  = PETSC_FALSE;
495e662bc50SJakub Kruzik   def->adaptiveconst = 1.0;
49637eeb815SJakub Kruzik 
49737eeb815SJakub Kruzik   /*
49837eeb815SJakub Kruzik       Set the pointers for the functions that are provided above.
49937eeb815SJakub Kruzik       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
50037eeb815SJakub Kruzik       are called, they will automatically call these functions.  Note we
50137eeb815SJakub Kruzik       choose not to provide a couple of these functions since they are
50237eeb815SJakub Kruzik       not needed.
50337eeb815SJakub Kruzik   */
50437eeb815SJakub Kruzik   pc->ops->apply               = PCApply_Deflation;
50537eeb815SJakub Kruzik   pc->ops->applytranspose      = PCApply_Deflation;
50637eeb815SJakub Kruzik   pc->ops->setup               = PCSetUp_Deflation;
50737eeb815SJakub Kruzik   pc->ops->reset               = PCReset_Deflation;
50837eeb815SJakub Kruzik   pc->ops->destroy             = PCDestroy_Deflation;
50937eeb815SJakub Kruzik   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
51037eeb815SJakub Kruzik   pc->ops->view                = 0;
51137eeb815SJakub Kruzik   pc->ops->applyrichardson     = 0;
51237eeb815SJakub Kruzik   pc->ops->applysymmetricleft  = 0;
51337eeb815SJakub Kruzik   pc->ops->applysymmetricright = 0;
51437eeb815SJakub Kruzik 
51537eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
51637eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
517e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
51837eeb815SJakub Kruzik   PetscFunctionReturn(0);
51937eeb815SJakub Kruzik }
52037eeb815SJakub Kruzik 
521