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