xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 8513960b4b54188a0d25269fe55126c0a8b7a2d0)
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 PCDeflationTypes[]    = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0};
54 
55 const char *const PCDeflationSpaceTypes[] = {
56   "haar",
57   "jacket-haar",
58   "db2",
59   "db4",
60   "db8",
61   "db16",
62   "biorth22",
63   "meyer",
64   "aggregation",
65   "slepc",
66   "slepc-cheap",
67   "user",
68   "DdefSpaceType",
69   "Ddef_SPACE_",
70   0
71 };
72 
73 
74 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
75 {
76   PC_Deflation   *def = (PC_Deflation*)pc->data;
77   PetscErrorCode ierr;
78 
79   PetscFunctionBegin;
80   if (transpose) {
81     def->Wt = W;
82     def->W = NULL;
83   } else {
84     def->W = W;
85   }
86   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 /* TODO create PCDeflationSetSpaceTranspose? */
91 /*@
92    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
93 
94    Logically Collective on PC
95 
96    Input Parameters:
97 +  pc - the preconditioner context
98 .  W  - deflation matrix
99 -  tranpose - indicates that W is an explicit transpose of the deflation matrix
100 
101    Level: intermediate
102 
103 .seealso: PCDEFLATION
104 @*/
105 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
106 {
107   PetscErrorCode ierr;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
112   PetscValidLogicalCollectiveBool(pc,transpose,3);
113   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
118 {
119   PC_Deflation   *def = (PC_Deflation*)pc->data;
120 
121   PetscFunctionBegin;
122   def->nestedlvl = current;
123   def->maxnestedlvl = max;
124   PetscFunctionReturn(0);
125 }
126 
127 /*@
128    PCDeflationSetMaxLvl - Set maximum level of deflation.
129 
130    Logically Collective on PC
131 
132    Input Parameters:
133 +  pc  - the preconditioner context
134 .  max - maximum deflation level
135 
136    Level: intermediate
137 
138 .seealso: PCDEFLATION
139 @*/
140 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
141 {
142   PetscErrorCode ierr;
143 
144   PetscFunctionBegin;
145   PetscValidLogicalCollectiveInt(pc,max,2);
146   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 /*
151   TODO CP corection?
152   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
153 */
154 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
155 {
156   PC_Deflation     *def = (PC_Deflation*)pc->data;
157   Mat              A;
158   Vec              r,w1,w2;
159   PetscBool        nonzero;
160   PetscErrorCode   ierr;
161 
162   PetscFunctionBegin;
163   w1 = def->workcoarse[0];
164   w2 = def->workcoarse[1];
165   r  = def->work;
166   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
167 
168   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
169   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
170   if (nonzero) {
171     ierr = MatMult(A,x,r);CHKERRQ(ierr);               /*    r  <- b - Ax              */
172     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
173   } else {
174     ierr = VecCopy(b,r);CHKERRQ(ierr);                 /*    r  <- b (x is 0)          */
175   }
176 
177   ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
178   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);   /*    w2 <- (W'*A*W)^{-1}*w1    */
179   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);           /*    r  <- W*w2                */
180   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 /*
185   if (def->correct) {
186     z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*r
187   } else {
188     z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*r
189   }
190 */
191 static PetscErrorCode PCApply_Deflation(PC pc,Vec r, Vec z)
192 {
193   PC_Deflation     *def = (PC_Deflation*)pc->data;
194   Mat              A;
195   Vec              u,w1,w2;
196   PetscErrorCode   ierr;
197 
198   PetscFunctionBegin;
199   w1 = def->workcoarse[0];
200   w2 = def->workcoarse[1];
201   u  = def->work;
202   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
203 
204   if (!def->AW) {
205     ierr = MatMult(A,r,u);CHKERRQ(ierr);                       /*    u  <- A*r                 */
206     if (def->correct) ierr = VecAXPY(u,-1.0,r);CHKERRQ(ierr);  /*    u  <- A*r -r              */
207     ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr);        /*    w1 <- W'*u                */
208   } else {
209     ierr = MatMultTranspose(def->AW,u,w1);CHKERRQ(ierr);       /*    u  <- A*r                 */
210     if (def->correct) {
211       ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr);      /*    w2 <- W'*u                */
212       ierr = VecAXPY(w1,-1.0,w2);CHKERRQ(ierr);                /*    w1 <- w1 - w2             */
213     }
214   }
215   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);           /*    w2 <- (W'*A*W)^{-1}*w1    */
216   ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                   /*    u  <- W*w2                */
217   ierr = VecWAXPY(z,-1.0,u,r);CHKERRQ(ierr);                   /*    z  <- r - u               */
218   PetscFunctionReturn(0);
219 }
220 
221 static PetscErrorCode  PCDeflationSetType_Deflation(PC pc,PCDeflationType type)
222 {
223   PC_Deflation *def = (PC_Deflation*)pc->data;
224 
225   PetscFunctionBegin;
226   def->init = PETSC_FALSE;
227   def->pre = PETSC_FALSE;
228   if (type == PC_DEFLATION_POST) {
229     //pc->ops->postsolve = PCPostSolve_Deflation;
230     pc->ops->presolve = 0;
231   } else {
232     pc->ops->presolve = PCPreSolve_Deflation;
233     pc->ops->postsolve = 0;
234     if (type == PC_DEFLATION_INIT) {
235       def->init = PETSC_TRUE;
236       pc->ops->apply = 0;
237     } else {
238       def->pre  = PETSC_TRUE;
239     }
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 /*@
245    PCDeflationSetType - Causes the deflation preconditioner to use only a special
246     initial gues or pre/post solve solution update
247 
248    Logically Collective on PC
249 
250    Input Parameters:
251 +  pc - the preconditioner context
252 -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
253 
254    Options Database Key:
255 .  -pc_deflation_type <pre,init,post>
256 
257    Level: intermediate
258 
259    Concepts: Deflation preconditioner
260 
261 .seealso: PCDeflationGetType()
262 @*/
263 PetscErrorCode  PCDeflationSetType(PC pc,PCDeflationType type)
264 {
265   PetscErrorCode ierr;
266 
267   PetscFunctionBegin;
268   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
269   ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr);
270   PetscFunctionReturn(0);
271 }
272 
273 static PetscErrorCode  PCDeflationGetType_Deflation(PC pc,PCDeflationType *type)
274 {
275   PC_Deflation *def = (PC_Deflation*)pc->data;
276 
277   PetscFunctionBegin;
278   if (def->init) {
279     *type = PC_DEFLATION_INIT;
280   } else if (def->pre) {
281     *type = PC_DEFLATION_PRE;
282   } else {
283     *type = PC_DEFLATION_POST;
284   }
285   PetscFunctionReturn(0);
286 }
287 
288 /*@
289    PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner
290 
291    Not Collective on PC
292 
293    Input Parameter:
294 .  pc - the preconditioner context
295 
296    Output Parameter:
297 -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
298 
299    Level: intermediate
300 
301    Concepts: Deflation preconditioner
302 
303 .seealso: PCDeflationSetType()
304 @*/
305 PetscErrorCode  PCDeflationGetType(PC pc,PCDeflationType *type)
306 {
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
311   ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr);
312   PetscFunctionReturn(0);
313 }
314 
315 static PetscErrorCode PCSetUp_Deflation(PC pc)
316 {
317   PC_Deflation     *def = (PC_Deflation*)pc->data;
318   KSP              innerksp;
319   PC               pcinner;
320   Mat              Amat,nextDef=NULL,*mats;
321   PetscInt         i,m,red,size,commsize;
322   PetscBool        match,flgspd,transp=PETSC_FALSE;
323   MatCompositeType ctype;
324   MPI_Comm         comm;
325   const char       *prefix;
326   PetscErrorCode   ierr;
327 
328   PetscFunctionBegin;
329   if (pc->setupcalled) PetscFunctionReturn(0);
330   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
331   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
332 
333   /* compute a deflation space */
334   if (def->W || def->Wt) {
335     def->spacetype = PC_DEFLATION_SPACE_USER;
336   } else {
337     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
338   }
339 
340   /* nested deflation */
341   if (def->W) {
342     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
343     if (match) {
344       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
345       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
346     }
347   } else {
348     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
349     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
350     if (match) {
351       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
352       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
353     }
354     transp = PETSC_TRUE;
355   }
356   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
357     ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
358     if (!transp) {
359       for (i=0; i<size; i++) {
360         ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
361         //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr);
362       }
363       if (def->nestedlvl < def->maxnestedlvl) {
364         size -= 1;
365         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
366         def->W = mats[size];
367         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
368         if (size > 1) {
369           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
370           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
371         } else {
372           nextDef = mats[0];
373           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
374         }
375       } else {
376         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
377         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
378       }
379     }
380     ierr = PetscFree(mats);CHKERRQ(ierr);
381   }
382 
383   /* setup coarse problem */
384   if (!def->WtAWinv) {
385     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */
386     if (!def->WtAW) {
387       /* TODO add implicit product version ? */
388       if (!def->AW) {
389         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
390       } else {
391         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
392       }
393       /* TODO create MatInheritOption(Mat,MatOption) */
394       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
395       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
396 #if defined(PETSC_USE_DEBUG)
397       /* Check WtAW is not sigular */
398       PetscReal *norms;
399       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
400       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
401       for (i=0; i<m; i++) {
402         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
403           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
404         }
405       }
406       ierr = PetscFree(norms);CHKERRQ(ierr);
407 #endif
408     } else {
409       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
410     }
411     /* TODO use MATINV */
412     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
413     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
414     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
415     ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
416     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
417     ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
418     /* ugly hack to not have overwritten PCTELESCOPE */
419     if (prefix) {
420       ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
421       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
422     } else {
423       ierr = KSPSetOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
424     }
425     ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
426     /* Reduction factor choice */
427     red = def->reductionfact;
428     if (red < 0) {
429       ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
430       red  = ceil((float)commsize/ceil((float)m/commsize));
431       ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
432       if (match) red = commsize;
433       ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
434     }
435     ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
436     ierr = PCSetUp(pcinner);CHKERRQ(ierr);
437     ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
438     if (innerksp) {
439       ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
440       /* Setup KSP and PC */
441       if (nextDef) { /* next level for multilevel deflation */
442         /* set default KSPtype */
443         if (!def->ksptype) {
444           def->ksptype = KSPFGMRES;
445           if (flgspd) { /* SPD system */
446             def->ksptype = KSPFCG;
447           }
448         }
449         ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
450         ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
451         ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
452         ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
453         /* inherit options TODO if not set */
454         ((PC_Deflation*)(pcinner))->ksptype = def->ksptype;
455         ((PC_Deflation*)(pcinner))->correct = def->correct;
456         ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv;
457         ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst;
458         ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
459       } else { /* the last level */
460         ierr = KSPSetType(innerksp,KSPGMRES);CHKERRQ(ierr);
461         /* TODO Cholesky if flgspd? */
462         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
463         //TODO remove explicit matSolverPackage
464         if (commsize == red) {
465           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
466         } else {
467           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
468         }
469       }
470       /* TODO use def_[lvl]_ if lvl > 0? */
471       if (prefix) {
472         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
473         ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
474       } else {
475         ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
476       }
477     }
478     /* TODO: check if WtAWinv is KSP and move following from this if */
479     ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
480     //if (def->adaptiveconv) {
481     //  PetscReal *rnorm;
482     //  PetscNew(&rnorm);
483     //  ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr);
484     //}
485     ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
486   }
487 
488   /* create work vecs */
489   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
490   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 static PetscErrorCode PCReset_Deflation(PC pc)
495 {
496   PC_Deflation      *def = (PC_Deflation*)pc->data;
497   PetscErrorCode    ierr;
498 
499   PetscFunctionBegin;
500   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
501   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
502   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
503   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
504   ierr = MatDestroy(&def->AW);CHKERRQ(ierr);
505   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
506   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 /*
511    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
512    that was created with PCCreate_Deflation().
513 
514    Input Parameter:
515 .  pc - the preconditioner context
516 
517    Application Interface Routine: PCDestroy()
518 */
519 static PetscErrorCode PCDestroy_Deflation(PC pc)
520 {
521   PetscErrorCode ierr;
522 
523   PetscFunctionBegin;
524   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
525   ierr = PetscFree(pc->data);CHKERRQ(ierr);
526   PetscFunctionReturn(0);
527 }
528 
529 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
530 {
531   PC_Deflation      *def = (PC_Deflation*)pc->data;
532   PetscBool         iascii;
533   PetscErrorCode    ierr;
534 
535   PetscFunctionBegin;
536   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
537   if (iascii) {
538     //if (cg->adaptiveconv) {
539     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
540     //}
541     if (def->correct) {
542       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
543     }
544     if (!def->nestedlvl) {
545       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
546       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
547     }
548 
549     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
550     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
551     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
552     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
558 {
559   PC_Deflation      *jac = (PC_Deflation*)pc->data;
560   PetscErrorCode ierr;
561   PetscBool      flg;
562   PCDeflationType   deflt,type;
563 
564   PetscFunctionBegin;
565   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
566   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
567   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
568   if (flg) {
569     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
570   }
571   ierr = PetscOptionsTail();CHKERRQ(ierr);
572   PetscFunctionReturn(0);
573 }
574 
575 /*MC
576      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
577      or to a predefined value
578 
579    Options Database Key:
580 +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
581 -    -pc_jacobi_abs - use the absolute value of the diagonal entry
582 
583    Level: beginner
584 
585   Notes:
586     todo
587 
588 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
589            PCDeflationSetType(), PCDeflationSetSpace()
590 M*/
591 
592 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
593 {
594   PC_Deflation   *def;
595   PetscErrorCode ierr;
596 
597   PetscFunctionBegin;
598   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
599   pc->data = (void*)def;
600 
601   def->init          = PETSC_FALSE;
602   def->pre           = PETSC_TRUE;
603   def->correct       = PETSC_FALSE;
604   def->truenorm      = PETSC_TRUE;
605   def->reductionfact = -1;
606   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
607   def->spacesize     = 1;
608   def->extendsp      = PETSC_FALSE;
609   def->nestedlvl     = 0;
610   def->maxnestedlvl  = 0;
611   def->adaptiveconv  = PETSC_FALSE;
612   def->adaptiveconst = 1.0;
613 
614   /*
615       Set the pointers for the functions that are provided above.
616       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
617       are called, they will automatically call these functions.  Note we
618       choose not to provide a couple of these functions since they are
619       not needed.
620   */
621   pc->ops->apply               = PCApply_Deflation;
622   pc->ops->applytranspose      = PCApply_Deflation;
623   pc->ops->presolve            = PCPreSolve_Deflation;
624   pc->ops->postsolve           = 0;
625   pc->ops->setup               = PCSetUp_Deflation;
626   pc->ops->reset               = PCReset_Deflation;
627   pc->ops->destroy             = PCDestroy_Deflation;
628   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
629   pc->ops->view                = PCView_Deflation;
630   pc->ops->applyrichardson     = 0;
631   pc->ops->applysymmetricleft  = 0;
632   pc->ops->applysymmetricright = 0;
633 
634   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
635   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
636   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
637   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
638   PetscFunctionReturn(0);
639 }
640 
641