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