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