xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision e924b0029923c70f7ab6f0e2847aad89b52f3fdf)
1 
2 /*  --------------------------------------------------------------------
3 
4      This file implements a Deflation preconditioner in PETSc as part of PC.
5      You can use this as a starting point for implementing your own
6      preconditioner that is not provided with PETSc. (You might also consider
7      just using PCSHELL)
8 
9      The following basic routines are required for each preconditioner.
10           PCCreate_XXX()          - Creates a preconditioner context
11           PCSetFromOptions_XXX()  - Sets runtime options
12           PCApply_XXX()           - Applies the preconditioner
13           PCDestroy_XXX()         - Destroys the preconditioner context
14      where the suffix "_XXX" denotes a particular implementation, in
15      this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation).
16      These routines are actually called via the common user interface
17      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
18      so the application code interface remains identical for all
19      preconditioners.
20 
21      Another key routine is:
22           PCSetUp_XXX()           - Prepares for the use of a preconditioner
23      by setting data structures and options.   The interface routine PCSetUp()
24      is not usually called directly by the user, but instead is called by
25      PCApply() if necessary.
26 
27      Additional basic routines are:
28           PCView_XXX()            - Prints details of runtime options that
29                                     have actually been used.
30      These are called by application codes via the interface routines
31      PCView().
32 
33      The various types of solvers (preconditioners, Krylov subspace methods,
34      nonlinear solvers, timesteppers) are all organized similarly, so the
35      above description applies to these categories also.  One exception is
36      that the analogues of PCApply() for these components are KSPSolve(),
37      SNESSolve(), and TSSolve().
38 
39      Additional optional functionality unique to preconditioners is left and
40      right symmetric preconditioner application via PCApplySymmetricLeft()
41      and PCApplySymmetricRight().  The Deflation implementation is
42      PCApplySymmetricLeftOrRight_Deflation().
43 
44     -------------------------------------------------------------------- */
45 
46 /*
47    Include files needed for the Deflation preconditioner:
48      pcimpl.h - private include file intended for use by all preconditioners
49 */
50 
51 #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
52 
53 const char *const PCDeflationSpaceTypes[] = {
54   "haar",
55   "jacket-haar",
56   "db2",
57   "db4",
58   "db8",
59   "db16",
60   "biorth22",
61   "meyer",
62   "aggregation",
63   "slepc",
64   "slepc-cheap",
65   "user",
66   "DdefSpaceType",
67   "Ddef_SPACE_",
68   0
69 };
70 
71 
72 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
73 {
74   PC_Deflation   *def = (PC_Deflation*)pc->data;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   if (transpose) {
79     def->Wt = W;
80     def->W = NULL;
81   } else {
82     def->W = W;
83   }
84   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 /* TODO create PCDeflationSetSpaceTranspose? */
89 /*@
90    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
91 
92    Logically Collective on PC
93 
94    Input Parameters:
95 +  pc - the preconditioner context
96 .  W  - deflation matrix
97 -  tranpose - indicates that W is an explicit transpose of the deflation matrix
98 
99    Level: intermediate
100 
101 .seealso: PCDEFLATION
102 @*/
103 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
104 {
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
109   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
110   PetscValidLogicalCollectiveBool(pc,transpose,3);
111   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
112   PetscFunctionReturn(0);
113 }
114 
115 static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
116 {
117   PC_Deflation     *def = (PC_Deflation*)pc->data;
118   PetscErrorCode   ierr;
119 
120   PetscFunctionBegin;
121   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
122   def->WtAW = mat;
123   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
124   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 /*@
129    PCDeflationSetCoarseMat - Set coarse problem Mat.
130 
131    Not Collective
132 
133    Input Parameters:
134 +  pc - preconditioner context
135 -  mat - coarse problem mat
136 
137    Level: developer
138 
139 .seealso: PCDEFLATION
140 @*/
141 PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
142 {
143   PetscErrorCode ierr;
144 
145   PetscFunctionBegin;
146   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
147   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
148   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 
152 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
153 {
154   PC_Deflation   *def = (PC_Deflation*)pc->data;
155 
156   PetscFunctionBegin;
157   if (current) def->nestedlvl = current;
158   def->maxnestedlvl = max;
159   PetscFunctionReturn(0);
160 }
161 
162 /*@
163    PCDeflationSetMaxLvl - Set maximum level of deflation.
164 
165    Logically Collective on PC
166 
167    Input Parameters:
168 +  pc  - the preconditioner context
169 -  max - maximum deflation level
170 
171    Level: intermediate
172 
173 .seealso: PCDEFLATION
174 @*/
175 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
176 {
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181   PetscValidLogicalCollectiveInt(pc,max,2);
182   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
187 {
188   PC_Deflation   *def = (PC_Deflation*)pc->data;
189 
190   PetscFunctionBegin;
191   *apc = def->pc;
192   PetscFunctionReturn(0);
193 }
194 
195 /*@
196    PCDeflationGetPC - Returns a pointer to additional preconditioner.
197 
198    Not Collective
199 
200    Input Parameters:
201 .  pc  - the preconditioner context
202 
203    Output Parameter:
204 .  apc - additional preconditioner
205 
206    Level: advanced
207 
208 .seealso: PCDeflationSetPC(), PCDEFLATION
209 @*/
210 PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
211 {
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
216   PetscValidPointer(pc,2);
217   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 
221 static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
222 {
223   PC_Deflation   *def = (PC_Deflation*)pc->data;
224   PetscErrorCode ierr;
225 
226   PetscFunctionBegin;
227   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
228   def->pc = apc;
229   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
230   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
231   PetscFunctionReturn(0);
232 }
233 
234 /*@
235    PCDeflationSetPC - Set additional preconditioner.
236 
237    Collective on PC
238 
239    Input Parameters:
240 +  pc  - the preconditioner context
241 -  apc - additional preconditioner
242 
243    Level: developer
244 
245 .seealso: PCDeflationGetPC(), PCDEFLATION
246 @*/
247 PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
248 {
249   PetscErrorCode ierr;
250 
251   PetscFunctionBegin;
252   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
253   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
254   PetscCheckSameComm(pc,1,apc,2);
255   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
260 {
261   PC_Deflation     *def = (PC_Deflation*)pc->data;
262 
263   PetscFunctionBegin;
264   *ksp = def->WtAWinv;
265   PetscFunctionReturn(0);
266 }
267 
268 /*@
269    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
270 
271    Not Collective
272 
273    Input Parameters:
274 .  pc - preconditioner context
275 
276    Output Parameter:
277 .  ksp - coarse problem KSP context
278 
279    Level: developer
280 
281 .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION
282 @*/
283 PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
284 {
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
289   PetscValidPointer(ksp,2);
290   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp)
295 {
296   PC_Deflation     *def = (PC_Deflation*)pc->data;
297   PetscErrorCode   ierr;
298 
299   PetscFunctionBegin;
300   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
301   def->WtAWinv = ksp;
302   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
303   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 /*@
308    PCDeflationSetCoarseKSP - Set coarse problem KSP.
309 
310    Collective on PC
311 
312    Input Parameters:
313 +  pc - preconditioner context
314 -  ksp - coarse problem KSP context
315 
316    Level: developer
317 
318 .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION
319 @*/
320 PetscErrorCode  PCDeflationSetCoarseKSP(PC pc,KSP ksp)
321 {
322   PetscErrorCode ierr;
323 
324   PetscFunctionBegin;
325   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
326   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
327   PetscCheckSameComm(pc,1,ksp,2);
328   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr);
329   PetscFunctionReturn(0);
330 }
331 
332 /*
333   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
334 */
335 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
336 {
337   PC_Deflation     *def = (PC_Deflation*)pc->data;
338   Mat              A;
339   Vec              r,w1,w2;
340   PetscBool        nonzero;
341   PetscErrorCode   ierr;
342 
343   PetscFunctionBegin;
344   w1 = def->workcoarse[0];
345   w2 = def->workcoarse[1];
346   r  = def->work;
347   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
348 
349   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
350   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
351   if (nonzero) {
352     ierr = MatMult(A,x,r);CHKERRQ(ierr);               /*    r  <- b - Ax              */
353     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
354   } else {
355     ierr = VecCopy(b,r);CHKERRQ(ierr);                 /*    r  <- b (x is 0)          */
356   }
357 
358   ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
359   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);   /*    w2 <- (W'*A*W)^{-1}*w1    */
360   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);           /*    r  <- W*w2                */
361   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 /*
366   if (def->correct) {
367     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
368   } else {
369     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
370   }
371 */
372 static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
373 {
374   PC_Deflation     *def = (PC_Deflation*)pc->data;
375   Mat              A;
376   Vec              u,w1,w2;
377   PetscErrorCode   ierr;
378 
379   PetscFunctionBegin;
380   w1 = def->workcoarse[0];
381   w2 = def->workcoarse[1];
382   u  = def->work;
383   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
384 
385   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                      /*    z <- M^{-1}*r             */
386   if (!def->init) {
387     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                  /*    w1 <- W'*A*z              */
388     if (def->correct) {
389       if (def->Wt) {
390         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);               /*    w2 <- W'*r                */
391       } else {
392         ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr);       /*    w2 <- W'*r                */
393       }
394       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);  /*    w1 <- w1 - l*w2           */
395     }
396     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);            /*    w2 <- (W'*A*W)^{-1}*w1    */
397     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                    /*    u  <- W*w2                */
398     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                       /*    z  <- z - u               */
399   }
400   PetscFunctionReturn(0);
401 }
402 
403 static PetscErrorCode PCSetUp_Deflation(PC pc)
404 {
405   PC_Deflation     *def = (PC_Deflation*)pc->data;
406   KSP              innerksp;
407   PC               pcinner;
408   Mat              Amat,nextDef=NULL,*mats;
409   PetscInt         i,m,red,size,commsize;
410   PetscBool        match,flgspd,transp=PETSC_FALSE;
411   MatCompositeType ctype;
412   MPI_Comm         comm;
413   const char       *prefix;
414   PetscErrorCode   ierr;
415 
416   PetscFunctionBegin;
417   if (pc->setupcalled) PetscFunctionReturn(0);
418   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
419   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
420 
421   /* compute a deflation space */
422   if (def->W || def->Wt) {
423     def->spacetype = PC_DEFLATION_SPACE_USER;
424   } else {
425     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
426   }
427 
428   /* nested deflation */
429   if (def->W) {
430     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
431     if (match) {
432       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
433       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
434     }
435   } else {
436     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
437     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
438     if (match) {
439       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
440       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
441     }
442     transp = PETSC_TRUE;
443   }
444   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
445     if (!transp) {
446       if (def->nestedlvl < def->maxnestedlvl) {
447         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
448         for (i=0; i<size; i++) {
449           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
450         }
451         size -= 1;
452         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
453         def->W = mats[size];
454         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
455         if (size > 1) {
456           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
457           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
458         } else {
459           nextDef = mats[0];
460           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
461         }
462         ierr = PetscFree(mats);CHKERRQ(ierr);
463       } else {
464         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
465         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
466       }
467     } else {
468       if (def->nestedlvl < def->maxnestedlvl) {
469         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
470         for (i=0; i<size; i++) {
471           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
472         }
473         size -= 1;
474         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
475         def->Wt = mats[0];
476         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
477         if (size > 1) {
478           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
479           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
480         } else {
481           nextDef = mats[1];
482           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
483         }
484         ierr = PetscFree(mats);CHKERRQ(ierr);
485       } else {
486         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
487         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
488       }
489     }
490   }
491 
492   if (transp) {
493     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
494     ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
495   }
496 
497 
498   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
499 
500   /* assemble WtA */
501   if (!def->WtA) {
502     if (def->Wt) {
503       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
504     } else {
505       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
506     }
507   }
508   /* setup coarse problem */
509   if (!def->WtAWinv) {
510     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
511     if (!def->WtAW) {
512       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
513       /* TODO create MatInheritOption(Mat,MatOption) */
514       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
515       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
516 #if defined(PETSC_USE_DEBUG)
517       /* Check columns of W are not in kernel of A */
518       PetscReal *norms;
519       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
520       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
521       for (i=0; i<m; i++) {
522         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
523           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
524         }
525       }
526       ierr = PetscFree(norms);CHKERRQ(ierr);
527 #endif
528     } else {
529       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
530     }
531     /* TODO use MATINV */
532     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
533     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
534     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
535     /* Setup KSP and PC */
536     if (nextDef) { /* next level for multilevel deflation */
537       innerksp = def->WtAWinv;
538       /* set default KSPtype */
539       if (!def->ksptype) {
540         def->ksptype = KSPFGMRES;
541         if (flgspd) { /* SPD system */
542           def->ksptype = KSPFCG;
543         }
544       }
545       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
546       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
547       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
548       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
549       /* inherit options */
550       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
551       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
552       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
553     } else { /* the last level */
554       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
555       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
556       /* ugly hack to not have overwritten PCTELESCOPE */
557       if (prefix) {
558         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
559       }
560       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
561       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
562       /* Reduction factor choice */
563       red = def->reductionfact;
564       if (red < 0) {
565         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
566         red  = ceil((float)commsize/ceil((float)m/commsize));
567         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
568         if (match) red = commsize;
569         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
570       }
571       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
572       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
573       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
574       if (innerksp) {
575         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
576         /* TODO Cholesky if flgspd? */
577         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
578         //TODO remove explicit matSolverPackage
579         if (commsize == red) {
580           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
581         } else {
582           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
583         }
584       }
585     }
586 
587     if (innerksp) {
588       /* TODO use def_[lvl]_ if lvl > 0? */
589       if (prefix) {
590         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
591       }
592       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
593       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
594       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
595     }
596   }
597   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
598   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
599 
600   if (!def->pc) {
601     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
602     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
603     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
604     if (prefix) {
605       ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
606     }
607     ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr);
608     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
609     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
610   }
611 
612   /* create work vecs */
613   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
614   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
615   PetscFunctionReturn(0);
616 }
617 
618 static PetscErrorCode PCReset_Deflation(PC pc)
619 {
620   PC_Deflation      *def = (PC_Deflation*)pc->data;
621   PetscErrorCode    ierr;
622 
623   PetscFunctionBegin;
624   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
625   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
626   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
627   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
628   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
629   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
630   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
631   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
632   PetscFunctionReturn(0);
633 }
634 
635 /*
636    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
637    that was created with PCCreate_Deflation().
638 
639    Input Parameter:
640 .  pc - the preconditioner context
641 
642    Application Interface Routine: PCDestroy()
643 */
644 static PetscErrorCode PCDestroy_Deflation(PC pc)
645 {
646   PetscErrorCode ierr;
647 
648   PetscFunctionBegin;
649   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
650   ierr = PetscFree(pc->data);CHKERRQ(ierr);
651   PetscFunctionReturn(0);
652 }
653 
654 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
655 {
656   PC_Deflation      *def = (PC_Deflation*)pc->data;
657   PetscBool         iascii;
658   PetscErrorCode    ierr;
659 
660   PetscFunctionBegin;
661   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
662   if (iascii) {
663     //if (cg->adaptiveconv) {
664     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
665     //}
666     if (def->correct) {
667       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
668     }
669     if (!def->nestedlvl) {
670       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
671       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
672     }
673 
674     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr);
675     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
676     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
677     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
678 
679     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
680     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
681     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
682     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
683   }
684   PetscFunctionReturn(0);
685 }
686 
687 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
688 {
689   PC_Deflation      *def = (PC_Deflation*)pc->data;
690   PetscErrorCode    ierr;
691 
692   PetscFunctionBegin;
693   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
694   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
695   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
696   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
697 //TODO add set function and fix manpages
698   ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr);
699   ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr);
700   ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
701   ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
702   ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
703   ierr = PetscOptionsTail();CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 /*MC
708      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
709      or to a predefined value
710 
711    Options Database Key:
712 +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
713 -    -pc_jacobi_abs - use the absolute value of the diagonal entry
714 
715    Level: beginner
716 
717   Notes:
718     todo
719 
720 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
721            PCDeflationSetType(), PCDeflationSetSpace()
722 M*/
723 
724 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
725 {
726   PC_Deflation   *def;
727   PetscErrorCode ierr;
728 
729   PetscFunctionBegin;
730   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
731   pc->data = (void*)def;
732 
733   def->init          = PETSC_FALSE;
734   def->correct       = PETSC_FALSE;
735   def->correctfact   = 1.0;
736   def->reductionfact = -1;
737   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
738   def->spacesize     = 1;
739   def->extendsp      = PETSC_FALSE;
740   def->nestedlvl     = 0;
741   def->maxnestedlvl  = 0;
742 
743   /*
744       Set the pointers for the functions that are provided above.
745       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
746       are called, they will automatically call these functions.  Note we
747       choose not to provide a couple of these functions since they are
748       not needed.
749   */
750   pc->ops->apply               = PCApply_Deflation;
751   pc->ops->applytranspose      = PCApply_Deflation;
752   pc->ops->presolve            = PCPreSolve_Deflation;
753   pc->ops->postsolve           = 0;
754   pc->ops->setup               = PCSetUp_Deflation;
755   pc->ops->reset               = PCReset_Deflation;
756   pc->ops->destroy             = PCDestroy_Deflation;
757   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
758   pc->ops->view                = PCView_Deflation;
759   pc->ops->applyrichardson     = 0;
760   pc->ops->applysymmetricleft  = 0;
761   pc->ops->applysymmetricright = 0;
762 
763   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
764   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
765   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
766   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
767   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
768   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
769   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
770   PetscFunctionReturn(0);
771 }
772 
773