xref: /petsc/src/ksp/pc/interface/precon.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1 /*
2     The PC (preconditioner) interface routines, callable by users.
3 */
4 #include "src/ksp/pc/pcimpl.h"            /*I "petscksp.h" I*/
5 
6 /* Logging support */
7 PetscCookie PC_COOKIE = 0;
8 PetscEvent  PC_SetUp = 0, PC_SetUpOnBlocks = 0, PC_Apply = 0, PC_ApplyCoarse = 0, PC_ApplyMultiple = 0, PC_ApplySymmetricLeft = 0;
9 PetscEvent  PC_ApplySymmetricRight = 0, PC_ModifySubMatrices = 0;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "PCGetDefaultType_Private"
13 PetscErrorCode PCGetDefaultType_Private(PC pc,const char* type[])
14 {
15   PetscErrorCode ierr;
16   PetscMPIInt    size;
17   PetscTruth     flg1,flg2,set,flg3;
18 
19   PetscFunctionBegin;
20   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
21   if (pc->pmat) {
22     PetscErrorCode (*f)(Mat,PetscTruth*,MatReuse,Mat*);
23     ierr = PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
24     if (size == 1) {
25       ierr = MatHasOperation(pc->pmat,MATOP_ICCFACTOR_SYMBOLIC,&flg1);CHKERRQ(ierr);
26       ierr = MatHasOperation(pc->pmat,MATOP_ILUFACTOR_SYMBOLIC,&flg2);CHKERRQ(ierr);
27       ierr = MatIsSymmetricKnown(pc->pmat,&set,&flg3);CHKERRQ(ierr);
28       if (flg1 && (!flg2 || (set && flg3))) {
29 	*type = PCICC;
30       } else if (flg2) {
31 	*type = PCILU;
32       } else if (f) { /* likely is a parallel matrix run on one processor */
33 	*type = PCBJACOBI;
34       } else {
35 	*type = PCNONE;
36       }
37     } else {
38        if (f) {
39 	*type = PCBJACOBI;
40       } else {
41 	*type = PCNONE;
42       }
43     }
44   } else {
45     if (size == 1) {
46       *type = PCILU;
47     } else {
48       *type = PCBJACOBI;
49     }
50   }
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "PCDestroy"
56 /*@C
57    PCDestroy - Destroys PC context that was created with PCCreate().
58 
59    Collective on PC
60 
61    Input Parameter:
62 .  pc - the preconditioner context
63 
64    Level: developer
65 
66 .keywords: PC, destroy
67 
68 .seealso: PCCreate(), PCSetUp()
69 @*/
70 PetscErrorCode PCDestroy(PC pc)
71 {
72   PetscErrorCode ierr;
73 
74   PetscFunctionBegin;
75   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
76   if (--pc->refct > 0) PetscFunctionReturn(0);
77 
78   /* if memory was published with AMS then destroy it */
79   ierr = PetscObjectDepublish(pc);CHKERRQ(ierr);
80 
81   if (pc->ops->destroy)       {ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);}
82   if (pc->diagonalscaleright) {ierr = VecDestroy(pc->diagonalscaleright);CHKERRQ(ierr);}
83   if (pc->diagonalscaleleft)  {ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);}
84 
85   ierr = PetscHeaderDestroy(pc);CHKERRQ(ierr);
86   PetscFunctionReturn(0);
87 }
88 
89 #undef __FUNCT__
90 #define __FUNCT__ "PCDiagonalScale"
91 /*@C
92    PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
93       scaling as needed by certain time-stepping codes.
94 
95    Collective on PC
96 
97    Input Parameter:
98 .  pc - the preconditioner context
99 
100    Output Parameter:
101 .  flag - PETSC_TRUE if it applies the scaling
102 
103    Level: developer
104 
105    Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
106 $           D M A D^{-1} y = D M b  for left preconditioning or
107 $           D A M D^{-1} z = D b for right preconditioning
108 
109 .keywords: PC
110 
111 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
112 @*/
113 PetscErrorCode PCDiagonalScale(PC pc,PetscTruth *flag)
114 {
115   PetscFunctionBegin;
116   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
117   PetscValidPointer(flag,2);
118   *flag = pc->diagonalscale;
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "PCDiagonalScaleSet"
124 /*@
125    PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
126       scaling as needed by certain time-stepping codes.
127 
128    Collective on PC
129 
130    Input Parameters:
131 +  pc - the preconditioner context
132 -  s - scaling vector
133 
134    Level: intermediate
135 
136    Notes: The system solved via the Krylov method is
137 $           D M A D^{-1} y = D M b  for left preconditioning or
138 $           D A M D^{-1} z = D b for right preconditioning
139 
140    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
141 
142 .keywords: PC
143 
144 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
145 @*/
146 PetscErrorCode PCDiagonalScaleSet(PC pc,Vec s)
147 {
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
152   PetscValidHeaderSpecific(s,VEC_COOKIE,2);
153   pc->diagonalscale     = PETSC_TRUE;
154   if (pc->diagonalscaleleft) {
155     ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);
156   }
157   pc->diagonalscaleleft = s;
158   ierr                  = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
159   if (!pc->diagonalscaleright) {
160     ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr);
161   }
162   ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr);
163   ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "PCDiagonalScaleLeft"
169 /*@C
170    PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
171       scaling as needed by certain time-stepping codes.
172 
173    Collective on PC
174 
175    Input Parameters:
176 +  pc - the preconditioner context
177 .  in - input vector
178 +  out - scaled vector (maybe the same as in)
179 
180    Level: intermediate
181 
182    Notes: The system solved via the Krylov method is
183 $           D M A D^{-1} y = D M b  for left preconditioning or
184 $           D A M D^{-1} z = D b for right preconditioning
185 
186    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
187 
188    If diagonal scaling is turned off and in is not out then in is copied to out
189 
190 .keywords: PC
191 
192 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
193 @*/
194 PetscErrorCode PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
195 {
196   PetscErrorCode ierr;
197 
198   PetscFunctionBegin;
199   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
200   PetscValidHeaderSpecific(in,VEC_COOKIE,2);
201   PetscValidHeaderSpecific(out,VEC_COOKIE,3);
202   if (pc->diagonalscale) {
203     ierr = VecPointwiseMult(pc->diagonalscaleleft,in,out);CHKERRQ(ierr);
204   } else if (in != out) {
205     ierr = VecCopy(in,out);CHKERRQ(ierr);
206   }
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "PCDiagonalScaleRight"
212 /*@C
213    PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
214 
215    Collective on PC
216 
217    Input Parameters:
218 +  pc - the preconditioner context
219 .  in - input vector
220 +  out - scaled vector (maybe the same as in)
221 
222    Level: intermediate
223 
224    Notes: The system solved via the Krylov method is
225 $           D M A D^{-1} y = D M b  for left preconditioning or
226 $           D A M D^{-1} z = D b for right preconditioning
227 
228    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
229 
230    If diagonal scaling is turned off and in is not out then in is copied to out
231 
232 .keywords: PC
233 
234 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
235 @*/
236 PetscErrorCode PCDiagonalScaleRight(PC pc,Vec in,Vec out)
237 {
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
242   PetscValidHeaderSpecific(in,VEC_COOKIE,2);
243   PetscValidHeaderSpecific(out,VEC_COOKIE,3);
244   if (pc->diagonalscale) {
245     ierr = VecPointwiseMult(pc->diagonalscaleright,in,out);CHKERRQ(ierr);
246   } else if (in != out) {
247     ierr = VecCopy(in,out);CHKERRQ(ierr);
248   }
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "PCPublish_Petsc"
254 static PetscErrorCode PCPublish_Petsc(PetscObject obj)
255 {
256   PetscFunctionBegin;
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "PCCreate"
262 /*@C
263    PCCreate - Creates a preconditioner context.
264 
265    Collective on MPI_Comm
266 
267    Input Parameter:
268 .  comm - MPI communicator
269 
270    Output Parameter:
271 .  pc - location to put the preconditioner context
272 
273    Notes:
274    The default preconditioner on one processor is PCILU with 0 fill on more
275    then one it is PCBJACOBI with ILU() on each processor.
276 
277    Level: developer
278 
279 .keywords: PC, create, context
280 
281 .seealso: PCSetUp(), PCApply(), PCDestroy()
282 @*/
283 PetscErrorCode PCCreate(MPI_Comm comm,PC *newpc)
284 {
285   PC             pc;
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   PetscValidPointer(newpc,1)
290   *newpc = 0;
291 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
292   ierr = PCInitializePackage(PETSC_NULL);CHKERRQ(ierr);
293 #endif
294 
295   ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
296   pc->bops->publish      = PCPublish_Petsc;
297   pc->mat                = 0;
298   pc->setupcalled        = 0;
299   pc->data               = 0;
300   pc->diagonalscale      = PETSC_FALSE;
301   pc->diagonalscaleleft  = 0;
302   pc->diagonalscaleright = 0;
303 
304   pc->ops->destroy             = 0;
305   pc->ops->apply               = 0;
306   pc->ops->applytranspose      = 0;
307   pc->ops->applyBA             = 0;
308   pc->ops->applyBAtranspose    = 0;
309   pc->ops->applyrichardson     = 0;
310   pc->ops->view                = 0;
311   pc->ops->getfactoredmatrix   = 0;
312   pc->ops->applysymmetricright = 0;
313   pc->ops->applysymmetricleft  = 0;
314   pc->ops->setuponblocks       = 0;
315 
316   pc->modifysubmatrices   = 0;
317   pc->modifysubmatricesP  = 0;
318   *newpc                  = pc;
319   ierr = PetscPublishAll(pc);CHKERRQ(ierr);
320   PetscFunctionReturn(0);
321 
322 }
323 
324 /* -------------------------------------------------------------------------------*/
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "PCApply"
328 /*@
329    PCApply - Applies the preconditioner to a vector.
330 
331    Collective on PC and Vec
332 
333    Input Parameters:
334 +  pc - the preconditioner context
335 -  x - input vector
336 
337    Output Parameter:
338 .  y - output vector
339 
340    Level: developer
341 
342 .keywords: PC, apply
343 
344 .seealso: PCApplyTranspose(), PCApplyBAorAB()
345 @*/
346 PetscErrorCode PCApply(PC pc,Vec x,Vec y)
347 {
348   PetscErrorCode ierr;
349 
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
352   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
353   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
354   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
355 
356   if (pc->setupcalled < 2) {
357     ierr = PCSetUp(pc);CHKERRQ(ierr);
358   }
359   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
360   ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr);
361   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "PCApplySymmetricLeft"
367 /*@
368    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
369 
370    Collective on PC and Vec
371 
372    Input Parameters:
373 +  pc - the preconditioner context
374 -  x - input vector
375 
376    Output Parameter:
377 .  y - output vector
378 
379    Notes:
380    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
381 
382    Level: developer
383 
384 .keywords: PC, apply, symmetric, left
385 
386 .seealso: PCApply(), PCApplySymmetricRight()
387 @*/
388 PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y)
389 {
390   PetscErrorCode ierr;
391 
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
394   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
395   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
396   if (!pc->ops->applysymmetricleft) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
397 
398   if (pc->setupcalled < 2) {
399     ierr = PCSetUp(pc);CHKERRQ(ierr);
400   }
401 
402   ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
403   ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr);
404   ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCApplySymmetricRight"
410 /*@
411    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
412 
413    Collective on PC and Vec
414 
415    Input Parameters:
416 +  pc - the preconditioner context
417 -  x - input vector
418 
419    Output Parameter:
420 .  y - output vector
421 
422    Level: developer
423 
424    Notes:
425    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
426 
427 .keywords: PC, apply, symmetric, right
428 
429 .seealso: PCApply(), PCApplySymmetricLeft()
430 @*/
431 PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y)
432 {
433   PetscErrorCode ierr;
434 
435   PetscFunctionBegin;
436   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
437   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
438   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
439   if (!pc->ops->applysymmetricright) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
440 
441   if (pc->setupcalled < 2) {
442     ierr = PCSetUp(pc);CHKERRQ(ierr);
443   }
444 
445   ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
446   ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr);
447   ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "PCApplyTranspose"
453 /*@
454    PCApplyTranspose - Applies the transpose of preconditioner to a vector.
455 
456    Collective on PC and Vec
457 
458    Input Parameters:
459 +  pc - the preconditioner context
460 -  x - input vector
461 
462    Output Parameter:
463 .  y - output vector
464 
465    Level: developer
466 
467 .keywords: PC, apply, transpose
468 
469 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCHasApplyTranspose()
470 @*/
471 PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y)
472 {
473   PetscErrorCode ierr;
474 
475   PetscFunctionBegin;
476   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
477   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
478   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
479   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
480   if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP," ");
481 
482   if (pc->setupcalled < 2) {
483     ierr = PCSetUp(pc);CHKERRQ(ierr);
484   }
485 
486   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
487   ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr);
488   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 #undef __FUNCT__
493 #define __FUNCT__ "PCHasApplyTranspose"
494 /*@
495    PCHasApplyTranspose - Test whether the preconditioner has a transpose apply operation
496 
497    Collective on PC and Vec
498 
499    Input Parameters:
500 .  pc - the preconditioner context
501 
502    Output Parameter:
503 .  flg - PETSC_TRUE if a transpose operation is defined
504 
505    Level: developer
506 
507 .keywords: PC, apply, transpose
508 
509 .seealso: PCApplyTranspose()
510 @*/
511 PetscErrorCode PCHasApplyTranspose(PC pc,PetscTruth *flg)
512 {
513   PetscFunctionBegin;
514   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
515   PetscValidPointer(flg,2);
516   *flg = (PetscTruth) (pc->ops->applytranspose == 0);
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "PCApplyBAorAB"
522 /*@
523    PCApplyBAorAB - Applies the preconditioner and operator to a vector.
524 
525    Collective on PC and Vec
526 
527    Input Parameters:
528 +  pc - the preconditioner context
529 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
530 .  x - input vector
531 -  work - work vector
532 
533    Output Parameter:
534 .  y - output vector
535 
536    Level: developer
537 
538 .keywords: PC, apply, operator
539 
540 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
541 @*/
542 PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
543 {
544   PetscErrorCode ierr;
545 
546   PetscFunctionBegin;
547   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
548   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
549   PetscValidHeaderSpecific(y,VEC_COOKIE,4);
550   PetscValidHeaderSpecific(work,VEC_COOKIE,5);
551   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
552   if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
553     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
554   }
555   if (pc->diagonalscale && side == PC_SYMMETRIC) {
556     SETERRQ(PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
557   }
558 
559   if (pc->setupcalled < 2) {
560     ierr = PCSetUp(pc);CHKERRQ(ierr);
561   }
562 
563   if (pc->diagonalscale) {
564     if (pc->ops->applyBA) {
565       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
566       ierr = VecDuplicate(x,&work2);CHKERRQ(ierr);
567       ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr);
568       ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr);
569       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
570       ierr = VecDestroy(work2);CHKERRQ(ierr);
571     } else if (side == PC_RIGHT) {
572       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
573       ierr = PCApply(pc,y,work);CHKERRQ(ierr);
574       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
575       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
576     } else if (side == PC_LEFT) {
577       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
578       ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr);
579       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
580       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
581     } else if (side == PC_SYMMETRIC) {
582       SETERRQ(PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
583     }
584   } else {
585     if (pc->ops->applyBA) {
586       ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr);
587     } else if (side == PC_RIGHT) {
588       ierr = PCApply(pc,x,work);CHKERRQ(ierr);
589       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
590     } else if (side == PC_LEFT) {
591       ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr);
592       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
593     } else if (side == PC_SYMMETRIC) {
594       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
595       ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr);
596       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
597       ierr = VecCopy(y,work);CHKERRQ(ierr);
598       ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr);
599     }
600   }
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "PCApplyBAorABTranspose"
606 /*@
607    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
608    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
609    not tr(B*A) = tr(A)*tr(B).
610 
611    Collective on PC and Vec
612 
613    Input Parameters:
614 +  pc - the preconditioner context
615 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
616 .  x - input vector
617 -  work - work vector
618 
619    Output Parameter:
620 .  y - output vector
621 
622    Level: developer
623 
624 .keywords: PC, apply, operator, transpose
625 
626 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
627 @*/
628 PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
629 {
630   PetscErrorCode ierr;
631 
632   PetscFunctionBegin;
633   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
634   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
635   PetscValidHeaderSpecific(y,VEC_COOKIE,4);
636   PetscValidHeaderSpecific(work,VEC_COOKIE,5);
637   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
638   if (pc->ops->applyBAtranspose) {
639     ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr);
640     PetscFunctionReturn(0);
641   }
642   if (side != PC_LEFT && side != PC_RIGHT) {
643     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
644   }
645 
646   if (pc->setupcalled < 2) {
647     ierr = PCSetUp(pc);CHKERRQ(ierr);
648   }
649 
650   if (side == PC_RIGHT) {
651     ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr);
652     ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr);
653   } else if (side == PC_LEFT) {
654     ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr);
655     ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr);
656   }
657   /* add support for PC_SYMMETRIC */
658   PetscFunctionReturn(0); /* actually will never get here */
659 }
660 
661 /* -------------------------------------------------------------------------------*/
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "PCApplyRichardsonExists"
665 /*@
666    PCApplyRichardsonExists - Determines whether a particular preconditioner has a
667    built-in fast application of Richardson's method.
668 
669    Not Collective
670 
671    Input Parameter:
672 .  pc - the preconditioner
673 
674    Output Parameter:
675 .  exists - PETSC_TRUE or PETSC_FALSE
676 
677    Level: developer
678 
679 .keywords: PC, apply, Richardson, exists
680 
681 .seealso: PCApplyRichardson()
682 @*/
683 PetscErrorCode PCApplyRichardsonExists(PC pc,PetscTruth *exists)
684 {
685   PetscFunctionBegin;
686   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
687   PetscValidIntPointer(exists,2);
688   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
689   else                    *exists = PETSC_FALSE;
690   PetscFunctionReturn(0);
691 }
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "PCApplyRichardson"
695 /*@
696    PCApplyRichardson - Applies several steps of Richardson iteration with
697    the particular preconditioner. This routine is usually used by the
698    Krylov solvers and not the application code directly.
699 
700    Collective on PC
701 
702    Input Parameters:
703 +  pc  - the preconditioner context
704 .  x   - the initial guess
705 .  w   - one work vector
706 .  rtol - relative decrease in residual norm convergence criteria
707 .  abstol - absolute residual norm convergence criteria
708 .  dtol - divergence residual norm increase criteria
709 -  its - the number of iterations to apply.
710 
711    Output Parameter:
712 .  y - the solution
713 
714    Notes:
715    Most preconditioners do not support this function. Use the command
716    PCApplyRichardsonExists() to determine if one does.
717 
718    Except for the multigrid PC this routine ignores the convergence tolerances
719    and always runs for the number of iterations
720 
721    Level: developer
722 
723 .keywords: PC, apply, Richardson
724 
725 .seealso: PCApplyRichardsonExists()
726 @*/
727 PetscErrorCode PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
728 {
729   PetscErrorCode ierr;
730 
731   PetscFunctionBegin;
732   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
733   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
734   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
735   PetscValidHeaderSpecific(w,VEC_COOKIE,4);
736   if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP," ");
737 
738   if (pc->setupcalled < 2) {
739     ierr = PCSetUp(pc);CHKERRQ(ierr);
740   }
741 
742   ierr = (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its);CHKERRQ(ierr);
743   PetscFunctionReturn(0);
744 }
745 
746 /*
747       a setupcall of 0 indicates never setup,
748                      1 needs to be resetup,
749                      2 does not need any changes.
750 */
751 #undef __FUNCT__
752 #define __FUNCT__ "PCSetUp"
753 /*@
754    PCSetUp - Prepares for the use of a preconditioner.
755 
756    Collective on PC
757 
758    Input Parameter:
759 .  pc - the preconditioner context
760 
761    Level: developer
762 
763 .keywords: PC, setup
764 
765 .seealso: PCCreate(), PCApply(), PCDestroy()
766 @*/
767 PetscErrorCode PCSetUp(PC pc)
768 {
769   PetscErrorCode ierr;
770   const char     *def;
771 
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
774 
775   if (pc->setupcalled > 1) {
776     ierr = PetscLogInfo((pc,"PCSetUp:Setting PC with identical preconditioner\n"));CHKERRQ(ierr);
777     PetscFunctionReturn(0);
778   } else if (!pc->setupcalled) {
779     ierr = PetscLogInfo((pc,"PCSetUp:Setting up new PC\n"));CHKERRQ(ierr);
780   } else if (pc->flag == SAME_NONZERO_PATTERN) {
781     ierr = PetscLogInfo((pc,"PCSetUp:Setting up PC with same nonzero pattern\n"));CHKERRQ(ierr);
782   } else {
783     ierr = PetscLogInfo((pc,"PCSetUp:Setting up PC with different nonzero pattern\n"));CHKERRQ(ierr);
784   }
785 
786   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
787   if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
788 
789   if (!pc->type_name) {
790     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
791     ierr = PCSetType(pc,def);CHKERRQ(ierr);
792   }
793 
794   if (pc->ops->setup) {
795     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
796   }
797   pc->setupcalled = 2;
798   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "PCSetUpOnBlocks"
804 /*@
805    PCSetUpOnBlocks - Sets up the preconditioner for each block in
806    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
807    methods.
808 
809    Collective on PC
810 
811    Input Parameters:
812 .  pc - the preconditioner context
813 
814    Level: developer
815 
816 .keywords: PC, setup, blocks
817 
818 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
819 @*/
820 PetscErrorCode PCSetUpOnBlocks(PC pc)
821 {
822   PetscErrorCode ierr;
823 
824   PetscFunctionBegin;
825   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
826   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
827   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
828   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
829   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
830   PetscFunctionReturn(0);
831 }
832 
833 #undef __FUNCT__
834 #define __FUNCT__ "PCSetModifySubMatrices"
835 /*@C
836    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
837    submatrices that arise within certain subdomain-based preconditioners.
838    The basic submatrices are extracted from the preconditioner matrix as
839    usual; the user can then alter these (for example, to set different boundary
840    conditions for each submatrix) before they are used for the local solves.
841 
842    Collective on PC
843 
844    Input Parameters:
845 +  pc - the preconditioner context
846 .  func - routine for modifying the submatrices
847 -  ctx - optional user-defined context (may be null)
848 
849    Calling sequence of func:
850 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
851 
852 .  row - an array of index sets that contain the global row numbers
853          that comprise each local submatrix
854 .  col - an array of index sets that contain the global column numbers
855          that comprise each local submatrix
856 .  submat - array of local submatrices
857 -  ctx - optional user-defined context for private data for the
858          user-defined func routine (may be null)
859 
860    Notes:
861    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
862    KSPSolve().
863 
864    A routine set by PCSetModifySubMatrices() is currently called within
865    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
866    preconditioners.  All other preconditioners ignore this routine.
867 
868    Level: advanced
869 
870 .keywords: PC, set, modify, submatrices
871 
872 .seealso: PCModifySubMatrices()
873 @*/
874 PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
875 {
876   PetscFunctionBegin;
877   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
878   pc->modifysubmatrices  = func;
879   pc->modifysubmatricesP = ctx;
880   PetscFunctionReturn(0);
881 }
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "PCModifySubMatrices"
885 /*@C
886    PCModifySubMatrices - Calls an optional user-defined routine within
887    certain preconditioners if one has been set with PCSetModifySubMarices().
888 
889    Collective on PC
890 
891    Input Parameters:
892 +  pc - the preconditioner context
893 .  nsub - the number of local submatrices
894 .  row - an array of index sets that contain the global row numbers
895          that comprise each local submatrix
896 .  col - an array of index sets that contain the global column numbers
897          that comprise each local submatrix
898 .  submat - array of local submatrices
899 -  ctx - optional user-defined context for private data for the
900          user-defined routine (may be null)
901 
902    Output Parameter:
903 .  submat - array of local submatrices (the entries of which may
904             have been modified)
905 
906    Notes:
907    The user should NOT generally call this routine, as it will
908    automatically be called within certain preconditioners (currently
909    block Jacobi, additive Schwarz) if set.
910 
911    The basic submatrices are extracted from the preconditioner matrix
912    as usual; the user can then alter these (for example, to set different
913    boundary conditions for each submatrix) before they are used for the
914    local solves.
915 
916    Level: developer
917 
918 .keywords: PC, modify, submatrices
919 
920 .seealso: PCSetModifySubMatrices()
921 @*/
922 PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
923 {
924   PetscErrorCode ierr;
925 
926   PetscFunctionBegin;
927   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
928   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
929   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
930   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
931   PetscFunctionReturn(0);
932 }
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "PCSetOperators"
936 /*@
937    PCSetOperators - Sets the matrix associated with the linear system and
938    a (possibly) different one associated with the preconditioner.
939 
940    Collective on PC and Mat
941 
942    Input Parameters:
943 +  pc - the preconditioner context
944 .  Amat - the matrix associated with the linear system
945 .  Pmat - the matrix to be used in constructing the preconditioner, usually
946           the same as Amat.
947 -  flag - flag indicating information about the preconditioner matrix structure
948    during successive linear solves.  This flag is ignored the first time a
949    linear system is solved, and thus is irrelevant when solving just one linear
950    system.
951 
952    Notes:
953    The flag can be used to eliminate unnecessary work in the preconditioner
954    during the repeated solution of linear systems of the same size.  The
955    available options are
956 +    SAME_PRECONDITIONER -
957        Pmat is identical during successive linear solves.
958        This option is intended for folks who are using
959        different Amat and Pmat matrices and wish to reuse the
960        same preconditioner matrix.  For example, this option
961        saves work by not recomputing incomplete factorization
962        for ILU/ICC preconditioners.
963 .     SAME_NONZERO_PATTERN -
964        Pmat has the same nonzero structure during
965        successive linear solves.
966 -     DIFFERENT_NONZERO_PATTERN -
967        Pmat does not have the same nonzero structure.
968 
969    Caution:
970    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
971    and does not check the structure of the matrix.  If you erroneously
972    claim that the structure is the same when it actually is not, the new
973    preconditioner will not function correctly.  Thus, use this optimization
974    feature carefully!
975 
976    If in doubt about whether your preconditioner matrix has changed
977    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
978 
979    More Notes about Repeated Solution of Linear Systems:
980    PETSc does NOT reset the matrix entries of either Amat or Pmat
981    to zero after a linear solve; the user is completely responsible for
982    matrix assembly.  See the routine MatZeroEntries() if desiring to
983    zero all elements of a matrix.
984 
985    Level: intermediate
986 
987 .keywords: PC, set, operators, matrix, linear system
988 
989 .seealso: PCGetOperators(), MatZeroEntries()
990  @*/
991 PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
992 {
993   PetscErrorCode ierr;
994   PetscTruth     isbjacobi,isrowbs;
995 
996   PetscFunctionBegin;
997   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
998   if (Amat) PetscValidHeaderSpecific(Amat,MAT_COOKIE,2);
999   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_COOKIE,3);
1000 
1001   /*
1002       BlockSolve95 cannot use default BJacobi preconditioning
1003   */
1004   if (Amat) {
1005     ierr = PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);CHKERRQ(ierr);
1006     if (isrowbs) {
1007       ierr = PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);CHKERRQ(ierr);
1008       if (isbjacobi) {
1009         ierr = PCSetType(pc,PCILU);CHKERRQ(ierr);
1010         ierr = PetscLogInfo((pc,"PCSetOperators:Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n"));CHKERRQ(ierr);
1011       }
1012     }
1013   }
1014 
1015   pc->mat  = Amat;
1016   pc->pmat = Pmat;
1017   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1018     pc->setupcalled = 1;
1019   }
1020   pc->flag = flag;
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 #undef __FUNCT__
1025 #define __FUNCT__ "PCGetOperators"
1026 /*@C
1027    PCGetOperators - Gets the matrix associated with the linear system and
1028    possibly a different one associated with the preconditioner.
1029 
1030    Not collective, though parallel Mats are returned if the PC is parallel
1031 
1032    Input Parameter:
1033 .  pc - the preconditioner context
1034 
1035    Output Parameters:
1036 +  mat - the matrix associated with the linear system
1037 .  pmat - matrix associated with the preconditioner, usually the same
1038           as mat.
1039 -  flag - flag indicating information about the preconditioner
1040           matrix structure.  See PCSetOperators() for details.
1041 
1042    Level: intermediate
1043 
1044 .keywords: PC, get, operators, matrix, linear system
1045 
1046 .seealso: PCSetOperators()
1047 @*/
1048 PetscErrorCode PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1049 {
1050   PetscFunctionBegin;
1051   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1052   if (mat)  *mat  = pc->mat;
1053   if (pmat) *pmat = pc->pmat;
1054   if (flag) *flag = pc->flag;
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "PCGetFactoredMatrix"
1060 /*@C
1061    PCGetFactoredMatrix - Gets the factored matrix from the
1062    preconditioner context.  This routine is valid only for the LU,
1063    incomplete LU, Cholesky, and incomplete Cholesky methods.
1064 
1065    Not Collective on PC though Mat is parallel if PC is parallel
1066 
1067    Input Parameters:
1068 .  pc - the preconditioner context
1069 
1070    Output parameters:
1071 .  mat - the factored matrix
1072 
1073    Level: advanced
1074 
1075    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1076 
1077 .keywords: PC, get, factored, matrix
1078 @*/
1079 PetscErrorCode PCGetFactoredMatrix(PC pc,Mat *mat)
1080 {
1081   PetscErrorCode ierr;
1082 
1083   PetscFunctionBegin;
1084   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1085   PetscValidPointer(mat,2);
1086   if (pc->ops->getfactoredmatrix) {
1087     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1088   }
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "PCSetOptionsPrefix"
1094 /*@C
1095    PCSetOptionsPrefix - Sets the prefix used for searching for all
1096    PC options in the database.
1097 
1098    Collective on PC
1099 
1100    Input Parameters:
1101 +  pc - the preconditioner context
1102 -  prefix - the prefix string to prepend to all PC option requests
1103 
1104    Notes:
1105    A hyphen (-) must NOT be given at the beginning of the prefix name.
1106    The first character of all runtime options is AUTOMATICALLY the
1107    hyphen.
1108 
1109    Level: advanced
1110 
1111 .keywords: PC, set, options, prefix, database
1112 
1113 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1114 @*/
1115 PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[])
1116 {
1117   PetscErrorCode ierr;
1118 
1119   PetscFunctionBegin;
1120   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1121   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "PCAppendOptionsPrefix"
1127 /*@C
1128    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1129    PC options in the database.
1130 
1131    Collective on PC
1132 
1133    Input Parameters:
1134 +  pc - the preconditioner context
1135 -  prefix - the prefix string to prepend to all PC option requests
1136 
1137    Notes:
1138    A hyphen (-) must NOT be given at the beginning of the prefix name.
1139    The first character of all runtime options is AUTOMATICALLY the
1140    hyphen.
1141 
1142    Level: advanced
1143 
1144 .keywords: PC, append, options, prefix, database
1145 
1146 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1147 @*/
1148 PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[])
1149 {
1150   PetscErrorCode ierr;
1151 
1152   PetscFunctionBegin;
1153   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1154   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1155   PetscFunctionReturn(0);
1156 }
1157 
1158 #undef __FUNCT__
1159 #define __FUNCT__ "PCGetOptionsPrefix"
1160 /*@C
1161    PCGetOptionsPrefix - Gets the prefix used for searching for all
1162    PC options in the database.
1163 
1164    Not Collective
1165 
1166    Input Parameters:
1167 .  pc - the preconditioner context
1168 
1169    Output Parameters:
1170 .  prefix - pointer to the prefix string used, is returned
1171 
1172    Notes: On the fortran side, the user should pass in a string 'prifix' of
1173    sufficient length to hold the prefix.
1174 
1175    Level: advanced
1176 
1177 .keywords: PC, get, options, prefix, database
1178 
1179 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1180 @*/
1181 PetscErrorCode PCGetOptionsPrefix(PC pc,char *prefix[])
1182 {
1183   PetscErrorCode ierr;
1184 
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1187   PetscValidPointer(prefix,2);
1188   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 #undef __FUNCT__
1193 #define __FUNCT__ "PCPreSolve"
1194 /*@
1195    PCPreSolve - Optional pre-solve phase, intended for any
1196    preconditioner-specific actions that must be performed before
1197    the iterative solve itself.
1198 
1199    Collective on PC
1200 
1201    Input Parameters:
1202 +  pc - the preconditioner context
1203 -  ksp - the Krylov subspace context
1204 
1205    Level: developer
1206 
1207    Sample of Usage:
1208 .vb
1209     PCPreSolve(pc,ksp);
1210     KSPSolve(ksp,b,x);
1211     PCPostSolve(pc,ksp);
1212 .ve
1213 
1214    Notes:
1215    The pre-solve phase is distinct from the PCSetUp() phase.
1216 
1217    KSPSolve() calls this directly, so is rarely called by the user.
1218 
1219 .keywords: PC, pre-solve
1220 
1221 .seealso: PCPostSolve()
1222 @*/
1223 PetscErrorCode PCPreSolve(PC pc,KSP ksp)
1224 {
1225   PetscErrorCode ierr;
1226   Vec            x,rhs;
1227   Mat            A,B;
1228 
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1231   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
1232   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1233   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1234   /*
1235       Scale the system and have the matrices use the scaled form
1236     only if the two matrices are actually the same (and hence
1237     have the same scaling
1238   */
1239   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1240   if (A == B) {
1241     ierr = MatScaleSystem(pc->mat,x,rhs);CHKERRQ(ierr);
1242     ierr = MatUseScaledForm(pc->mat,PETSC_TRUE);CHKERRQ(ierr);
1243   }
1244 
1245   if (pc->ops->presolve) {
1246     ierr = (*pc->ops->presolve)(pc,ksp,x,rhs);CHKERRQ(ierr);
1247   }
1248   PetscFunctionReturn(0);
1249 }
1250 
1251 #undef __FUNCT__
1252 #define __FUNCT__ "PCPostSolve"
1253 /*@
1254    PCPostSolve - Optional post-solve phase, intended for any
1255    preconditioner-specific actions that must be performed after
1256    the iterative solve itself.
1257 
1258    Collective on PC
1259 
1260    Input Parameters:
1261 +  pc - the preconditioner context
1262 -  ksp - the Krylov subspace context
1263 
1264    Sample of Usage:
1265 .vb
1266     PCPreSolve(pc,ksp);
1267     KSPSolve(ksp,b,x);
1268     PCPostSolve(pc,ksp);
1269 .ve
1270 
1271    Note:
1272    KSPSolve() calls this routine directly, so it is rarely called by the user.
1273 
1274    Level: developer
1275 
1276 .keywords: PC, post-solve
1277 
1278 .seealso: PCPreSolve(), KSPSolve()
1279 @*/
1280 PetscErrorCode PCPostSolve(PC pc,KSP ksp)
1281 {
1282   PetscErrorCode ierr;
1283   Vec            x,rhs;
1284   Mat            A,B;
1285 
1286   PetscFunctionBegin;
1287   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1288   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
1289   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1290   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1291   if (pc->ops->postsolve) {
1292     ierr =  (*pc->ops->postsolve)(pc,ksp,x,rhs);CHKERRQ(ierr);
1293   }
1294 
1295   /*
1296       Scale the system and have the matrices use the scaled form
1297     only if the two matrices are actually the same (and hence
1298     have the same scaling
1299   */
1300   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1301   if (A == B) {
1302     ierr = MatUnScaleSystem(pc->mat,x,rhs);CHKERRQ(ierr);
1303     ierr = MatUseScaledForm(pc->mat,PETSC_FALSE);CHKERRQ(ierr);
1304   }
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 #undef __FUNCT__
1309 #define __FUNCT__ "PCView"
1310 /*@C
1311    PCView - Prints the PC data structure.
1312 
1313    Collective on PC
1314 
1315    Input Parameters:
1316 +  PC - the PC context
1317 -  viewer - optional visualization context
1318 
1319    Note:
1320    The available visualization contexts include
1321 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1322 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1323          output where only the first processor opens
1324          the file.  All other processors send their
1325          data to the first processor to print.
1326 
1327    The user can open an alternative visualization contexts with
1328    PetscViewerASCIIOpen() (output to a specified file).
1329 
1330    Level: developer
1331 
1332 .keywords: PC, view
1333 
1334 .seealso: KSPView(), PetscViewerASCIIOpen()
1335 @*/
1336 PetscErrorCode PCView(PC pc,PetscViewer viewer)
1337 {
1338   PCType            cstr;
1339   PetscErrorCode    ierr;
1340   PetscTruth        mat_exists,iascii,isstring;
1341   PetscViewerFormat format;
1342 
1343   PetscFunctionBegin;
1344   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1345   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(pc->comm);
1346   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
1347   PetscCheckSameComm(pc,1,viewer,2);
1348 
1349   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1350   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1351   if (iascii) {
1352     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1353     if (pc->prefix) {
1354       ierr = PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",pc->prefix);CHKERRQ(ierr);
1355     } else {
1356       ierr = PetscViewerASCIIPrintf(viewer,"PC Object:\n");CHKERRQ(ierr);
1357     }
1358     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1359     if (cstr) {
1360       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",cstr);CHKERRQ(ierr);
1361     } else {
1362       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
1363     }
1364     if (pc->ops->view) {
1365       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1366       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1367       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1368     }
1369     ierr = PetscObjectExists((PetscObject)pc->mat,&mat_exists);CHKERRQ(ierr);
1370     if (mat_exists) {
1371       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1372       if (pc->pmat == pc->mat) {
1373         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1374         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1375         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1376         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1377       } else {
1378         ierr = PetscObjectExists((PetscObject)pc->pmat,&mat_exists);CHKERRQ(ierr);
1379         if (mat_exists) {
1380           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1381         } else {
1382           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1383         }
1384         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1385         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1386         if (mat_exists) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1387         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1388       }
1389       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1390     }
1391   } else if (isstring) {
1392     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1393     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1394     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1395   } else {
1396     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1397   }
1398   PetscFunctionReturn(0);
1399 }
1400 
1401 #undef __FUNCT__
1402 #define __FUNCT__ "PCRegister"
1403 /*@C
1404   PCRegister - See PCRegisterDynamic()
1405 
1406   Level: advanced
1407 @*/
1408 PetscErrorCode PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1409 {
1410   PetscErrorCode ierr;
1411   char           fullname[PETSC_MAX_PATH_LEN];
1412 
1413   PetscFunctionBegin;
1414 
1415   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1416   ierr = PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 #undef __FUNCT__
1421 #define __FUNCT__ "PCComputeExplicitOperator"
1422 /*@
1423     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1424 
1425     Collective on PC
1426 
1427     Input Parameter:
1428 .   pc - the preconditioner object
1429 
1430     Output Parameter:
1431 .   mat - the explict preconditioned operator
1432 
1433     Notes:
1434     This computation is done by applying the operators to columns of the
1435     identity matrix.
1436 
1437     Currently, this routine uses a dense matrix format when 1 processor
1438     is used and a sparse format otherwise.  This routine is costly in general,
1439     and is recommended for use only with relatively small systems.
1440 
1441     Level: advanced
1442 
1443 .keywords: PC, compute, explicit, operator
1444 
1445 @*/
1446 PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat)
1447 {
1448   Vec            in,out;
1449   PetscErrorCode ierr;
1450   PetscInt       i,M,m,*rows,start,end;
1451   PetscMPIInt    size;
1452   MPI_Comm       comm;
1453   PetscScalar    *array,zero = 0.0,one = 1.0;
1454 
1455   PetscFunctionBegin;
1456   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1457   PetscValidPointer(mat,2);
1458 
1459   comm = pc->comm;
1460   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1461 
1462   if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1463   ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1464   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1465   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1466   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1467   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1468   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr);
1469   for (i=0; i<m; i++) {rows[i] = start + i;}
1470 
1471   ierr = MatCreate(comm,m,m,M,M,mat);CHKERRQ(ierr);
1472   if (size == 1) {
1473     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1474     ierr = MatSeqDenseSetPreallocation(*mat,PETSC_NULL);CHKERRQ(ierr);
1475   } else {
1476     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1477     ierr = MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1478   }
1479 
1480   for (i=0; i<M; i++) {
1481 
1482     ierr = VecSet(&zero,in);CHKERRQ(ierr);
1483     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1484     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1485     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1486 
1487     /* should fix, allowing user to choose side */
1488     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1489 
1490     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1491     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1492     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1493 
1494   }
1495   ierr = PetscFree(rows);CHKERRQ(ierr);
1496   ierr = VecDestroy(out);CHKERRQ(ierr);
1497   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1498   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1499   PetscFunctionReturn(0);
1500 }
1501 
1502