xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision f3eaa4f06f3a57e68cf6844593c2de26424fb800)
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 PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
116 {
117   PC_Deflation   *def = (PC_Deflation*)pc->data;
118 
119   PetscFunctionBegin;
120   if (current) def->nestedlvl = current;
121   def->maxnestedlvl = max;
122   PetscFunctionReturn(0);
123 }
124 
125 /*@
126    PCDeflationSetMaxLvl - Set maximum level of deflation.
127 
128    Logically Collective on PC
129 
130    Input Parameters:
131 +  pc  - the preconditioner context
132 -  max - maximum deflation level
133 
134    Level: intermediate
135 
136 .seealso: PCDEFLATION
137 @*/
138 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
139 {
140   PetscErrorCode ierr;
141 
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
144   PetscValidLogicalCollectiveInt(pc,max,2);
145   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 }
148 
149 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
150 {
151   PC_Deflation   *def = (PC_Deflation*)pc->data;
152 
153   PetscFunctionBegin;
154   *apc = def->pc;
155   PetscFunctionReturn(0);
156 }
157 
158 /*@
159    PCDeflationGetPC - Returns a pointer to additional preconditioner.
160 
161    Not Collective
162 
163    Input Parameters:
164 .  pc  - the preconditioner context
165 
166    Output Parameter:
167 .  apc - additional preconditioner
168 
169    Level: advanced
170 
171 .seealso: PCDeflationSetPC(), PCDEFLATION
172 @*/
173 PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
174 {
175   PetscErrorCode ierr;
176 
177   PetscFunctionBegin;
178   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
179   PetscValidPointer(pc,2);
180   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
185 {
186   PC_Deflation   *def = (PC_Deflation*)pc->data;
187   PetscErrorCode ierr;
188 
189   PetscFunctionBegin;
190   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
191   def->pc = apc;
192   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
193   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 /*@
198    PCDeflationSetPC - Set additional preconditioner.
199 
200    Collective on PC
201 
202    Input Parameters:
203 +  pc  - the preconditioner context
204 -  apc - additional preconditioner
205 
206    Level: developer
207 
208 .seealso: PCDeflationGetPC(), PCDEFLATION
209 @*/
210 PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
211 {
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
216   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
217   PetscCheckSameComm(pc,1,apc,2);
218   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
223 {
224   PC_Deflation     *def = (PC_Deflation*)pc->data;
225 
226   PetscFunctionBegin;
227   *ksp = def->WtAWinv;
228   PetscFunctionReturn(0);
229 }
230 
231 /*@
232    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
233 
234    Not Collective
235 
236    Input Parameters:
237 .  pc - preconditioner context
238 
239    Output Parameter:
240 .  ksp - coarse problem KSP context
241 
242    Level: developer
243 
244 .seealso: PCDeflationSetCoarseKSP(), PCDEFLATION
245 @*/
246 PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
247 {
248   PetscErrorCode ierr;
249 
250   PetscFunctionBegin;
251   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
252   PetscValidPointer(ksp,2);
253   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
254   PetscFunctionReturn(0);
255 }
256 
257 static PetscErrorCode PCDeflationSetCoarseKSP_Deflation(PC pc,KSP ksp)
258 {
259   PC_Deflation     *def = (PC_Deflation*)pc->data;
260   PetscErrorCode   ierr;
261 
262   PetscFunctionBegin;
263   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
264   def->WtAWinv = ksp;
265   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
266   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAWinv);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 /*@
271    PCDeflationSetCoarseKSP - Set coarse problem KSP.
272 
273    Collective on PC
274 
275    Input Parameters:
276 +  pc - preconditioner context
277 -  ksp - coarse problem KSP context
278 
279    Level: developer
280 
281 .seealso: PCDeflationGetCoarseKSP(), PCDEFLATION
282 @*/
283 PetscErrorCode  PCDeflationSetCoarseKSP(PC pc,KSP ksp)
284 {
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
289   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
290   PetscCheckSameComm(pc,1,ksp,2);
291   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseKSP_C",(PC,KSP),(pc,ksp));CHKERRQ(ierr);
292   PetscFunctionReturn(0);
293 }
294 
295 /*
296   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
297 */
298 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
299 {
300   PC_Deflation     *def = (PC_Deflation*)pc->data;
301   Mat              A;
302   Vec              r,w1,w2;
303   PetscBool        nonzero;
304   PetscErrorCode   ierr;
305 
306   PetscFunctionBegin;
307   w1 = def->workcoarse[0];
308   w2 = def->workcoarse[1];
309   r  = def->work;
310   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
311 
312   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
313   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
314   if (nonzero) {
315     ierr = MatMult(A,x,r);CHKERRQ(ierr);               /*    r  <- b - Ax              */
316     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
317   } else {
318     ierr = VecCopy(b,r);CHKERRQ(ierr);                 /*    r  <- b (x is 0)          */
319   }
320 
321   ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
322   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);   /*    w2 <- (W'*A*W)^{-1}*w1    */
323   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);           /*    r  <- W*w2                */
324   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
325   PetscFunctionReturn(0);
326 }
327 
328 /*
329   if (def->correct) {
330     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
331   } else {
332     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
333   }
334 */
335 static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
336 {
337   PC_Deflation     *def = (PC_Deflation*)pc->data;
338   Mat              A;
339   Vec              u,w1,w2;
340   PetscErrorCode   ierr;
341 
342   PetscFunctionBegin;
343   w1 = def->workcoarse[0];
344   w2 = def->workcoarse[1];
345   u  = def->work;
346   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
347 
348   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                      /*    z <- M^{-1}*r             */
349   if (!def->init) {
350     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                  /*    w1 <- W'*A*z              */
351     if (def->correct) {
352       if (def->Wt) {
353         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);               /*    w2 <- W'*r                */
354       } else {
355         ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr);       /*    w2 <- W'*r                */
356       }
357       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);  /*    w1 <- w1 - l*w2           */
358     }
359     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);            /*    w2 <- (W'*A*W)^{-1}*w1    */
360     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                    /*    u  <- W*w2                */
361     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                       /*    z  <- z - u               */
362   }
363   PetscFunctionReturn(0);
364 }
365 
366 static PetscErrorCode PCSetUp_Deflation(PC pc)
367 {
368   PC_Deflation     *def = (PC_Deflation*)pc->data;
369   KSP              innerksp;
370   PC               pcinner;
371   Mat              Amat,nextDef=NULL,*mats;
372   PetscInt         i,m,red,size,commsize;
373   PetscBool        match,flgspd,transp=PETSC_FALSE;
374   MatCompositeType ctype;
375   MPI_Comm         comm;
376   const char       *prefix;
377   PetscErrorCode   ierr;
378 
379   PetscFunctionBegin;
380   if (pc->setupcalled) PetscFunctionReturn(0);
381   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
382   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
383 
384   /* compute a deflation space */
385   if (def->W || def->Wt) {
386     def->spacetype = PC_DEFLATION_SPACE_USER;
387   } else {
388     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
389   }
390 
391   /* nested deflation */
392   if (def->W) {
393     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
394     if (match) {
395       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
396       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
397     }
398   } else {
399     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
400     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
401     if (match) {
402       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
403       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
404     }
405     transp = PETSC_TRUE;
406   }
407   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
408     if (!transp) {
409       if (def->nestedlvl < def->maxnestedlvl) {
410         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
411         for (i=0; i<size; i++) {
412           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
413         }
414         size -= 1;
415         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
416         def->W = mats[size];
417         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
418         if (size > 1) {
419           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
420           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
421         } else {
422           nextDef = mats[0];
423           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
424         }
425         ierr = PetscFree(mats);CHKERRQ(ierr);
426       } else {
427         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
428         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
429       }
430     } else {
431       if (def->nestedlvl < def->maxnestedlvl) {
432         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
433         for (i=0; i<size; i++) {
434           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
435         }
436         size -= 1;
437         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
438         def->Wt = mats[0];
439         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
440         if (size > 1) {
441           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
442           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
443         } else {
444           nextDef = mats[1];
445           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
446         }
447         ierr = PetscFree(mats);CHKERRQ(ierr);
448       } else {
449         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
450         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
451       }
452     }
453   }
454 
455   if (transp) {
456     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
457     ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
458   }
459 
460 
461   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
462 
463   /* assemble WtA */
464   if (!def->WtA) {
465     if (def->Wt) {
466       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
467     } else {
468       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
469     }
470   }
471   /* setup coarse problem */
472   if (!def->WtAWinv) {
473     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
474     if (!def->WtAW) {
475       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
476       /* TODO create MatInheritOption(Mat,MatOption) */
477       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
478       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
479 #if defined(PETSC_USE_DEBUG)
480       /* Check columns of W are not in kernel of A */
481       PetscReal *norms;
482       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
483       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
484       for (i=0; i<m; i++) {
485         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
486           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
487         }
488       }
489       ierr = PetscFree(norms);CHKERRQ(ierr);
490 #endif
491     } else {
492       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
493     }
494     /* TODO use MATINV */
495     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
496     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
497     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
498     /* Setup KSP and PC */
499     if (nextDef) { /* next level for multilevel deflation */
500       innerksp = def->WtAWinv;
501       /* set default KSPtype */
502       if (!def->ksptype) {
503         def->ksptype = KSPFGMRES;
504         if (flgspd) { /* SPD system */
505           def->ksptype = KSPFCG;
506         }
507       }
508       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
509       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
510       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
511       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
512       /* inherit options */
513       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
514       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
515       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
516     } else { /* the last level */
517       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
518       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
519       /* ugly hack to not have overwritten PCTELESCOPE */
520       if (prefix) {
521         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
522       }
523       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
524       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
525       /* Reduction factor choice */
526       red = def->reductionfact;
527       if (red < 0) {
528         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
529         red  = ceil((float)commsize/ceil((float)m/commsize));
530         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
531         if (match) red = commsize;
532         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
533       }
534       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
535       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
536       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
537       if (innerksp) {
538         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
539         /* TODO Cholesky if flgspd? */
540         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
541         //TODO remove explicit matSolverPackage
542         if (commsize == red) {
543           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
544         } else {
545           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
546         }
547       }
548     }
549 
550     if (innerksp) {
551       /* TODO use def_[lvl]_ if lvl > 0? */
552       if (prefix) {
553         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
554       }
555       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
556       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
557       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
558     }
559   }
560   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
561   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
562 
563   if (!def->pc) {
564     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
565     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
566     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
567     if (prefix) {
568       ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
569     }
570     ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr);
571     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
572     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
573   }
574 
575   /* create work vecs */
576   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
577   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580 
581 static PetscErrorCode PCReset_Deflation(PC pc)
582 {
583   PC_Deflation      *def = (PC_Deflation*)pc->data;
584   PetscErrorCode    ierr;
585 
586   PetscFunctionBegin;
587   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
588   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
589   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
590   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
591   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
592   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
593   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
594   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 /*
599    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
600    that was created with PCCreate_Deflation().
601 
602    Input Parameter:
603 .  pc - the preconditioner context
604 
605    Application Interface Routine: PCDestroy()
606 */
607 static PetscErrorCode PCDestroy_Deflation(PC pc)
608 {
609   PetscErrorCode ierr;
610 
611   PetscFunctionBegin;
612   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
613   ierr = PetscFree(pc->data);CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 
617 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
618 {
619   PC_Deflation      *def = (PC_Deflation*)pc->data;
620   PetscBool         iascii;
621   PetscErrorCode    ierr;
622 
623   PetscFunctionBegin;
624   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
625   if (iascii) {
626     //if (cg->adaptiveconv) {
627     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
628     //}
629     if (def->correct) {
630       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
631     }
632     if (!def->nestedlvl) {
633       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
634       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
635     }
636 
637     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr);
638     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
639     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
640     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
641 
642     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
643     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
644     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
645     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
646   }
647   PetscFunctionReturn(0);
648 }
649 
650 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
651 {
652   PC_Deflation      *def = (PC_Deflation*)pc->data;
653   PetscErrorCode    ierr;
654 
655   PetscFunctionBegin;
656   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
657   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
658   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
659   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
660 //TODO add set function and fix manpages
661   ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr);
662   ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr);
663   ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
664   ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
665   ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
666   ierr = PetscOptionsTail();CHKERRQ(ierr);
667   PetscFunctionReturn(0);
668 }
669 
670 /*MC
671      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
672      or to a predefined value
673 
674    Options Database Key:
675 +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
676 -    -pc_jacobi_abs - use the absolute value of the diagonal entry
677 
678    Level: beginner
679 
680   Notes:
681     todo
682 
683 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
684            PCDeflationSetType(), PCDeflationSetSpace()
685 M*/
686 
687 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
688 {
689   PC_Deflation   *def;
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
694   pc->data = (void*)def;
695 
696   def->init          = PETSC_FALSE;
697   def->correct       = PETSC_FALSE;
698   def->correctfact   = 1.0;
699   def->reductionfact = -1;
700   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
701   def->spacesize     = 1;
702   def->extendsp      = PETSC_FALSE;
703   def->nestedlvl     = 0;
704   def->maxnestedlvl  = 0;
705 
706   /*
707       Set the pointers for the functions that are provided above.
708       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
709       are called, they will automatically call these functions.  Note we
710       choose not to provide a couple of these functions since they are
711       not needed.
712   */
713   pc->ops->apply               = PCApply_Deflation;
714   pc->ops->applytranspose      = PCApply_Deflation;
715   pc->ops->presolve            = PCPreSolve_Deflation;
716   pc->ops->postsolve           = 0;
717   pc->ops->setup               = PCSetUp_Deflation;
718   pc->ops->reset               = PCReset_Deflation;
719   pc->ops->destroy             = PCDestroy_Deflation;
720   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
721   pc->ops->view                = PCView_Deflation;
722   pc->ops->applyrichardson     = 0;
723   pc->ops->applysymmetricleft  = 0;
724   pc->ops->applysymmetricright = 0;
725 
726   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
727   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
728   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
729   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
730   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
731   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
732   PetscFunctionReturn(0);
733 }
734 
735