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