xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 5c0c31c5171c0ed0d0162135c9816bb625821204)
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 static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
209 {
210   PC_Deflation   *def = (PC_Deflation*)pc->data;
211 
212   PetscFunctionBegin;
213   def->nestedlvl = current;
214   def->maxnestedlvl = max;
215   PetscFunctionReturn(0);
216 }
217 
218 /*@
219    PCDeflationSetMaxLvl - Set maximum level of deflation.
220 
221    Logically Collective on PC
222 
223    Input Parameters:
224 +  pc  - the preconditioner context
225 .  max - maximum deflation level
226 
227    Level: intermediate
228 
229 .seealso: PCDEFLATION
230 @*/
231 PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
232 {
233   PetscErrorCode ierr;
234 
235   PetscFunctionBegin;
236   PetscValidLogicalCollectiveInt(pc,max,2);
237   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
238   PetscFunctionReturn(0);
239 }
240 
241 /* -------------------------------------------------------------------------- */
242 /*
243    PCSetUp_Deflation - Prepares for the use of the Deflation preconditioner
244                     by setting data structures and options.
245 
246    Input Parameter:
247 .  pc - the preconditioner context
248 
249    Application Interface Routine: PCSetUp()
250 
251    Notes:
252    The interface routine PCSetUp() is not usually called directly by
253    the user, but instead is called by PCApply() if necessary.
254 */
255 static PetscErrorCode PCSetUp_Deflation(PC pc)
256 {
257   PC_Deflation     *def = (PC_Deflation*)pc->data;
258   KSP              innerksp;
259   PC               pcinner;
260   Mat              Amat,nextDef=NULL,*mats;
261   PetscInt         i,m,red,size,commsize;
262   PetscBool        match,flgspd,transp=PETSC_FALSE;
263   MatCompositeType ctype;
264   MPI_Comm         comm;
265   const char       *prefix;
266   PetscErrorCode   ierr;
267 
268   PetscFunctionBegin;
269   if (pc->setupcalled) PetscFunctionReturn(0);
270   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
271   if (def->W || def->Wt) {
272     def->spacetype = PC_DEFLATION_SPACE_USER;
273   } else {
274     //ierr = KSPDCGComputeDeflationSpace(ksp);CHKERRQ(ierr);
275   }
276 
277   /* nested deflation */
278   if (def->W) {
279     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
280     if (match) {
281       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
282       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
283     }
284   } else {
285     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
286     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
287     if (match) {
288       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
289       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
290     }
291     transp = PETSC_TRUE;
292   }
293   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
294     ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
295     if (!transp) {
296       for (i=0; i<size; i++) {
297         ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
298         //ierr = PetscObjectReference((PetscObject)mats[i]);CHKERRQ(ierr);
299       }
300       if (def->nestedlvl < def->maxnestedlvl) {
301         size -= 1;
302         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
303         def->W = mats[size];
304         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
305         if (size > 1) {
306           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
307           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
308         } else {
309           nextDef = mats[0];
310           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
311         }
312       } else {
313         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
314         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
315       }
316     }
317     ierr = PetscFree(mats);CHKERRQ(ierr);
318   }
319 
320   /* setup coarse problem */
321   if (!def->WtAWinv) {
322     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr); /* TODO works for W MatTranspose? */
323     if (!def->WtAW) {
324       ierr = PCGetOperators(pc,&Amat,NULL);CHKERRQ(ierr); /* using Amat! */
325       /* TODO add implicit product version ? */
326       if (!def->AW) {
327         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
328       } else {
329         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
330       }
331       /* TODO create MatInheritOption(Mat,MatOption) */
332       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
333       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
334 #if defined(PETSC_USE_DEBUG)
335       /* Check WtAW is not sigular */
336       PetscReal *norms;
337       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
338       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
339       for (i=0; i<m; i++) {
340         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
341           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
342         }
343       }
344       ierr = PetscFree(norms);CHKERRQ(ierr);
345 #endif
346     } else {
347       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
348     }
349     /* TODO use MATINV */
350     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
351     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
352     ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
353     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
354     ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
355     /* Reduction factor choice */
356     red = def->reductionfact;
357     if (red < 0) {
358       ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
359       red  = ceil((float)commsize/ceil((float)m/commsize));
360       ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
361       if (match) red = commsize;
362       ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
363     }
364     ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
365     ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
366     ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
367     /* Setup KSP and PC */
368     if (nextDef) { /* next level for multilevel deflation */
369       /* set default KSPtype */
370       if (!def->ksptype) {
371         def->ksptype = KSPFGMRES;
372         if (flgspd) { /* SPD system */
373           def->ksptype = KSPFCG;
374         }
375       }
376       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
377       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
378       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
379       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
380       /* inherit options TODO if not set */
381       ((PC_Deflation*)(pcinner))->ksptype = def->ksptype;
382       ((PC_Deflation*)(pcinner))->correct = def->correct;
383       ((PC_Deflation*)(pcinner))->adaptiveconv = def->adaptiveconv;
384       ((PC_Deflation*)(pcinner))->adaptiveconst = def->adaptiveconst;
385       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
386     } else { /* the last level */
387       ierr = KSPSetType(innerksp,KSPPREONLY);CHKERRQ(ierr);
388       /* TODO Cholesky if flgspd? */
389       ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
390       //TODO remove explicit matSolverPackage
391       if (commsize == red) {
392         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);CHKERRQ(ierr);
393       } else {
394         ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
395       }
396     }
397     /* TODO use def_[lvl]_ if lvl > 0? */
398     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
399     if (prefix) {
400       ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
401       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
402     } else {
403       ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
404     }
405     /* TODO: check if WtAWinv is KSP and move following from this if */
406     ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
407     //if (def->adaptiveconv) {
408     //  PetscReal *rnorm;
409     //  PetscNew(&rnorm);
410     //  ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr);
411     //}
412     ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
413   }
414 
415   ierr = KSPCreateVecs(def->WtAWinv,2,&def->work,0,NULL);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 /* -------------------------------------------------------------------------- */
419 /*
420    PCApply_Deflation - Applies the Deflation preconditioner to a vector.
421 
422    Input Parameters:
423 .  pc - the preconditioner context
424 .  x - input vector
425 
426    Output Parameter:
427 .  y - output vector
428 
429    Application Interface Routine: PCApply()
430  */
431 static PetscErrorCode PCApply_Deflation(PC pc,Vec x,Vec y)
432 {
433   PC_Deflation      *jac = (PC_Deflation*)pc->data;
434   PetscErrorCode ierr;
435 
436   PetscFunctionBegin;
437   PetscFunctionReturn(0);
438 }
439 /* -------------------------------------------------------------------------- */
440 static PetscErrorCode PCReset_Deflation(PC pc)
441 {
442   PC_Deflation      *jac = (PC_Deflation*)pc->data;
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   PetscFunctionReturn(0);
447 }
448 
449 /*
450    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
451    that was created with PCCreate_Deflation().
452 
453    Input Parameter:
454 .  pc - the preconditioner context
455 
456    Application Interface Routine: PCDestroy()
457 */
458 static PetscErrorCode PCDestroy_Deflation(PC pc)
459 {
460   PetscErrorCode ierr;
461 
462   PetscFunctionBegin;
463   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
464 
465   /*
466       Free the private data structure that was hanging off the PC
467   */
468   ierr = PetscFree(pc->data);CHKERRQ(ierr);
469   PetscFunctionReturn(0);
470 }
471 
472 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
473 {
474   PC_Deflation      *jac = (PC_Deflation*)pc->data;
475   PetscErrorCode ierr;
476   PetscBool      flg;
477   PCDeflationType   deflt,type;
478 
479   PetscFunctionBegin;
480   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
481   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
482   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
483   if (flg) {
484     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
485   }
486   ierr = PetscOptionsTail();CHKERRQ(ierr);
487   PetscFunctionReturn(0);
488 }
489 
490 /*MC
491      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
492      or to a predefined value
493 
494    Options Database Key:
495 +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
496 -    -pc_jacobi_abs - use the absolute value of the diagonal entry
497 
498    Level: beginner
499 
500   Notes:
501     todo
502 
503 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
504            PCDeflationSetType(), PCDeflationSetSpace()
505 M*/
506 
507 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
508 {
509   PC_Deflation   *def;
510   PetscErrorCode ierr;
511 
512   PetscFunctionBegin;
513   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
514   pc->data = (void*)def;
515 
516   def->init          = PETSC_FALSE;
517   def->pre           = PETSC_TRUE;
518   def->correct       = PETSC_FALSE;
519   def->truenorm      = PETSC_TRUE;
520   def->reductionfact = -1;
521   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
522   def->spacesize     = 1;
523   def->extendsp      = PETSC_FALSE;
524   def->nestedlvl     = 0;
525   def->maxnestedlvl  = 0;
526   def->adaptiveconv  = PETSC_FALSE;
527   def->adaptiveconst = 1.0;
528 
529   /*
530       Set the pointers for the functions that are provided above.
531       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
532       are called, they will automatically call these functions.  Note we
533       choose not to provide a couple of these functions since they are
534       not needed.
535   */
536   pc->ops->apply               = PCApply_Deflation;
537   pc->ops->applytranspose      = PCApply_Deflation;
538   pc->ops->setup               = PCSetUp_Deflation;
539   pc->ops->reset               = PCReset_Deflation;
540   pc->ops->destroy             = PCDestroy_Deflation;
541   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
542   pc->ops->view                = 0;
543   pc->ops->applyrichardson     = 0;
544   pc->ops->applysymmetricleft  = 0;
545   pc->ops->applysymmetricright = 0;
546 
547   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
548   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
549   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
550   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
551   PetscFunctionReturn(0);
552 }
553 
554