xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 117ef88edefbfc12e7c19efe87a19a2e1b0acd4f)
1 #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscksp.h" I*/  /* includes for fortran wrappers */
2 
3 const char *const PCDeflationSpaceTypes[] = {
4   "haar",
5   "db2",
6   "db4",
7   "db8",
8   "db16",
9   "biorth22",
10   "meyer",
11   "aggregation",
12   "user",
13   "PCDeflationSpaceType",
14   "PC_DEFLATION_SPACE_",
15   NULL
16 };
17 
18 static PetscErrorCode PCDeflationSetInitOnly_Deflation(PC pc,PetscBool flg)
19 {
20   PC_Deflation   *def = (PC_Deflation*)pc->data;
21 
22   PetscFunctionBegin;
23   def->init = flg;
24   PetscFunctionReturn(0);
25 }
26 
27 /*@
28    PCDeflationSetInitOnly - Do only initialization step.
29     Sets initial guess to the solution on the deflation space but does not apply
30     the deflation preconditioner. The additional preconditioner is still applied.
31 
32    Logically Collective
33 
34    Input Parameters:
35 +  pc  - the preconditioner context
36 -  flg - default PETSC_FALSE
37 
38    Options Database Keys:
39 .    -pc_deflation_init_only <false> - if true computes only the special guess
40 
41    Level: intermediate
42 
43 .seealso: PCDEFLATION
44 @*/
45 PetscErrorCode PCDeflationSetInitOnly(PC pc,PetscBool flg)
46 {
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
51   PetscValidLogicalCollectiveBool(pc,flg,2);
52   ierr = PetscTryMethod(pc,"PCDeflationSetInitOnly_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
53   PetscFunctionReturn(0);
54 }
55 
56 
57 static PetscErrorCode PCDeflationSetLevels_Deflation(PC pc,PetscInt current,PetscInt max)
58 {
59   PC_Deflation   *def = (PC_Deflation*)pc->data;
60 
61   PetscFunctionBegin;
62   if (current) def->lvl = current;
63   def->maxlvl = max;
64   PetscFunctionReturn(0);
65 }
66 
67 /*@
68    PCDeflationSetLevels - Set the maximum level of deflation nesting.
69 
70    Logically Collective
71 
72    Input Parameters:
73 +  pc  - the preconditioner context
74 -  max - maximum deflation level
75 
76    Options Database Keys:
77 .    -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation
78 
79    Level: intermediate
80 
81 .seealso: PCDeflationSetSpaceToCompute(), PCDeflationSetSpace(), PCDEFLATION
82 @*/
83 PetscErrorCode PCDeflationSetLevels(PC pc,PetscInt max)
84 {
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
89   PetscValidLogicalCollectiveInt(pc,max,2);
90   ierr = PetscTryMethod(pc,"PCDeflationSetLevels_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
91   PetscFunctionReturn(0);
92 }
93 
94 static PetscErrorCode PCDeflationSetReductionFactor_Deflation(PC pc,PetscInt red)
95 {
96   PC_Deflation   *def = (PC_Deflation*)pc->data;
97 
98   PetscFunctionBegin;
99   def->reductionfact = red;
100   PetscFunctionReturn(0);
101 }
102 
103 /*@
104    PCDeflationSetReductionFactor - Set reduction factor for the bottom PCTELESCOPE coarse problem solver.
105 
106    Logically Collective
107 
108    Input Parameters:
109 +  pc  - the preconditioner context
110 -  red - reduction factor (or PETSC_DETERMINE)
111 
112    Options Database Keys:
113 .    -pc_deflation_reduction_factor <\-1> - reduction factor on bottom level coarse problem for PCTELESCOPE
114 
115    Notes:
116      Default is computed based on the size of the coarse problem.
117 
118    Level: intermediate
119 
120 .seealso: PCTELESCOPE, PCDEFLATION
121 @*/
122 PetscErrorCode PCDeflationSetReductionFactor(PC pc,PetscInt red)
123 {
124   PetscErrorCode ierr;
125 
126   PetscFunctionBegin;
127   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
128   PetscValidLogicalCollectiveInt(pc,red,2);
129   ierr = PetscTryMethod(pc,"PCDeflationSetReductionFactor_C",(PC,PetscInt),(pc,red));CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 static PetscErrorCode PCDeflationSetCorrectionFactor_Deflation(PC pc,PetscScalar fact)
134 {
135   PC_Deflation   *def = (PC_Deflation*)pc->data;
136 
137   PetscFunctionBegin;
138   /* TODO PETSC_DETERMINE -> compute max eigenvalue with power method */
139   def->correct = PETSC_TRUE;
140   def->correctfact = fact;
141   if (def->correct == 0.0) {
142     def->correct = PETSC_FALSE;
143   }
144   PetscFunctionReturn(0);
145 }
146 
147 /*@
148    PCDeflationSetCorrectionFactor - Set coarse problem correction factor.
149     The Preconditioner becomes P*M^{-1} + fact*Q.
150 
151    Logically Collective
152 
153    Input Parameters:
154 +  pc   - the preconditioner context
155 -  fact - correction factor
156 
157    Options Database Keys:
158 +    -pc_deflation_correction        <false> - if true apply coarse problem correction
159 -    -pc_deflation_correction_factor <1.0>   - sets coarse problem correction factor
160 
161    Notes:
162     Any non-zero fact enables the coarse problem correction.
163 
164    Level: intermediate
165 
166 .seealso: PCDEFLATION
167 @*/
168 PetscErrorCode PCDeflationSetCorrectionFactor(PC pc,PetscScalar fact)
169 {
170   PetscErrorCode ierr;
171 
172   PetscFunctionBegin;
173   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
174   PetscValidLogicalCollectiveScalar(pc,fact,2);
175   ierr = PetscTryMethod(pc,"PCDeflationSetCorrectionFactor_C",(PC,PetscScalar),(pc,fact));CHKERRQ(ierr);
176   PetscFunctionReturn(0);
177 }
178 
179 static PetscErrorCode PCDeflationSetSpaceToCompute_Deflation(PC pc,PCDeflationSpaceType type,PetscInt size)
180 {
181   PC_Deflation   *def = (PC_Deflation*)pc->data;
182 
183   PetscFunctionBegin;
184   if (type) def->spacetype = type;
185   if (size > 0) def->spacesize = size;
186   PetscFunctionReturn(0);
187 }
188 
189 /*@
190    PCDeflationSetSpaceToCompute - Set deflation space type and size to compute.
191 
192    Logically Collective
193 
194    Input Parameters:
195 +  pc   - the preconditioner context
196 .  type - deflation space type to compute (or PETSC_IGNORE)
197 -  size - size of the space to compute (or PETSC_DEFAULT)
198 
199    Options Database Keys:
200 +    -pc_deflation_compute_space      <haar> - compute PCDeflationSpaceType deflation space
201 -    -pc_deflation_compute_space_size <1>    - size of the deflation space
202 
203    Notes:
204     For wavelet-based deflation, size represents number of levels.
205 
206     The deflation space is computed in PCSetUp().
207 
208    Level: intermediate
209 
210 .seealso: PCDeflationSetLevels(), PCDEFLATION
211 @*/
212 PetscErrorCode PCDeflationSetSpaceToCompute(PC pc,PCDeflationSpaceType type,PetscInt size)
213 {
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
218   if (type) PetscValidLogicalCollectiveEnum(pc,type,2);
219   if (size > 0) PetscValidLogicalCollectiveInt(pc,size,3);
220   ierr = PetscTryMethod(pc,"PCDeflationSetSpaceToCompute_C",(PC,PCDeflationSpaceType,PetscInt),(pc,type,size));CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
225 {
226   PC_Deflation   *def = (PC_Deflation*)pc->data;
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   /* possibly allows W' = Wt (which is valid but not tested) */
231   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
232   if (transpose) {
233     ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
234     def->Wt = W;
235   } else {
236     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
237     def->W = W;
238   }
239   PetscFunctionReturn(0);
240 }
241 
242 /*@
243    PCDeflationSetSpace - Set the deflation space matrix (or its (Hermitian) transpose).
244 
245    Logically Collective
246 
247    Input Parameters:
248 +  pc        - the preconditioner context
249 .  W         - deflation matrix
250 -  transpose - indicates that W is an explicit transpose of the deflation matrix
251 
252    Notes:
253     Setting W as a multipliplicative MATCOMPOSITE enables use of the multilevel
254     deflation. If W = W0*W1*W2*...*Wn, W0 is taken as the first deflation space and
255     the coarse problem (W0'*A*W0)^{-1} is again preconditioned by deflation with
256     W1 as the deflation matrix. This repeats until the maximum level set by
257     PCDeflationSetLevels() is reached or there are no more matrices available.
258     If there are matrices left after reaching the maximum level,
259     they are merged into a deflation matrix ...*W{n-1}*Wn.
260 
261    Level: intermediate
262 
263 .seealso: PCDeflationSetLevels(), PCDEFLATION
264 @*/
265 PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
266 {
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
271   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
272   PetscValidLogicalCollectiveBool(pc,transpose,3);
273   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
274   PetscFunctionReturn(0);
275 }
276 
277 static PetscErrorCode PCDeflationSetProjectionNullSpaceMat_Deflation(PC pc,Mat mat)
278 {
279   PC_Deflation     *def = (PC_Deflation*)pc->data;
280   PetscErrorCode   ierr;
281 
282   PetscFunctionBegin;
283   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
284   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
285   def->WtA = mat;
286   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtA);CHKERRQ(ierr);
287   PetscFunctionReturn(0);
288 }
289 
290 /*@
291    PCDeflationSetProjectionNullSpaceMat - Set the projection null space matrix (W'*A).
292 
293    Collective
294 
295    Input Parameters:
296 +  pc  - preconditioner context
297 -  mat - projection null space matrix
298 
299    Level: developer
300 
301 .seealso: PCDEFLATION
302 @*/
303 PetscErrorCode  PCDeflationSetProjectionNullSpaceMat(PC pc,Mat mat)
304 {
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
309   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
310   ierr = PetscTryMethod(pc,"PCDeflationSetProjectionNullSpaceMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
311   PetscFunctionReturn(0);
312 }
313 
314 static PetscErrorCode PCDeflationSetCoarseMat_Deflation(PC pc,Mat mat)
315 {
316   PC_Deflation     *def = (PC_Deflation*)pc->data;
317   PetscErrorCode   ierr;
318 
319   PetscFunctionBegin;
320   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
321   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
322   def->WtAW = mat;
323   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->WtAW);CHKERRQ(ierr);
324   PetscFunctionReturn(0);
325 }
326 
327 /*@
328    PCDeflationSetCoarseMat - Set the coarse problem Mat.
329 
330    Collective
331 
332    Input Parameters:
333 +  pc  - preconditioner context
334 -  mat - coarse problem mat
335 
336    Level: developer
337 
338 .seealso: PCDEFLATION
339 @*/
340 PetscErrorCode  PCDeflationSetCoarseMat(PC pc,Mat mat)
341 {
342   PetscErrorCode ierr;
343 
344   PetscFunctionBegin;
345   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
346   PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
347   ierr = PetscTryMethod(pc,"PCDeflationSetCoarseMat_C",(PC,Mat),(pc,mat));CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
352 {
353   PC_Deflation     *def = (PC_Deflation*)pc->data;
354 
355   PetscFunctionBegin;
356   *ksp = def->WtAWinv;
357   PetscFunctionReturn(0);
358 }
359 
360 /*@
361    PCDeflationGetCoarseKSP - Returns the coarse problem KSP.
362 
363    Not Collective
364 
365    Input Parameters:
366 .  pc - preconditioner context
367 
368    Output Parameters:
369 .  ksp - coarse problem KSP context
370 
371    Level: advanced
372 
373 .seealso: PCDEFLATION
374 @*/
375 PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
376 {
377   PetscErrorCode ierr;
378 
379   PetscFunctionBegin;
380   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
381   PetscValidPointer(ksp,2);
382   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
383   PetscFunctionReturn(0);
384 }
385 
386 static PetscErrorCode PCDeflationGetPC_Deflation(PC pc,PC *apc)
387 {
388   PC_Deflation   *def = (PC_Deflation*)pc->data;
389 
390   PetscFunctionBegin;
391   *apc = def->pc;
392   PetscFunctionReturn(0);
393 }
394 
395 /*@
396    PCDeflationGetPC - Returns the additional preconditioner M^{-1}.
397 
398    Not Collective
399 
400    Input Parameters:
401 .  pc  - the preconditioner context
402 
403    Output Parameters:
404 .  apc - additional preconditioner
405 
406    Level: advanced
407 
408 .seealso: PCDEFLATION
409 @*/
410 PetscErrorCode PCDeflationGetPC(PC pc,PC *apc)
411 {
412   PetscErrorCode ierr;
413 
414   PetscFunctionBegin;
415   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
416   PetscValidPointer(pc,2);
417   ierr = PetscTryMethod(pc,"PCDeflationGetPC_C",(PC,PC*),(pc,apc));CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 
421 /*
422   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
423 */
424 static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
425 {
426   PC_Deflation     *def = (PC_Deflation*)pc->data;
427   Mat              A;
428   Vec              r,w1,w2;
429   PetscBool        nonzero;
430   PetscErrorCode   ierr;
431 
432   PetscFunctionBegin;
433   w1 = def->workcoarse[0];
434   w2 = def->workcoarse[1];
435   r  = def->work;
436   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
437 
438   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
439   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
440   if (nonzero) {
441     ierr = MatMult(A,x,r);CHKERRQ(ierr);                          /*    r  <- b - Ax              */
442     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
443   } else {
444     ierr = VecCopy(b,r);CHKERRQ(ierr);                            /*    r  <- b (x is 0)          */
445   }
446 
447   if (def->Wt) {
448     ierr = MatMult(def->Wt,r,w1);CHKERRQ(ierr);                   /*    w1 <- W'*r                */
449   } else {
450     ierr = MatMultHermitianTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
451   }
452   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);              /*    w2 <- (W'*A*W)^{-1}*w1    */
453   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);                      /*    r  <- W*w2                */
454   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 /*
459   if (def->correct) {
460     z <- M^{-1}r - W*(W'*A*W)^{-1}*(W'*A*M^{-1}r - l*W'*r) = (P*M^{-1} + l*Q)*r
461   } else {
462     z <- M^{-1}*r - W*(W'*A*W)^{-1}*W'*A*M{-1}*r = P*M^{-1}*r
463   }
464 */
465 static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
466 {
467   PC_Deflation     *def = (PC_Deflation*)pc->data;
468   Mat              A;
469   Vec              u,w1,w2;
470   PetscErrorCode   ierr;
471 
472   PetscFunctionBegin;
473   w1 = def->workcoarse[0];
474   w2 = def->workcoarse[1];
475   u  = def->work;
476   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
477 
478   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);                          /*    z <- M^{-1}*r             */
479   if (!def->init) {
480     ierr = MatMult(def->WtA,z,w1);CHKERRQ(ierr);                      /*    w1 <- W'*A*z              */
481     if (def->correct) {
482       if (def->Wt) {
483         ierr = MatMult(def->Wt,r,w2);CHKERRQ(ierr);                   /*    w2 <- W'*r                */
484       } else {
485         ierr = MatMultHermitianTranspose(def->W,r,w2);CHKERRQ(ierr);  /*    w2 <- W'*r                */
486       }
487       ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);      /*    w1 <- w1 - l*w2           */
488     }
489     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);                /*    w2 <- (W'*A*W)^{-1}*w1    */
490     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                        /*    u  <- W*w2                */
491     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                           /*    z  <- z - u               */
492   }
493   PetscFunctionReturn(0);
494 }
495 
496 static PetscErrorCode PCSetUp_Deflation(PC pc)
497 {
498   PC_Deflation     *def = (PC_Deflation*)pc->data;
499   KSP              innerksp;
500   PC               pcinner;
501   Mat              Amat,nextDef=NULL,*mats;
502   PetscInt         i,m,red,size;
503   PetscMPIInt      commsize;
504   PetscBool        match,flgspd,transp=PETSC_FALSE;
505   MatCompositeType ctype;
506   MPI_Comm         comm;
507   char             prefix[128]="";
508   PetscErrorCode   ierr;
509 
510   PetscFunctionBegin;
511   if (pc->setupcalled) PetscFunctionReturn(0);
512   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
513   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
514   if (!def->lvl && !def->prefix) {
515     ierr = PCGetOptionsPrefix(pc,&def->prefix);CHKERRQ(ierr);
516   }
517   if (def->lvl) {
518     ierr = PetscSNPrintf(prefix,sizeof(prefix),"%d_",(int)def->lvl);CHKERRQ(ierr);
519   }
520 
521   /* compute a deflation space */
522   if (def->W || def->Wt) {
523     def->spacetype = PC_DEFLATION_SPACE_USER;
524   } else {
525     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
526   }
527 
528   /* nested deflation */
529   if (def->W) {
530     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
531     if (match) {
532       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
533       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
534     }
535   } else {
536     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
537     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
538     if (match) {
539       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
540       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
541     }
542     transp = PETSC_TRUE;
543   }
544   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
545     if (!transp) {
546       if (def->lvl < def->maxlvl) {
547         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
548         for (i=0; i<size; i++) {
549           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
550         }
551         size -= 1;
552         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
553         def->W = mats[size];
554         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
555         if (size > 1) {
556           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
557           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
558         } else {
559           nextDef = mats[0];
560           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
561         }
562         ierr = PetscFree(mats);CHKERRQ(ierr);
563       } else {
564         /* TODO test merge side performance */
565         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
566         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
567       }
568     } else {
569       if (def->lvl < def->maxlvl) {
570         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
571         for (i=0; i<size; i++) {
572           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
573         }
574         size -= 1;
575         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
576         def->Wt = mats[0];
577         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
578         if (size > 1) {
579           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
580           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
581         } else {
582           nextDef = mats[1];
583           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
584         }
585         ierr = PetscFree(mats);CHKERRQ(ierr);
586       } else {
587         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
588         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
589       }
590     }
591   }
592 
593   if (transp) {
594     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
595     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
596   }
597 
598   /* assemble WtA */
599   if (!def->WtA) {
600     if (def->Wt) {
601       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
602     } else {
603 #if defined(PETSC_USE_COMPLEX)
604       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
605       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
606 #else
607       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
608 #endif
609     }
610   }
611   /* setup coarse problem */
612   if (!def->WtAWinv) {
613     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
614     if (!def->WtAW) {
615       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
616       /* TODO create MatInheritOption(Mat,MatOption) */
617       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
618       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
619       if (PetscDefined(USE_DEBUG)) {
620         /* Check columns of W are not in kernel of A */
621         PetscReal *norms;
622         ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
623         ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
624         for (i=0; i<m; i++) {
625           if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
626             SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
627           }
628         }
629         ierr = PetscFree(norms);CHKERRQ(ierr);
630       }
631     } else {
632       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
633     }
634     /* TODO use MATINV ? */
635     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
636     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
637     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
638     /* Setup KSP and PC */
639     if (nextDef) { /* next level for multilevel deflation */
640       innerksp = def->WtAWinv;
641       /* set default KSPtype */
642       if (!def->ksptype) {
643         def->ksptype = KSPFGMRES;
644         if (flgspd) { /* SPD system */
645           def->ksptype = KSPFCG;
646         }
647       }
648       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
649       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
650       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
651       ierr = PCDeflationSetLevels_Deflation(pcinner,def->lvl+1,def->maxlvl);CHKERRQ(ierr);
652       /* inherit options */
653       if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix;
654       ((PC_Deflation*)(pcinner->data))->init          = def->init;
655       ((PC_Deflation*)(pcinner->data))->ksptype       = def->ksptype;
656       ((PC_Deflation*)(pcinner->data))->correct       = def->correct;
657       ((PC_Deflation*)(pcinner->data))->correctfact   = def->correctfact;
658       ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact;
659       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
660     } else { /* the last level */
661       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
662       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
663       /* do not overwrite PCTELESCOPE */
664       if (def->prefix) {
665         ierr = KSPSetOptionsPrefix(def->WtAWinv,def->prefix);CHKERRQ(ierr);
666       }
667       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"deflation_tel_");CHKERRQ(ierr);
668       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
669       ierr = PetscObjectTypeCompare((PetscObject)pcinner,PCTELESCOPE,&match);CHKERRQ(ierr);
670       if (!match) SETERRQ(comm,PETSC_ERR_SUP,"User can not owerwrite PCTELESCOPE on bottom level, use reduction factor = 1 instead.");
671       /* Reduction factor choice */
672       red = def->reductionfact;
673       if (red < 0) {
674         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
675         red  = ceil((float)commsize/ceil((float)m/commsize));
676         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
677         if (match) red = commsize;
678         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr);
679       }
680       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
681       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
682       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
683       if (innerksp) {
684         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
685         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
686 #if defined(PETSC_HAVE_SUPERLU)
687         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
688         if (match) {
689           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
690         }
691 #endif
692 #if defined(PETSC_HAVE_SUPERLU_DIST)
693         ierr = MatGetFactorAvailable(def->WtAW,MATSOLVERSUPERLU_DIST,MAT_FACTOR_LU,&match);CHKERRQ(ierr);
694         if (match) {
695           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
696         }
697 #endif
698       }
699     }
700 
701     if (innerksp) {
702       if (def->prefix) {
703         ierr = KSPSetOptionsPrefix(innerksp,def->prefix);CHKERRQ(ierr);
704         ierr = KSPAppendOptionsPrefix(innerksp,"deflation_");CHKERRQ(ierr);
705       } else {
706         ierr = KSPSetOptionsPrefix(innerksp,"deflation_");CHKERRQ(ierr);
707       }
708       ierr = KSPAppendOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
709       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
710       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
711     }
712   }
713   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
714   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
715 
716   /* create preconditioner */
717   if (!def->pc) {
718     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
719     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
720     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
721     if (def->prefix) {
722       ierr = PCSetOptionsPrefix(def->pc,def->prefix);CHKERRQ(ierr);
723     }
724     ierr = PCAppendOptionsPrefix(def->pc,"deflation_");CHKERRQ(ierr);
725     ierr = PCAppendOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
726     ierr = PCAppendOptionsPrefix(def->pc,"pc_");CHKERRQ(ierr);
727     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
728     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
729   }
730 
731   /* create work vecs */
732   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
733   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
734   PetscFunctionReturn(0);
735 }
736 
737 static PetscErrorCode PCReset_Deflation(PC pc)
738 {
739   PC_Deflation      *def = (PC_Deflation*)pc->data;
740   PetscErrorCode    ierr;
741 
742   PetscFunctionBegin;
743   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
744   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
745   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
746   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
747   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
748   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
749   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
750   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
751   PetscFunctionReturn(0);
752 }
753 
754 static PetscErrorCode PCDestroy_Deflation(PC pc)
755 {
756   PetscErrorCode ierr;
757 
758   PetscFunctionBegin;
759   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
760   ierr = PetscFree(pc->data);CHKERRQ(ierr);
761   PetscFunctionReturn(0);
762 }
763 
764 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
765 {
766   PC_Deflation      *def = (PC_Deflation*)pc->data;
767   PetscInt          its;
768   PetscBool         iascii;
769   PetscErrorCode    ierr;
770 
771   PetscFunctionBegin;
772   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
773   if (iascii) {
774     if (def->correct) {
775       ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n",
776                                     (double)PetscRealPart(def->correctfact),
777                                     (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr);
778     }
779     if (!def->lvl) {
780       ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
781     }
782 
783     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr);
784     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
785     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
786     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
787 
788     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr);
789     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
790     ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr);
791     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr);
792     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
793     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
794   }
795   PetscFunctionReturn(0);
796 }
797 
798 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
799 {
800   PC_Deflation      *def = (PC_Deflation*)pc->data;
801   PetscErrorCode    ierr;
802 
803   PetscFunctionBegin;
804   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
805   ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr);
806   ierr = PetscOptionsInt("-pc_deflation_levels","Maximum of deflation levels","PCDeflationSetLevels",def->maxlvl,&def->maxlvl,NULL);CHKERRQ(ierr);
807   ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
808   ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr);
809   ierr = PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
810   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
811   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
812   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
813   ierr = PetscOptionsTail();CHKERRQ(ierr);
814   PetscFunctionReturn(0);
815 }
816 
817 /*MC
818    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
819 
820    Options Database Keys:
821 +    -pc_deflation_init_only          <false> - if true computes only the special guess
822 .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
823 .    -pc_deflation_reduction_factor <\-1>     - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
824 .    -pc_deflation_correction         <false> - if true apply coarse problem correction
825 .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
826 .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
827 -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
828 
829    Notes:
830     Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
831     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
832 
833     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
834     If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the
835     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
836     application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
837     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
838     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
839     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
840 
841     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
842     be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
843     User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix
844     is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the
845     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned
846     by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum
847     level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged
848     (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by
849     PCDeflationSetLevels() or by -pc_deflation_levels.
850 
851     The coarse problem KSP can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1]
852     from the second level onward. You can also use
853     PCDeflationGetCoarseKSP() to control it from code. The bottom level KSP defaults to
854     KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE.
855     For convenience, the reduction factor can be set by PCDeflationSetReductionFactor()
856     or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
857 
858     The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for
859     coarse problem KSP apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use
860     PCDeflationGetPC() to control the additional preconditioner from code. It defaults to PCNONE.
861 
862     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
863     be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can
864     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
865     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
866     an isolated eigenvalue.
867 
868     The options are automatically inherited from the previous deflation level.
869 
870     The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also
871     recommend limiting the number of iterations for the coarse problems.
872 
873     See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
874     Section 4 describes some possible choices for the deflation space.
875 
876    Developer Notes:
877      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
878      Academy of Sciences and VSB - TU Ostrava.
879 
880      Developed from PERMON code used in [4] while on a research stay with
881      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
882 
883    References:
884 +    [1] - A. Nicolaides. “Deflation of conjugate gradients with applications to boundary valueproblems”, SIAM J. Numer. Anal. 24.2, 1987.
885 .    [2] - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
886 .    [3] - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
887 -    [4] - J. Kruzik "Implementation of the Deflated Variants of the Conjugate Gradient Method", Master's thesis, VSB-TUO, 2018 - http://dspace5.vsb.cz/bitstream/handle/10084/130303/KRU0097_USP_N2658_2612T078_2018.pdf
888 
889    Level: intermediate
890 
891 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
892            PCDeflationSetInitOnly(), PCDeflationSetLevels(), PCDeflationSetReductionFactor(),
893            PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompute(),
894            PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(),
895            PCDeflationSetCoarseMat(), PCDeflationGetCoarseKSP(), PCDeflationGetPC()
896 M*/
897 
898 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
899 {
900   PC_Deflation   *def;
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
905   pc->data = (void*)def;
906 
907   def->init          = PETSC_FALSE;
908   def->correct       = PETSC_FALSE;
909   def->correctfact   = 1.0;
910   def->reductionfact = -1;
911   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
912   def->spacesize     = 1;
913   def->extendsp      = PETSC_FALSE;
914   def->lvl           = 0;
915   def->maxlvl        = 0;
916   def->W             = NULL;
917   def->Wt            = NULL;
918 
919   pc->ops->apply          = PCApply_Deflation;
920   pc->ops->presolve       = PCPreSolve_Deflation;
921   pc->ops->setup          = PCSetUp_Deflation;
922   pc->ops->reset          = PCReset_Deflation;
923   pc->ops->destroy        = PCDestroy_Deflation;
924   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
925   pc->ops->view           = PCView_Deflation;
926 
927   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr);
928   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLevels_C",PCDeflationSetLevels_Deflation);CHKERRQ(ierr);
929   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr);
930   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr);
931   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr);
932   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
933   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr);
934   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
935   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
936   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
937   PetscFunctionReturn(0);
938 }
939 
940