xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 409a357bae74e46200d26aa0a2f21a6fc06036bf)
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 <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
52 
53 const char *const PCDeflationTypes[]    = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0};
54 
55 /*
56    Private context (data structure) for the deflation preconditioner.
57 */
58 typedef struct {
59   PetscBool init;            /* do only init step - error correction of direction is omitted */
60   PetscBool pre;             /* start with x0 being the solution in the deflation space */
61   PetscBool correct;         /* add CP (Qr) correction to descent direction */
62   PetscBool truenorm;
63   PetscBool adaptiveconv;
64   PetscReal adaptiveconst;
65   PetscInt  reductionfact;
66   Mat       W,Wt,AW,WtAW;    /* deflation space, coarse problem mats */
67   KSP       WtAWinv;         /* deflation coarse problem */
68   KSPType   ksptype;
69   Vec       *work;
70 
71   PCDeflationSpaceType spacetype;
72   PetscInt             spacesize;
73   PetscInt             nestedlvl;
74   PetscInt             maxnestedlvl;
75   PetscBool            extendsp;
76 } PC_Deflation;
77 
78 static PetscErrorCode  PCDeflationSetType_Deflation(PC pc,PCDeflationType type)
79 {
80   PC_Deflation *def = (PC_Deflation*)pc->data;
81 
82   PetscFunctionBegin;
83   def->init = PETSC_FALSE;
84   def->pre = PETSC_FALSE;
85   if (type == PC_DEFLATION_INIT) {
86     def->init = PETSC_TRUE;
87     def->pre  = PETSC_TRUE;
88   } else if (type == PC_DEFLATION_PRE) {
89     def->pre  = PETSC_TRUE;
90   }
91   PetscFunctionReturn(0);
92 }
93 
94 /*@
95    PCDeflationSetType - Causes the deflation preconditioner to use only a special
96     initial gues or pre/post solve solution update
97 
98    Logically Collective on PC
99 
100    Input Parameters:
101 +  pc - the preconditioner context
102 -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
103 
104    Options Database Key:
105 .  -pc_deflation_type <pre,init,post>
106 
107    Level: intermediate
108 
109    Concepts: Deflation preconditioner
110 
111 .seealso: PCDeflationGetType()
112 @*/
113 PetscErrorCode  PCDeflationSetType(PC pc,PCDeflationType type)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119   ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode  PCDeflationGetType_Deflation(PC pc,PCDeflationType *type)
124 {
125   PC_Deflation *def = (PC_Deflation*)pc->data;
126 
127   PetscFunctionBegin;
128   if (def->init) {
129     *type = PC_DEFLATION_INIT;
130   } else if (def->pre) {
131     *type = PC_DEFLATION_PRE;
132   } else {
133     *type = PC_DEFLATION_POST;
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 /*@
139    PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner
140 
141    Not Collective on PC
142 
143    Input Parameter:
144 .  pc - the preconditioner context
145 
146    Output Parameter:
147 -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
148 
149    Level: intermediate
150 
151    Concepts: Deflation preconditioner
152 
153 .seealso: PCDeflationSetType()
154 @*/
155 PetscErrorCode  PCDeflationGetType(PC pc,PCDeflationType *type)
156 {
157   PetscErrorCode ierr;
158 
159   PetscFunctionBegin;
160   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
161   ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr);
162   PetscFunctionReturn(0);
163 }
164 
165 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
166 {
167   PC_Deflation   *def = (PC_Deflation*)pc->data;
168   PetscErrorCode ierr;
169 
170   PetscFunctionBegin;
171   if (transpose) {
172     def->Wt = W;
173     def->W = NULL;
174   } else {
175     def->W = W;
176   }
177   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 /* TODO create PCDeflationSetSpaceTranspose? */
182 /*@
183    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
184 
185    Logically Collective on PC
186 
187    Input Parameters:
188 +  pc - the preconditioner context
189 .  W  - deflation matrix
190 -  tranpose - indicates that W is an explicit transpose of the deflation matrix
191 
192    Level: intermediate
193 
194 .seealso: PCDEFLATION
195 @*/
196 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
197 {
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
202   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
203   PetscValidLogicalCollectiveBool(pc,transpose,3);
204   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207 
208 
209 /* -------------------------------------------------------------------------- */
210 /*
211    PCSetUp_Deflation - Prepares for the use of the Deflation preconditioner
212                     by setting data structures and options.
213 
214    Input Parameter:
215 .  pc - the preconditioner context
216 
217    Application Interface Routine: PCSetUp()
218 
219    Notes:
220    The interface routine PCSetUp() is not usually called directly by
221    the user, but instead is called by PCApply() if necessary.
222 */
223 static PetscErrorCode PCSetUp_Deflation(PC pc)
224 {
225   PC_Deflation     *def = (PC_Deflation*)pc->data;
226   KSP              innerksp;
227   PC               pcinner;
228   Mat              Amat,nextDef=NULL,*mats;
229   PetscInt         i,m,red,size,commsize;
230   PetscBool        match,flgspd,transp=PETSC_FALSE;
231   MatCompositeType ctype;
232   MPI_Comm         comm;
233   const char       *prefix;
234   PetscErrorCode   ierr;
235 
236   PetscFunctionBegin;
237   if (pc->setupcalled) PetscFunctionReturn(0);
238   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
239   if (def->W || def->Wt) {
240     def->spacetype = PC_DEFLATION_SPACE_USER;
241   } else {
242     //ierr = KSPDCGComputeDeflationSpace(ksp);CHKERRQ(ierr);
243   }
244 
245   /* nested deflation */
246   if (def->W) {
247     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
248     if (match) {
249       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
250       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
251     }
252   } else {
253     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
254     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
255     if (match) {
256       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
257       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
258     }
259     transp = PETSC_TRUE;
260   }
261   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
262     ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
263     if (!transp) {
264       for (i=0; i<size; i++) {
265         ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
266         //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr);
267       }
268       if (def->nestedlvl < def->maxnestedlvl) {
269         size -= 1;
270         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
271         def->W = mats[size];
272         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
273         if (size > 1) {
274           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
275           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
276         } else {
277           nextDef = mats[0];
278           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
279         }
280       } else {
281         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
282         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
283       }
284     }
285     ierr = PetscFree(mats);CHKERRQ(ierr);
286   }
287 
288   /* setup coarse problem */
289   if (!def->WtAWinv) {
290     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */
291     if (!def->WtAW) {
292       ierr = PCGetOperators(pc,&Amat,NULL);CHKERRQ(ierr); /* using Amat! */
293       /* TODO add implicit product version ? */
294       if (!def->AW) {
295         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
296       } else {
297         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
298       }
299       /* TODO create MatInheritOption(Mat,MatOption) */
300       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
301       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
302 #if defined(PETSC_USE_DEBUG)
303       /* Check WtAW is not sigular */
304       PetscReal *norms;
305       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
306       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
307       for (i=0; i<m; i++) {
308         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
309           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
310         }
311       }
312       ierr = PetscFree(norms);CHKERRQ(ierr);
313 #endif
314     } else {
315       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
316     }
317     /* TODO use MATINV */
318     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
319     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
320     ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
321     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
322     ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
323     /* Reduction factor choice */
324     red = def->reductionfact;
325     if (red < 0) {
326       ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
327       red  = ceil((float)commsize/ceil((float)m/commsize));
328       ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
329       if (match) red = commsize;
330       ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
331     }
332     ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
333     ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
334     ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
335     /* Setup KSP and PC */
336     if (nextDef) { /* next level for multilevel deflation */
337       /* set default KSPtype */
338       if (!def->ksptype) {
339         def->ksptype = KSPFGMRES;
340         if (flgspd) { /* SPD system */
341           def->ksptype = KSPFCG;
342         }
343       }
344       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
345       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
346       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
347       ierr = PCDeflationSetNestLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
348       /* inherit options TODO if not set */
349       ((PC_Deflation*)(pcinner))->ksptype = def->ksptype;
350       ((PC_Deflation*)(pcinner))->correct = def->correct;
351       ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv;
352       ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst;
353       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
354     } else { /* the last level */
355       ierr = KSPSetType(innerksp,KSPPREONLY);CHKERRQ(ierr);
356       /* TODO Cholesky if flgspd? */
357       ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
358       //TODO remove explicit matSolverPackage
359       if (commsize == red) {
360         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);CHKERRQ(ierr);
361       } else {
362         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
363       }
364     }
365     /* TODO use def_[lvl]_ if lvl > 0? */
366     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
367     if (prefix) {
368       ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
369       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
370     } else {
371       ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
372     }
373     /* TODO: check if WtAWinv is KSP and move following from this if */
374     ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
375     //if (def->adaptiveconv) {
376     //  PetscReal *rnorm;
377     //  PetscNew(&rnorm);
378     //  ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr);
379     //}
380     ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
381   }
382 
383   ierr = KSPCreateVecs(def->WtAWinv,2,&def->work,0,NULL);CHKERRQ(ierr);
384   PetscFunctionReturn(0);
385 }
386 /* -------------------------------------------------------------------------- */
387 /*
388    PCApply_Deflation - Applies the Deflation preconditioner to a vector.
389 
390    Input Parameters:
391 .  pc - the preconditioner context
392 .  x - input vector
393 
394    Output Parameter:
395 .  y - output vector
396 
397    Application Interface Routine: PCApply()
398  */
399 static PetscErrorCode PCApply_Deflation(PC pc,Vec x,Vec y)
400 {
401   PC_Deflation      *jac = (PC_Deflation*)pc->data;
402   PetscErrorCode ierr;
403 
404   PetscFunctionBegin;
405   PetscFunctionReturn(0);
406 }
407 /* -------------------------------------------------------------------------- */
408 static PetscErrorCode PCReset_Deflation(PC pc)
409 {
410   PC_Deflation      *jac = (PC_Deflation*)pc->data;
411   PetscErrorCode ierr;
412 
413   PetscFunctionBegin;
414   PetscFunctionReturn(0);
415 }
416 
417 /*
418    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
419    that was created with PCCreate_Deflation().
420 
421    Input Parameter:
422 .  pc - the preconditioner context
423 
424    Application Interface Routine: PCDestroy()
425 */
426 static PetscErrorCode PCDestroy_Deflation(PC pc)
427 {
428   PetscErrorCode ierr;
429 
430   PetscFunctionBegin;
431   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
432 
433   /*
434       Free the private data structure that was hanging off the PC
435   */
436   ierr = PetscFree(pc->data);CHKERRQ(ierr);
437   PetscFunctionReturn(0);
438 }
439 
440 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
441 {
442   PC_Deflation      *jac = (PC_Deflation*)pc->data;
443   PetscErrorCode ierr;
444   PetscBool      flg;
445   PCDeflationType   deflt,type;
446 
447   PetscFunctionBegin;
448   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
449   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
450   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
451   if (flg) {
452     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
453   }
454   ierr = PetscOptionsTail();CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 /*MC
459      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
460      or to a predefined value
461 
462    Options Database Key:
463 +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
464 -    -pc_jacobi_abs - use the absolute value of the diagonal entry
465 
466    Level: beginner
467 
468   Notes:
469     todo
470 
471 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
472            PCDeflationSetType(), PCDeflationSetSpace()
473 M*/
474 
475 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
476 {
477   PC_Deflation   *def;
478   PetscErrorCode ierr;
479 
480   PetscFunctionBegin;
481   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
482   pc->data = (void*)def;
483 
484   def->init          = PETSC_FALSE;
485   def->pre           = PETSC_TRUE;
486   def->correct       = PETSC_FALSE;
487   def->truenorm      = PETSC_TRUE;
488   def->reductionfact = -1;
489   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
490   def->spacesize     = 1;
491   def->extendsp      = PETSC_FALSE;
492   def->nestedlvl     = 0;
493   def->maxnestedlvl  = 0;
494   def->adaptiveconv  = PETSC_FALSE;
495   def->adaptiveconst = 1.0;
496 
497   /*
498       Set the pointers for the functions that are provided above.
499       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
500       are called, they will automatically call these functions.  Note we
501       choose not to provide a couple of these functions since they are
502       not needed.
503   */
504   pc->ops->apply               = PCApply_Deflation;
505   pc->ops->applytranspose      = PCApply_Deflation;
506   pc->ops->setup               = PCSetUp_Deflation;
507   pc->ops->reset               = PCReset_Deflation;
508   pc->ops->destroy             = PCDestroy_Deflation;
509   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
510   pc->ops->view                = 0;
511   pc->ops->applyrichardson     = 0;
512   pc->ops->applysymmetricleft  = 0;
513   pc->ops->applysymmetricright = 0;
514 
515   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
516   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
517   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
518   PetscFunctionReturn(0);
519 }
520 
521