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