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