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