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