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