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