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