xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 9f604af8b6d29596cb0cec5adfb84567a79999b9)
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->lvl = current;
66   def->maxlvl = 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   char             prefix[128]="";
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   if (!def->lvl && !def->prefix) {
591     ierr = PCGetOptionsPrefix(pc,&def->prefix);CHKERRQ(ierr);
592   }
593   if (def->lvl) {
594     sprintf(prefix,"%d_",(int)def->lvl);
595   }
596 
597   /* compute a deflation space */
598   if (def->W || def->Wt) {
599     def->spacetype = PC_DEFLATION_SPACE_USER;
600   } else {
601     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
602   }
603 
604   /* nested deflation */
605   if (def->W) {
606     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
607     if (match) {
608       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
609       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
610     }
611   } else {
612     ierr = MatCreateHermitianTranspose(def->Wt,&def->W);CHKERRQ(ierr);
613     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
614     if (match) {
615       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
616       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
617     }
618     transp = PETSC_TRUE;
619   }
620   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
621     if (!transp) {
622       if (def->lvl < def->maxlvl) {
623         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
624         for (i=0; i<size; i++) {
625           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
626         }
627         size -= 1;
628         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
629         def->W = mats[size];
630         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
631         if (size > 1) {
632           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
633           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
634         } else {
635           nextDef = mats[0];
636           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
637         }
638         ierr = PetscFree(mats);CHKERRQ(ierr);
639       } else {
640         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
641         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
642       }
643     } else {
644       if (def->lvl < def->maxlvl) {
645         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
646         for (i=0; i<size; i++) {
647           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
648         }
649         size -= 1;
650         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
651         def->Wt = mats[0];
652         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
653         if (size > 1) {
654           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
655           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
656         } else {
657           nextDef = mats[1];
658           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
659         }
660         ierr = PetscFree(mats);CHKERRQ(ierr);
661       } else {
662         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
663         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
664       }
665     }
666   }
667 
668   if (transp) {
669     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
670     ierr = MatHermitianTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
671   }
672 
673   /* assemble WtA */
674   if (!def->WtA) {
675     if (def->Wt) {
676       ierr = MatMatMult(def->Wt,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
677     } else {
678 #if defined(PETSC_USE_COMPLEX)
679       ierr = MatHermitianTranspose(def->W,MAT_INITIAL_MATRIX,&def->Wt);CHKERRQ(ierr);
680       ierr = MatMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
681 #else
682       ierr = MatTransposeMatMult(def->W,Amat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtA);CHKERRQ(ierr);
683 #endif
684     }
685   }
686   /* setup coarse problem */
687   if (!def->WtAWinv) {
688     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
689     if (!def->WtAW) {
690       ierr = MatMatMult(def->WtA,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
691       /* TODO create MatInheritOption(Mat,MatOption) */
692       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
693       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
694 #if defined(PETSC_USE_DEBUG)
695       /* Check columns of W are not in kernel of A */
696       PetscReal *norms;
697       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
698       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
699       for (i=0; i<m; i++) {
700         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
701           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
702         }
703       }
704       ierr = PetscFree(norms);CHKERRQ(ierr);
705 #endif
706     } else {
707       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
708     }
709     /* TODO use MATINV */
710     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
711     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
712     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
713     /* Setup KSP and PC */
714     if (nextDef) { /* next level for multilevel deflation */
715       innerksp = def->WtAWinv;
716       /* set default KSPtype */
717       if (!def->ksptype) {
718         def->ksptype = KSPFGMRES;
719         if (flgspd) { /* SPD system */
720           def->ksptype = KSPFCG;
721         }
722       }
723       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP + tolerances */
724       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
725       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
726       ierr = PCDeflationSetLvl_Deflation(pcinner,def->lvl+1,def->maxlvl);CHKERRQ(ierr);
727       /* inherit options */
728       if (def->prefix) ((PC_Deflation*)(pcinner->data))->prefix = def->prefix;
729       ((PC_Deflation*)(pcinner->data))->init          = def->init;
730       ((PC_Deflation*)(pcinner->data))->ksptype       = def->ksptype;
731       ((PC_Deflation*)(pcinner->data))->correct       = def->correct;
732       ((PC_Deflation*)(pcinner->data))->correctfact   = def->correctfact;
733       ((PC_Deflation*)(pcinner->data))->reductionfact = def->reductionfact;
734       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
735     } else { /* the last level */
736       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
737       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
738       /* do not overwrite PCTELESCOPE */
739       if (def->prefix) {
740         ierr = KSPSetOptionsPrefix(def->WtAWinv,def->prefix);CHKERRQ(ierr);
741       }
742       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"def_tel_");CHKERRQ(ierr);
743       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
744       /* Reduction factor choice */
745       red = def->reductionfact;
746       if (red < 0) {
747         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
748         red  = ceil((float)commsize/ceil((float)m/commsize));
749         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
750         if (match) red = commsize;
751         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr);
752       }
753       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
754       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
755       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
756       if (innerksp) {
757         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
758         /* TODO Cholesky if flgspd? */
759         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
760         //TODO remove explicit matSolverPackage
761         if (commsize == red) {
762           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
763         } else {
764           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
765         }
766       }
767     }
768 
769     if (innerksp) {
770       if (def->prefix) {
771         ierr = KSPSetOptionsPrefix(innerksp,def->prefix);CHKERRQ(ierr);
772         ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
773       } else {
774         ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
775       }
776       ierr = KSPAppendOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
777       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
778       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
779     }
780   }
781   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
782   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
783 
784   if (!def->pc) {
785     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
786     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
787     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
788     if (def->prefix) {
789       ierr = PCSetOptionsPrefix(def->pc,def->prefix);CHKERRQ(ierr);
790     }
791     ierr = PCAppendOptionsPrefix(def->pc,"def_");CHKERRQ(ierr);
792     ierr = PCAppendOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
793     ierr = PCAppendOptionsPrefix(def->pc,"pc_");CHKERRQ(ierr);
794     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
795     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
796   }
797 
798   /* create work vecs */
799   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
800   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
801   PetscFunctionReturn(0);
802 }
803 
804 static PetscErrorCode PCReset_Deflation(PC pc)
805 {
806   PC_Deflation      *def = (PC_Deflation*)pc->data;
807   PetscErrorCode    ierr;
808 
809   PetscFunctionBegin;
810   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
811   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
812   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
813   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
814   ierr = MatDestroy(&def->WtA);CHKERRQ(ierr);
815   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
816   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
817   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
818   PetscFunctionReturn(0);
819 }
820 
821 static PetscErrorCode PCDestroy_Deflation(PC pc)
822 {
823   PetscErrorCode ierr;
824 
825   PetscFunctionBegin;
826   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
827   ierr = PetscFree(pc->data);CHKERRQ(ierr);
828   PetscFunctionReturn(0);
829 }
830 
831 static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
832 {
833   PC_Deflation      *def = (PC_Deflation*)pc->data;
834   PetscInt          its;
835   PetscBool         iascii;
836   PetscErrorCode    ierr;
837 
838   PetscFunctionBegin;
839   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
840   if (iascii) {
841     if (def->correct) {
842       ierr = PetscViewerASCIIPrintf(viewer,"using CP correction, factor = %g+%gi\n",
843                                     (double)PetscRealPart(def->correctfact),
844                                     (double)PetscImaginaryPart(def->correctfact));CHKERRQ(ierr);
845     }
846     if (!def->lvl) {
847       ierr = PetscViewerASCIIPrintf(viewer,"deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
848     }
849 
850     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC:\n");CHKERRQ(ierr);
851     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
852     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
853     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
854 
855     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse problem solver:\n");CHKERRQ(ierr);
856     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
857     ierr = KSPGetTotalIterations(def->WtAWinv,&its);CHKERRQ(ierr);
858     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations: %D\n",its);CHKERRQ(ierr);
859     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
860     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
861   }
862   PetscFunctionReturn(0);
863 }
864 
865 static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
866 {
867   PC_Deflation      *def = (PC_Deflation*)pc->data;
868   PetscErrorCode    ierr;
869 
870   PetscFunctionBegin;
871   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
872   ierr = PetscOptionsBool("-pc_deflation_init_only","Use only initialization step - Initdef","PCDeflationSetInitOnly",def->init,&def->init,NULL);CHKERRQ(ierr);
873   ierr = PetscOptionsInt("-pc_deflation_max_lvl","Maximum of deflation levels","PCDeflationSetMaxLvl",def->maxlvl,&def->maxlvl,NULL);CHKERRQ(ierr);
874   ierr = PetscOptionsInt("-pc_deflation_reduction_factor","Reduction factor for coarse problem solution using PCTELESCOPE","PCDeflationSetReductionFactor",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
875   ierr = PetscOptionsBool("-pc_deflation_correction","Add coarse problem correction Q to P","PCDeflationSetCorrectionFactor",def->correct,&def->correct,NULL);CHKERRQ(ierr);
876   ierr = PetscOptionsScalar("-pc_deflation_correction_factor","Set multiple of Q to use as coarse problem correction","PCDeflationSetCorrectionFactor",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
877   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
878   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
879   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
880 //TODO add set function and fix manpages
881   ierr = PetscOptionsTail();CHKERRQ(ierr);
882   PetscFunctionReturn(0);
883 }
884 
885 /*MC
886    PCDEFLATION - Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
887 
888    Options Database Keys:
889 +    -pc_deflation_init_only          <false> - if true computes only the special guess
890 .    -pc_deflation_max_lvl            <0>     - maximum number of levels for multilevel deflation
891 .    -pc_deflation_reduction_factor   <-1>    - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)
892 .    -pc_deflation_correction         <false> - if true apply coarse problem correction
893 .    -pc_deflation_correction_factor  <1.0>   - sets coarse problem correction factor
894 .    -pc_deflation_compute_space      <haar>  - compute PCDeflationSpaceType deflation space
895 -    -pc_deflation_compute_space_size <1>     - size of the deflation space (corresponds to number of levels for wavelet-based deflation)
896 
897    Notes:
898     Given a (complex - transpose is always hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
899     preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
900 
901     The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
902     If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the
903     preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
904     application consists of P*M{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
905     to zero while the addition of the coarse problem correction (factor*Q) makes the preconditioner to shift some
906     eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates
907     of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.
908 
909     The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
910     be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
911     User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix
912     is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the
913     deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by FCG (if A is MAT_SPD) or FGMRES preconditioned
914     by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum
915     level is reached or there are no more matrices. If the maximal level is reached, the remaining matrices are merged
916     (multiplied) to create the last deflation matrix. The maximal level defaults to 0 and can be set by
917     PCDeflationSetMaxLvl() or by -pc_deflation_max_lvl.
918 
919     The coarse problem KSP can be controlled from the command line with prefix -def_ for the first level and -def_[lvl-1]
920     from the second level onward. You can also use
921     PCDeflationGetCoarseKSP() or PCDeflationSetCoarseKSP() to control it from code. The bottom level KSP defaults to
922     KSPREONLY with PCLU direct solver wrapped into PCTELESCOPE. For convenience, the reduction factor can be set by
923     PCDeflationSetReductionFactor() or -pc_deflation_recduction_factor. The default is chosen heuristically based on the
924     coarse problem size.
925 
926     The additional preconditioner can be controlled from command line with prefix -def_[lvl]_pc (same rules used for
927     coarse problem KSP apply for [lvl]_ part of prefix), e.g., -def_1_pc_pc_type bjacobi. You can also use
928     PCDeflationGetPC() or PCDeflationSetPC() to control the additional preconditioner from code. It defaults to PCNONE.
929 
930     The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can
931     be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can
932     significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We
933     recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create
934     an isolated eigenvalue.
935 
936     The options are automatically inherited from previous deflation level.
937 
938     The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also
939     recommend limiting the number of iterations for the coarse problem.
940 
941     See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
942     Section 4 describes some possible choices for the deflation space.
943 
944    Developer Notes:
945      Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
946      Academy of Sciences and VSB - TU Ostrava.
947 
948      Developed from PERMON code used in [4] while on a research stay with
949      Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
950 
951    References:
952 +    [1] - A. Nicolaides. “Deflation of conjugate gradients with applications to boundary valueproblems”, SIAM J. Numer. Anal. 24.2, 1987.
953 .    [2] - Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
954 .    [3] - Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
955 .    [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
956 
957    Level: intermediate
958 
959 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
960            PCDeflationSetInitOnly(), PCDeflationSetMaxLvl(), PCDeflationSetReductionFactor(),
961            PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompoute(),
962            PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(),
963            PCDeflationSetCoarseMat(), PCDeflationSetCoarseKSP(), PCDeflationGetCoarseKSP(),
964            PCDeflationGetPC(), PCDeflationSetPC()
965 M*/
966 
967 PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
968 {
969   PC_Deflation   *def;
970   PetscErrorCode ierr;
971 
972   PetscFunctionBegin;
973   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
974   pc->data = (void*)def;
975 
976   def->init          = PETSC_FALSE;
977   def->correct       = PETSC_FALSE;
978   def->correctfact   = 1.0;
979   def->reductionfact = -1;
980   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
981   def->spacesize     = 1;
982   def->extendsp      = PETSC_FALSE;
983   def->lvl           = 0;
984   def->maxlvl        = 0;
985 
986   pc->ops->apply          = PCApply_Deflation;
987   pc->ops->presolve       = PCPreSolve_Deflation;
988   pc->ops->setup          = PCSetUp_Deflation;
989   pc->ops->reset          = PCReset_Deflation;
990   pc->ops->destroy        = PCDestroy_Deflation;
991   pc->ops->setfromoptions = PCSetFromOptions_Deflation;
992   pc->ops->view           = PCView_Deflation;
993 
994   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetInitOnly_C",PCDeflationSetInitOnly_Deflation);CHKERRQ(ierr);
995   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
996   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetReductionFactor_C",PCDeflationSetReductionFactor_Deflation);CHKERRQ(ierr);
997   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCorrectionFactor_C",PCDeflationSetCorrectionFactor_Deflation);CHKERRQ(ierr);
998   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpaceToCompute_C",PCDeflationSetSpaceToCompute_Deflation);CHKERRQ(ierr);
999   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
1000   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetProjectionNullSpaceMat_C",PCDeflationSetProjectionNullSpaceMat_Deflation);CHKERRQ(ierr);
1001   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseMat_C",PCDeflationSetCoarseMat_Deflation);CHKERRQ(ierr);
1002   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
1003   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetCoarseKSP_C",PCDeflationSetCoarseKSP_Deflation);CHKERRQ(ierr);
1004   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetPC_C",PCDeflationGetPC_Deflation);CHKERRQ(ierr);
1005   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
1006   PetscFunctionReturn(0);
1007 }
1008 
1009