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