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