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