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