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