xref: /petsc/src/ksp/pc/interface/precon.c (revision 6ec9520ff15c6c6c97b00a55e49c5b80be990a1d)
1 
2 /*
3     The PC (preconditioner) interface routines, callable by users.
4 */
5 #include <petsc-private/pcimpl.h>            /*I "petscksp.h" I*/
6 #include <petscdm.h>
7 
8 /* Logging support */
9 PetscClassId  PC_CLASSID;
10 PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
11 PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks, PC_ApplyOnMproc;
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   PetscBool      flg1,flg2,set,flg3;
20 
21   PetscFunctionBegin;
22   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
23   if (pc->pmat) {
24     PetscErrorCode (*f)(Mat,MatReuse,Mat*);
25     ierr = PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr);
26     if (size == 1) {
27       ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);CHKERRQ(ierr);
28       ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&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__ "PCReset"
58 /*@
59    PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated Vecs and Mats
60 
61    Collective on PC
62 
63    Input Parameter:
64 .  pc - the preconditioner context
65 
66    Level: developer
67 
68    Notes: This allows a PC to be reused for a different sized linear system but using the same options that have been previously set in the PC
69 
70 .keywords: PC, destroy
71 
72 .seealso: PCCreate(), PCSetUp()
73 @*/
74 PetscErrorCode  PCReset(PC pc)
75 {
76   PetscErrorCode ierr;
77 
78   PetscFunctionBegin;
79   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
80   if (pc->ops->reset) {
81     ierr = (*pc->ops->reset)(pc);CHKERRQ(ierr);
82   }
83   ierr = VecDestroy(&pc->diagonalscaleright);CHKERRQ(ierr);
84   ierr = VecDestroy(&pc->diagonalscaleleft);CHKERRQ(ierr);
85   ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr);
86   ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
87 
88   pc->setupcalled = 0;
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "PCDestroy"
94 /*@
95    PCDestroy - Destroys PC context that was created with PCCreate().
96 
97    Collective on PC
98 
99    Input Parameter:
100 .  pc - the preconditioner context
101 
102    Level: developer
103 
104 .keywords: PC, destroy
105 
106 .seealso: PCCreate(), PCSetUp()
107 @*/
108 PetscErrorCode  PCDestroy(PC *pc)
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (!*pc) PetscFunctionReturn(0);
114   PetscValidHeaderSpecific((*pc),PC_CLASSID,1);
115   if (--((PetscObject)(*pc))->refct > 0) {*pc = 0; PetscFunctionReturn(0);}
116 
117   ierr = PCReset(*pc);CHKERRQ(ierr);
118 
119   /* if memory was published with SAWs then destroy it */
120   ierr = PetscObjectSAWsViewOff((PetscObject)*pc);CHKERRQ(ierr);
121   if ((*pc)->ops->destroy) {ierr = (*(*pc)->ops->destroy)((*pc));CHKERRQ(ierr);}
122   ierr = DMDestroy(&(*pc)->dm);CHKERRQ(ierr);
123   ierr = PetscHeaderDestroy(pc);CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "PCGetDiagonalScale"
129 /*@C
130    PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
131       scaling as needed by certain time-stepping codes.
132 
133    Logically Collective on PC
134 
135    Input Parameter:
136 .  pc - the preconditioner context
137 
138    Output Parameter:
139 .  flag - PETSC_TRUE if it applies the scaling
140 
141    Level: developer
142 
143    Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
144 $           D M A D^{-1} y = D M b  for left preconditioning or
145 $           D A M D^{-1} z = D b for right preconditioning
146 
147 .keywords: PC
148 
149 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale()
150 @*/
151 PetscErrorCode  PCGetDiagonalScale(PC pc,PetscBool  *flag)
152 {
153   PetscFunctionBegin;
154   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
155   PetscValidPointer(flag,2);
156   *flag = pc->diagonalscale;
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "PCSetDiagonalScale"
162 /*@
163    PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
164       scaling as needed by certain time-stepping codes.
165 
166    Logically Collective on PC
167 
168    Input Parameters:
169 +  pc - the preconditioner context
170 -  s - scaling vector
171 
172    Level: intermediate
173 
174    Notes: The system solved via the Krylov method is
175 $           D M A D^{-1} y = D M b  for left preconditioning or
176 $           D A M D^{-1} z = D b for right preconditioning
177 
178    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
179 
180 .keywords: PC
181 
182 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale()
183 @*/
184 PetscErrorCode  PCSetDiagonalScale(PC pc,Vec s)
185 {
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
190   PetscValidHeaderSpecific(s,VEC_CLASSID,2);
191   pc->diagonalscale     = PETSC_TRUE;
192 
193   ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
194   ierr = VecDestroy(&pc->diagonalscaleleft);CHKERRQ(ierr);
195 
196   pc->diagonalscaleleft = s;
197 
198   ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr);
199   ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr);
200   ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr);
201   PetscFunctionReturn(0);
202 }
203 
204 #undef __FUNCT__
205 #define __FUNCT__ "PCDiagonalScaleLeft"
206 /*@
207    PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
208 
209    Logically Collective on PC
210 
211    Input Parameters:
212 +  pc - the preconditioner context
213 .  in - input vector
214 +  out - scaled vector (maybe the same as in)
215 
216    Level: intermediate
217 
218    Notes: The system solved via the Krylov method is
219 $           D M A D^{-1} y = D M b  for left preconditioning or
220 $           D A M D^{-1} z = D b for right preconditioning
221 
222    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
223 
224    If diagonal scaling is turned off and in is not out then in is copied to out
225 
226 .keywords: PC
227 
228 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
229 @*/
230 PetscErrorCode  PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
231 {
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
236   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
237   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
238   if (pc->diagonalscale) {
239     ierr = VecPointwiseMult(out,pc->diagonalscaleleft,in);CHKERRQ(ierr);
240   } else if (in != out) {
241     ierr = VecCopy(in,out);CHKERRQ(ierr);
242   }
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "PCDiagonalScaleRight"
248 /*@
249    PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
250 
251    Logically Collective on PC
252 
253    Input Parameters:
254 +  pc - the preconditioner context
255 .  in - input vector
256 +  out - scaled vector (maybe the same as in)
257 
258    Level: intermediate
259 
260    Notes: The system solved via the Krylov method is
261 $           D M A D^{-1} y = D M b  for left preconditioning or
262 $           D A M D^{-1} z = D b for right preconditioning
263 
264    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
265 
266    If diagonal scaling is turned off and in is not out then in is copied to out
267 
268 .keywords: PC
269 
270 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
271 @*/
272 PetscErrorCode  PCDiagonalScaleRight(PC pc,Vec in,Vec out)
273 {
274   PetscErrorCode ierr;
275 
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
278   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
279   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
280   if (pc->diagonalscale) {
281     ierr = VecPointwiseMult(out,pc->diagonalscaleright,in);CHKERRQ(ierr);
282   } else if (in != out) {
283     ierr = VecCopy(in,out);CHKERRQ(ierr);
284   }
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "PCSetUseAmat"
290 /*@
291    PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
292    operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
293    TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
294 
295    Logically Collective on PC
296 
297    Input Parameters:
298 +  pc - the preconditioner context
299 -  flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
300 
301    Options Database Key:
302 .  -pc_use_amat <true,false>
303 
304    Notes:
305    For the common case in which the linear system matrix and the matrix used to construct the
306    preconditioner are identical, this routine is does nothing.
307 
308    Level: intermediate
309 
310 .seealso: PCGetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
311 @*/
312 PetscErrorCode  PCSetUseAmat(PC pc,PetscBool flg)
313 {
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
316   pc->useAmat = flg;
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "PCGetUseAmat"
322 /*@
323    PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
324    operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
325    TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
326 
327    Logically Collective on PC
328 
329    Input Parameter:
330 .  pc - the preconditioner context
331 
332    Output Parameter:
333 .  flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
334 
335    Notes:
336    For the common case in which the linear system matrix and the matrix used to construct the
337    preconditioner are identical, this routine is does nothing.
338 
339    Level: intermediate
340 
341 .seealso: PCSetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
342 @*/
343 PetscErrorCode  PCGetUseAmat(PC pc,PetscBool *flg)
344 {
345   PetscFunctionBegin;
346   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
347   *flg = pc->useAmat;
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "PCCreate"
353 /*@
354    PCCreate - Creates a preconditioner context.
355 
356    Collective on MPI_Comm
357 
358    Input Parameter:
359 .  comm - MPI communicator
360 
361    Output Parameter:
362 .  pc - location to put the preconditioner context
363 
364    Notes:
365    The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
366    in parallel. For dense matrices it is always PCNONE.
367 
368    Level: developer
369 
370 .keywords: PC, create, context
371 
372 .seealso: PCSetUp(), PCApply(), PCDestroy()
373 @*/
374 PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
375 {
376   PC             pc;
377   PetscErrorCode ierr;
378 
379   PetscFunctionBegin;
380   PetscValidPointer(newpc,1);
381   *newpc = 0;
382   ierr = PCInitializePackage();CHKERRQ(ierr);
383 
384   ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,"PC","Preconditioner","PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
385 
386   pc->mat                  = 0;
387   pc->pmat                 = 0;
388   pc->setupcalled          = 0;
389   pc->setfromoptionscalled = 0;
390   pc->data                 = 0;
391   pc->diagonalscale        = PETSC_FALSE;
392   pc->diagonalscaleleft    = 0;
393   pc->diagonalscaleright   = 0;
394 
395   pc->modifysubmatrices  = 0;
396   pc->modifysubmatricesP = 0;
397 
398   *newpc = pc;
399   PetscFunctionReturn(0);
400 
401 }
402 
403 /* -------------------------------------------------------------------------------*/
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "PCApply"
407 /*@
408    PCApply - Applies the preconditioner to a vector.
409 
410    Collective on PC and Vec
411 
412    Input Parameters:
413 +  pc - the preconditioner context
414 -  x - input vector
415 
416    Output Parameter:
417 .  y - output vector
418 
419    Level: developer
420 
421 .keywords: PC, apply
422 
423 .seealso: PCApplyTranspose(), PCApplyBAorAB()
424 @*/
425 PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
426 {
427   PetscErrorCode ierr;
428   PetscInt       m,n,mv,nv;
429 
430   PetscFunctionBegin;
431   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
432   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
433   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
434   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
435   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
436   ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr);
437   ierr = VecGetLocalSize(x,&nv);CHKERRQ(ierr);
438   ierr = VecGetLocalSize(y,&mv);CHKERRQ(ierr);
439   if (mv != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Preconditioner number of local rows %D does not equal resulting vector number of rows %D",m,mv);CHKERRQ(ierr);
440   if (nv != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Preconditioner number of local columns %D does not equal resulting vector number of rows %D",n,nv);CHKERRQ(ierr);
441   VecLocked(y,3);
442 
443   if (pc->setupcalled < 2) {
444     ierr = PCSetUp(pc);CHKERRQ(ierr);
445   }
446   if (!pc->ops->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply");
447   ierr = VecLockPush(x);CHKERRQ(ierr);
448   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
449   ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr);
450   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
451   ierr = VecLockPop(x);CHKERRQ(ierr);
452   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 
456 #undef __FUNCT__
457 #define __FUNCT__ "PCApplySymmetricLeft"
458 /*@
459    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
460 
461    Collective on PC and Vec
462 
463    Input Parameters:
464 +  pc - the preconditioner context
465 -  x - input vector
466 
467    Output Parameter:
468 .  y - output vector
469 
470    Notes:
471    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
472 
473    Level: developer
474 
475 .keywords: PC, apply, symmetric, left
476 
477 .seealso: PCApply(), PCApplySymmetricRight()
478 @*/
479 PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
480 {
481   PetscErrorCode ierr;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
485   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
486   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
487   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
488   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
489   if (pc->setupcalled < 2) {
490     ierr = PCSetUp(pc);CHKERRQ(ierr);
491   }
492   if (!pc->ops->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
493   ierr = VecLockPush(x);CHKERRQ(ierr);
494   ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
495   ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr);
496   ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
497   ierr = VecLockPop(x);CHKERRQ(ierr);
498   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
499   PetscFunctionReturn(0);
500 }
501 
502 #undef __FUNCT__
503 #define __FUNCT__ "PCApplySymmetricRight"
504 /*@
505    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
506 
507    Collective on PC and Vec
508 
509    Input Parameters:
510 +  pc - the preconditioner context
511 -  x - input vector
512 
513    Output Parameter:
514 .  y - output vector
515 
516    Level: developer
517 
518    Notes:
519    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
520 
521 .keywords: PC, apply, symmetric, right
522 
523 .seealso: PCApply(), PCApplySymmetricLeft()
524 @*/
525 PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
526 {
527   PetscErrorCode ierr;
528 
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
531   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
532   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
533   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
534   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
535   if (pc->setupcalled < 2) {
536     ierr = PCSetUp(pc);CHKERRQ(ierr);
537   }
538   if (!pc->ops->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
539   ierr = VecLockPush(x);CHKERRQ(ierr);
540   ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
541   ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr);
542   ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
543   ierr = VecLockPop(x);CHKERRQ(ierr);
544   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
545   PetscFunctionReturn(0);
546 }
547 
548 #undef __FUNCT__
549 #define __FUNCT__ "PCApplyTranspose"
550 /*@
551    PCApplyTranspose - Applies the transpose of preconditioner to a vector.
552 
553    Collective on PC and Vec
554 
555    Input Parameters:
556 +  pc - the preconditioner context
557 -  x - input vector
558 
559    Output Parameter:
560 .  y - output vector
561 
562    Notes: For complex numbers this applies the non-Hermitian transpose.
563 
564    Developer Notes: We need to implement a PCApplyHermitianTranspose()
565 
566    Level: developer
567 
568 .keywords: PC, apply, transpose
569 
570 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
571 @*/
572 PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
573 {
574   PetscErrorCode ierr;
575 
576   PetscFunctionBegin;
577   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
578   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
579   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
580   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
581   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
582   if (pc->setupcalled < 2) {
583     ierr = PCSetUp(pc);CHKERRQ(ierr);
584   }
585   if (!pc->ops->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply transpose");
586   ierr = VecLockPush(x);CHKERRQ(ierr);
587   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
588   ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr);
589   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
590   ierr = VecLockPop(x);CHKERRQ(ierr);
591   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
592   PetscFunctionReturn(0);
593 }
594 
595 #undef __FUNCT__
596 #define __FUNCT__ "PCApplyTransposeExists"
597 /*@
598    PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
599 
600    Collective on PC and Vec
601 
602    Input Parameters:
603 .  pc - the preconditioner context
604 
605    Output Parameter:
606 .  flg - PETSC_TRUE if a transpose operation is defined
607 
608    Level: developer
609 
610 .keywords: PC, apply, transpose
611 
612 .seealso: PCApplyTranspose()
613 @*/
614 PetscErrorCode  PCApplyTransposeExists(PC pc,PetscBool  *flg)
615 {
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
618   PetscValidPointer(flg,2);
619   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
620   else *flg = PETSC_FALSE;
621   PetscFunctionReturn(0);
622 }
623 
624 #undef __FUNCT__
625 #define __FUNCT__ "PCApplyBAorAB"
626 /*@
627    PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
628 
629    Collective on PC and Vec
630 
631    Input Parameters:
632 +  pc - the preconditioner context
633 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
634 .  x - input vector
635 -  work - work vector
636 
637    Output Parameter:
638 .  y - output vector
639 
640    Level: developer
641 
642    Notes: If the PC has had PCSetDiagonalScale() set then D M A D^{-1} for left preconditioning or  D A M D^{-1} is actually applied. Note that the
643    specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling.
644 
645 .keywords: PC, apply, operator
646 
647 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
648 @*/
649 PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
650 {
651   PetscErrorCode ierr;
652 
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
655   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
656   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
657   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
658   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
659   ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);
660   if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
661   if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
662 
663   if (pc->setupcalled < 2) {
664     ierr = PCSetUp(pc);CHKERRQ(ierr);
665   }
666 
667   if (pc->diagonalscale) {
668     if (pc->ops->applyBA) {
669       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
670       ierr = VecDuplicate(x,&work2);CHKERRQ(ierr);
671       ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr);
672       ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr);
673       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
674       ierr = VecDestroy(&work2);CHKERRQ(ierr);
675     } else if (side == PC_RIGHT) {
676       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
677       ierr = PCApply(pc,y,work);CHKERRQ(ierr);
678       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
679       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
680     } else if (side == PC_LEFT) {
681       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
682       ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr);
683       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
684       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
685     } else if (side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
686   } else {
687     if (pc->ops->applyBA) {
688       ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr);
689     } else if (side == PC_RIGHT) {
690       ierr = PCApply(pc,x,work);CHKERRQ(ierr);
691       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
692     } else if (side == PC_LEFT) {
693       ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr);
694       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
695     } else if (side == PC_SYMMETRIC) {
696       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
697       ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr);
698       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
699       ierr = VecCopy(y,work);CHKERRQ(ierr);
700       ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr);
701     }
702   }
703   ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 #undef __FUNCT__
708 #define __FUNCT__ "PCApplyBAorABTranspose"
709 /*@
710    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
711    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
712    NOT tr(B*A) = tr(A)*tr(B).
713 
714    Collective on PC and Vec
715 
716    Input Parameters:
717 +  pc - the preconditioner context
718 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
719 .  x - input vector
720 -  work - work vector
721 
722    Output Parameter:
723 .  y - output vector
724 
725 
726    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
727       defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
728 
729     Level: developer
730 
731 .keywords: PC, apply, operator, transpose
732 
733 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
734 @*/
735 PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
736 {
737   PetscErrorCode ierr;
738 
739   PetscFunctionBegin;
740   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
741   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
742   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
743   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
744   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
745   ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);
746   if (pc->ops->applyBAtranspose) {
747     ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr);
748     ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
749     PetscFunctionReturn(0);
750   }
751   if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
752 
753   if (pc->setupcalled < 2) {
754     ierr = PCSetUp(pc);CHKERRQ(ierr);
755   }
756 
757   if (side == PC_RIGHT) {
758     ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr);
759     ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr);
760   } else if (side == PC_LEFT) {
761     ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr);
762     ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr);
763   }
764   /* add support for PC_SYMMETRIC */
765   ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
766   PetscFunctionReturn(0);
767 }
768 
769 /* -------------------------------------------------------------------------------*/
770 
771 #undef __FUNCT__
772 #define __FUNCT__ "PCApplyRichardsonExists"
773 /*@
774    PCApplyRichardsonExists - Determines whether a particular preconditioner has a
775    built-in fast application of Richardson's method.
776 
777    Not Collective
778 
779    Input Parameter:
780 .  pc - the preconditioner
781 
782    Output Parameter:
783 .  exists - PETSC_TRUE or PETSC_FALSE
784 
785    Level: developer
786 
787 .keywords: PC, apply, Richardson, exists
788 
789 .seealso: PCApplyRichardson()
790 @*/
791 PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscBool  *exists)
792 {
793   PetscFunctionBegin;
794   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
795   PetscValidIntPointer(exists,2);
796   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
797   else *exists = PETSC_FALSE;
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "PCApplyRichardson"
803 /*@
804    PCApplyRichardson - Applies several steps of Richardson iteration with
805    the particular preconditioner. This routine is usually used by the
806    Krylov solvers and not the application code directly.
807 
808    Collective on PC
809 
810    Input Parameters:
811 +  pc  - the preconditioner context
812 .  b   - the right hand side
813 .  w   - one work vector
814 .  rtol - relative decrease in residual norm convergence criteria
815 .  abstol - absolute residual norm convergence criteria
816 .  dtol - divergence residual norm increase criteria
817 .  its - the number of iterations to apply.
818 -  guesszero - if the input x contains nonzero initial guess
819 
820    Output Parameter:
821 +  outits - number of iterations actually used (for SOR this always equals its)
822 .  reason - the reason the apply terminated
823 -  y - the solution (also contains initial guess if guesszero is PETSC_FALSE
824 
825    Notes:
826    Most preconditioners do not support this function. Use the command
827    PCApplyRichardsonExists() to determine if one does.
828 
829    Except for the multigrid PC this routine ignores the convergence tolerances
830    and always runs for the number of iterations
831 
832    Level: developer
833 
834 .keywords: PC, apply, Richardson
835 
836 .seealso: PCApplyRichardsonExists()
837 @*/
838 PetscErrorCode  PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
839 {
840   PetscErrorCode ierr;
841 
842   PetscFunctionBegin;
843   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
844   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
845   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
846   PetscValidHeaderSpecific(w,VEC_CLASSID,4);
847   if (b == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"b and y must be different vectors");
848   if (pc->setupcalled < 2) {
849     ierr = PCSetUp(pc);CHKERRQ(ierr);
850   }
851   if (!pc->ops->applyrichardson) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply richardson");
852   ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr);
853   PetscFunctionReturn(0);
854 }
855 
856 /*
857       a setupcall of 0 indicates never setup,
858                      1 indicates has been previously setup
859 */
860 #undef __FUNCT__
861 #define __FUNCT__ "PCSetUp"
862 /*@
863    PCSetUp - Prepares for the use of a preconditioner.
864 
865    Collective on PC
866 
867    Input Parameter:
868 .  pc - the preconditioner context
869 
870    Level: developer
871 
872 .keywords: PC, setup
873 
874 .seealso: PCCreate(), PCApply(), PCDestroy()
875 @*/
876 PetscErrorCode  PCSetUp(PC pc)
877 {
878   PetscErrorCode   ierr;
879   const char       *def;
880   PetscObjectState matstate, matnonzerostate;
881 
882   PetscFunctionBegin;
883   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
884   if (!pc->mat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");
885 
886   if (pc->setupcalled && pc->reusepreconditioner) {
887     ierr = PetscInfo(pc,"Leaving PC with identical preconditioner since reuse preconditioner is set\n");CHKERRQ(ierr);
888     PetscFunctionReturn(0);
889   }
890 
891   ierr = PetscObjectStateGet((PetscObject)pc->pmat,&matstate);CHKERRQ(ierr);
892   ierr = MatGetNonzeroState(pc->pmat,&matnonzerostate);CHKERRQ(ierr);
893   if (!pc->setupcalled) {
894     ierr            = PetscInfo(pc,"Setting up PC for first time\n");CHKERRQ(ierr);
895     pc->flag        = DIFFERENT_NONZERO_PATTERN;
896   } else if (matstate == pc->matstate) {
897     ierr = PetscInfo(pc,"Leaving PC with identical preconditioner since operator is unchanged\n");CHKERRQ(ierr);
898     PetscFunctionReturn(0);
899   } else {
900     if (matnonzerostate > pc->matnonzerostate) {
901        ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
902        pc->flag            = DIFFERENT_NONZERO_PATTERN;
903     } else {
904       ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
905       pc->flag            = SAME_NONZERO_PATTERN;
906     }
907   }
908   pc->matstate        = matstate;
909   pc->matnonzerostate = matnonzerostate;
910 
911   if (!((PetscObject)pc)->type_name) {
912     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
913     ierr = PCSetType(pc,def);CHKERRQ(ierr);
914   }
915 
916   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
917   if (pc->ops->setup) {
918     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
919   }
920   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
921   pc->setupcalled = 1;
922   PetscFunctionReturn(0);
923 }
924 
925 #undef __FUNCT__
926 #define __FUNCT__ "PCSetUpOnBlocks"
927 /*@
928    PCSetUpOnBlocks - Sets up the preconditioner for each block in
929    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
930    methods.
931 
932    Collective on PC
933 
934    Input Parameters:
935 .  pc - the preconditioner context
936 
937    Level: developer
938 
939 .keywords: PC, setup, blocks
940 
941 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
942 @*/
943 PetscErrorCode  PCSetUpOnBlocks(PC pc)
944 {
945   PetscErrorCode ierr;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
949   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
950   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
951   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
952   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
953   PetscFunctionReturn(0);
954 }
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "PCSetModifySubMatrices"
958 /*@C
959    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
960    submatrices that arise within certain subdomain-based preconditioners.
961    The basic submatrices are extracted from the preconditioner matrix as
962    usual; the user can then alter these (for example, to set different boundary
963    conditions for each submatrix) before they are used for the local solves.
964 
965    Logically Collective on PC
966 
967    Input Parameters:
968 +  pc - the preconditioner context
969 .  func - routine for modifying the submatrices
970 -  ctx - optional user-defined context (may be null)
971 
972    Calling sequence of func:
973 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
974 
975 .  row - an array of index sets that contain the global row numbers
976          that comprise each local submatrix
977 .  col - an array of index sets that contain the global column numbers
978          that comprise each local submatrix
979 .  submat - array of local submatrices
980 -  ctx - optional user-defined context for private data for the
981          user-defined func routine (may be null)
982 
983    Notes:
984    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
985    KSPSolve().
986 
987    A routine set by PCSetModifySubMatrices() is currently called within
988    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
989    preconditioners.  All other preconditioners ignore this routine.
990 
991    Level: advanced
992 
993 .keywords: PC, set, modify, submatrices
994 
995 .seealso: PCModifySubMatrices(), PCASMGetSubMatrices()
996 @*/
997 PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
998 {
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1001   pc->modifysubmatrices  = func;
1002   pc->modifysubmatricesP = ctx;
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 #undef __FUNCT__
1007 #define __FUNCT__ "PCModifySubMatrices"
1008 /*@C
1009    PCModifySubMatrices - Calls an optional user-defined routine within
1010    certain preconditioners if one has been set with PCSetModifySubMarices().
1011 
1012    Collective on PC
1013 
1014    Input Parameters:
1015 +  pc - the preconditioner context
1016 .  nsub - the number of local submatrices
1017 .  row - an array of index sets that contain the global row numbers
1018          that comprise each local submatrix
1019 .  col - an array of index sets that contain the global column numbers
1020          that comprise each local submatrix
1021 .  submat - array of local submatrices
1022 -  ctx - optional user-defined context for private data for the
1023          user-defined routine (may be null)
1024 
1025    Output Parameter:
1026 .  submat - array of local submatrices (the entries of which may
1027             have been modified)
1028 
1029    Notes:
1030    The user should NOT generally call this routine, as it will
1031    automatically be called within certain preconditioners (currently
1032    block Jacobi, additive Schwarz) if set.
1033 
1034    The basic submatrices are extracted from the preconditioner matrix
1035    as usual; the user can then alter these (for example, to set different
1036    boundary conditions for each submatrix) before they are used for the
1037    local solves.
1038 
1039    Level: developer
1040 
1041 .keywords: PC, modify, submatrices
1042 
1043 .seealso: PCSetModifySubMatrices()
1044 @*/
1045 PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
1046 {
1047   PetscErrorCode ierr;
1048 
1049   PetscFunctionBegin;
1050   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1051   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
1052   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
1053   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
1054   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "PCSetOperators"
1060 /*@
1061    PCSetOperators - Sets the matrix associated with the linear system and
1062    a (possibly) different one associated with the preconditioner.
1063 
1064    Logically Collective on PC and Mat
1065 
1066    Input Parameters:
1067 +  pc - the preconditioner context
1068 .  Amat - the matrix that defines the linear system
1069 -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
1070 
1071    Notes:
1072     Passing a NULL for Amat or Pmat removes the matrix that is currently used.
1073 
1074     If you wish to replace either Amat or Pmat but leave the other one untouched then
1075     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
1076     on it and then pass it back in in your call to KSPSetOperators().
1077 
1078    More Notes about Repeated Solution of Linear Systems:
1079    PETSc does NOT reset the matrix entries of either Amat or Pmat
1080    to zero after a linear solve; the user is completely responsible for
1081    matrix assembly.  See the routine MatZeroEntries() if desiring to
1082    zero all elements of a matrix.
1083 
1084    Level: intermediate
1085 
1086 .keywords: PC, set, operators, matrix, linear system
1087 
1088 .seealso: PCGetOperators(), MatZeroEntries()
1089  @*/
1090 PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat)
1091 {
1092   PetscErrorCode   ierr;
1093   PetscInt         m1,n1,m2,n2;
1094 
1095   PetscFunctionBegin;
1096   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1097   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1098   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1099   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1100   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1101   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
1102     ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr);
1103     ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr);
1104     if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Amat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1);
1105     ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr);
1106     ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr);
1107     if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Pmat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1);
1108   }
1109 
1110   if (Pmat != pc->pmat) {
1111     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1112     pc->matnonzerostate = -1;
1113     pc->matstate        = -1;
1114   }
1115 
1116   /* reference first in case the matrices are the same */
1117   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1118   ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1119   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1120   ierr     = MatDestroy(&pc->pmat);CHKERRQ(ierr);
1121   pc->mat  = Amat;
1122   pc->pmat = Pmat;
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "PCSetReusePreconditioner"
1128 /*@
1129    PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
1130 
1131    Logically Collective on PC
1132 
1133    Input Parameters:
1134 +  pc - the preconditioner context
1135 -  flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1136 
1137     Level: intermediate
1138 
1139 .seealso: PCGetOperators(), MatZeroEntries(), PCGetReusePreconditioner(), KSPSetReusePreconditioner()
1140  @*/
1141 PetscErrorCode  PCSetReusePreconditioner(PC pc,PetscBool flag)
1142 {
1143   PetscFunctionBegin;
1144   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1145   pc->reusepreconditioner = flag;
1146   PetscFunctionReturn(0);
1147 }
1148 
1149 #undef __FUNCT__
1150 #define __FUNCT__ "PCGetReusePreconditioner"
1151 /*@
1152    PCGetReusePreconditioner - Determines if the PC reuses the current preconditioner even if the operator in the preconditioner has changed.
1153 
1154    Not Collective
1155 
1156    Input Parameter:
1157 .  pc - the preconditioner context
1158 
1159    Output Parameter:
1160 .  flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1161 
1162 .seealso: PCGetOperators(), MatZeroEntries(), PCSetReusePreconditioner()
1163  @*/
1164 PetscErrorCode  PCGetReusePreconditioner(PC pc,PetscBool *flag)
1165 {
1166   PetscFunctionBegin;
1167   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1168   *flag = pc->reusepreconditioner;
1169   PetscFunctionReturn(0);
1170 }
1171 
1172 #undef __FUNCT__
1173 #define __FUNCT__ "PCGetOperators"
1174 /*@C
1175    PCGetOperators - Gets the matrix associated with the linear system and
1176    possibly a different one associated with the preconditioner.
1177 
1178    Not collective, though parallel Mats are returned if the PC is parallel
1179 
1180    Input Parameter:
1181 .  pc - the preconditioner context
1182 
1183    Output Parameters:
1184 +  Amat - the matrix defining the linear system
1185 -  Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1186 
1187    Level: intermediate
1188 
1189    Notes: Does not increase the reference count of the matrices, so you should not destroy them
1190 
1191    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1192       are created in PC and returned to the user. In this case, if both operators
1193       mat and pmat are requested, two DIFFERENT operators will be returned. If
1194       only one is requested both operators in the PC will be the same (i.e. as
1195       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1196       The user must set the sizes of the returned matrices and their type etc just
1197       as if the user created them with MatCreate(). For example,
1198 
1199 $         KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1200 $           set size, type, etc of Amat
1201 
1202 $         MatCreate(comm,&mat);
1203 $         KSP/PCSetOperators(ksp/pc,Amat,Amat);
1204 $         PetscObjectDereference((PetscObject)mat);
1205 $           set size, type, etc of Amat
1206 
1207      and
1208 
1209 $         KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1210 $           set size, type, etc of Amat and Pmat
1211 
1212 $         MatCreate(comm,&Amat);
1213 $         MatCreate(comm,&Pmat);
1214 $         KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1215 $         PetscObjectDereference((PetscObject)Amat);
1216 $         PetscObjectDereference((PetscObject)Pmat);
1217 $           set size, type, etc of Amat and Pmat
1218 
1219     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1220     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1221     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1222     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1223     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1224     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1225     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1226     it can be created for you?
1227 
1228 
1229 .keywords: PC, get, operators, matrix, linear system
1230 
1231 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1232 @*/
1233 PetscErrorCode  PCGetOperators(PC pc,Mat *Amat,Mat *Pmat)
1234 {
1235   PetscErrorCode ierr;
1236 
1237   PetscFunctionBegin;
1238   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1239   if (Amat) {
1240     if (!pc->mat) {
1241       if (pc->pmat && !Pmat) {  /* Apmat has been set, but user did not request it, so use for Amat */
1242         pc->mat = pc->pmat;
1243         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1244       } else {                  /* both Amat and Pmat are empty */
1245         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr);
1246         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1247           pc->pmat = pc->mat;
1248           ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1249         }
1250       }
1251     }
1252     *Amat = pc->mat;
1253   }
1254   if (Pmat) {
1255     if (!pc->pmat) {
1256       if (pc->mat && !Amat) {    /* Amat has been set but was not requested, so use for pmat */
1257         pc->pmat = pc->mat;
1258         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1259       } else {
1260         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr);
1261         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1262           pc->mat = pc->pmat;
1263           ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1264         }
1265       }
1266     }
1267     *Pmat = pc->pmat;
1268   }
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNCT__
1273 #define __FUNCT__ "PCGetOperatorsSet"
1274 /*@C
1275    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1276    possibly a different one associated with the preconditioner have been set in the PC.
1277 
1278    Not collective, though the results on all processes should be the same
1279 
1280    Input Parameter:
1281 .  pc - the preconditioner context
1282 
1283    Output Parameters:
1284 +  mat - the matrix associated with the linear system was set
1285 -  pmat - matrix associated with the preconditioner was set, usually the same
1286 
1287    Level: intermediate
1288 
1289 .keywords: PC, get, operators, matrix, linear system
1290 
1291 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1292 @*/
1293 PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1294 {
1295   PetscFunctionBegin;
1296   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1297   if (mat) *mat = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1298   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1299   PetscFunctionReturn(0);
1300 }
1301 
1302 #undef __FUNCT__
1303 #define __FUNCT__ "PCFactorGetMatrix"
1304 /*@
1305    PCFactorGetMatrix - Gets the factored matrix from the
1306    preconditioner context.  This routine is valid only for the LU,
1307    incomplete LU, Cholesky, and incomplete Cholesky methods.
1308 
1309    Not Collective on PC though Mat is parallel if PC is parallel
1310 
1311    Input Parameters:
1312 .  pc - the preconditioner context
1313 
1314    Output parameters:
1315 .  mat - the factored matrix
1316 
1317    Level: advanced
1318 
1319    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1320 
1321 .keywords: PC, get, factored, matrix
1322 @*/
1323 PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1324 {
1325   PetscErrorCode ierr;
1326 
1327   PetscFunctionBegin;
1328   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1329   PetscValidPointer(mat,2);
1330   if (pc->ops->getfactoredmatrix) {
1331     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1332   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 #undef __FUNCT__
1337 #define __FUNCT__ "PCSetOptionsPrefix"
1338 /*@C
1339    PCSetOptionsPrefix - Sets the prefix used for searching for all
1340    PC options in the database.
1341 
1342    Logically Collective on PC
1343 
1344    Input Parameters:
1345 +  pc - the preconditioner context
1346 -  prefix - the prefix string to prepend to all PC option requests
1347 
1348    Notes:
1349    A hyphen (-) must NOT be given at the beginning of the prefix name.
1350    The first character of all runtime options is AUTOMATICALLY the
1351    hyphen.
1352 
1353    Level: advanced
1354 
1355 .keywords: PC, set, options, prefix, database
1356 
1357 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1358 @*/
1359 PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1360 {
1361   PetscErrorCode ierr;
1362 
1363   PetscFunctionBegin;
1364   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1365   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1366   PetscFunctionReturn(0);
1367 }
1368 
1369 #undef __FUNCT__
1370 #define __FUNCT__ "PCAppendOptionsPrefix"
1371 /*@C
1372    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1373    PC options in the database.
1374 
1375    Logically Collective on PC
1376 
1377    Input Parameters:
1378 +  pc - the preconditioner context
1379 -  prefix - the prefix string to prepend to all PC option requests
1380 
1381    Notes:
1382    A hyphen (-) must NOT be given at the beginning of the prefix name.
1383    The first character of all runtime options is AUTOMATICALLY the
1384    hyphen.
1385 
1386    Level: advanced
1387 
1388 .keywords: PC, append, options, prefix, database
1389 
1390 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1391 @*/
1392 PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1393 {
1394   PetscErrorCode ierr;
1395 
1396   PetscFunctionBegin;
1397   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1398   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1399   PetscFunctionReturn(0);
1400 }
1401 
1402 #undef __FUNCT__
1403 #define __FUNCT__ "PCGetOptionsPrefix"
1404 /*@C
1405    PCGetOptionsPrefix - Gets the prefix used for searching for all
1406    PC options in the database.
1407 
1408    Not Collective
1409 
1410    Input Parameters:
1411 .  pc - the preconditioner context
1412 
1413    Output Parameters:
1414 .  prefix - pointer to the prefix string used, is returned
1415 
1416    Notes: On the fortran side, the user should pass in a string 'prifix' of
1417    sufficient length to hold the prefix.
1418 
1419    Level: advanced
1420 
1421 .keywords: PC, get, options, prefix, database
1422 
1423 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1424 @*/
1425 PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1426 {
1427   PetscErrorCode ierr;
1428 
1429   PetscFunctionBegin;
1430   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1431   PetscValidPointer(prefix,2);
1432   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1433   PetscFunctionReturn(0);
1434 }
1435 
1436 #undef __FUNCT__
1437 #define __FUNCT__ "PCPreSolve"
1438 /*@
1439    PCPreSolve - Optional pre-solve phase, intended for any
1440    preconditioner-specific actions that must be performed before
1441    the iterative solve itself.
1442 
1443    Collective on PC
1444 
1445    Input Parameters:
1446 +  pc - the preconditioner context
1447 -  ksp - the Krylov subspace context
1448 
1449    Level: developer
1450 
1451    Sample of Usage:
1452 .vb
1453     PCPreSolve(pc,ksp);
1454     KSPSolve(ksp,b,x);
1455     PCPostSolve(pc,ksp);
1456 .ve
1457 
1458    Notes:
1459    The pre-solve phase is distinct from the PCSetUp() phase.
1460 
1461    KSPSolve() calls this directly, so is rarely called by the user.
1462 
1463 .keywords: PC, pre-solve
1464 
1465 .seealso: PCPostSolve()
1466 @*/
1467 PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1468 {
1469   PetscErrorCode ierr;
1470   Vec            x,rhs;
1471 
1472   PetscFunctionBegin;
1473   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1474   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1475   pc->presolvedone++;
1476   if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice");
1477   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1478   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1479 
1480   if (pc->ops->presolve) {
1481     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1482   }
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "PCPostSolve"
1488 /*@
1489    PCPostSolve - Optional post-solve phase, intended for any
1490    preconditioner-specific actions that must be performed after
1491    the iterative solve itself.
1492 
1493    Collective on PC
1494 
1495    Input Parameters:
1496 +  pc - the preconditioner context
1497 -  ksp - the Krylov subspace context
1498 
1499    Sample of Usage:
1500 .vb
1501     PCPreSolve(pc,ksp);
1502     KSPSolve(ksp,b,x);
1503     PCPostSolve(pc,ksp);
1504 .ve
1505 
1506    Note:
1507    KSPSolve() calls this routine directly, so it is rarely called by the user.
1508 
1509    Level: developer
1510 
1511 .keywords: PC, post-solve
1512 
1513 .seealso: PCPreSolve(), KSPSolve()
1514 @*/
1515 PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1516 {
1517   PetscErrorCode ierr;
1518   Vec            x,rhs;
1519 
1520   PetscFunctionBegin;
1521   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1522   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1523   pc->presolvedone--;
1524   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1525   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1526   if (pc->ops->postsolve) {
1527     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1528   }
1529   PetscFunctionReturn(0);
1530 }
1531 
1532 #undef __FUNCT__
1533 #define __FUNCT__ "PCLoad"
1534 /*@C
1535   PCLoad - Loads a PC that has been stored in binary  with PCView().
1536 
1537   Collective on PetscViewer
1538 
1539   Input Parameters:
1540 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or
1541            some related function before a call to PCLoad().
1542 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1543 
1544    Level: intermediate
1545 
1546   Notes:
1547    The type is determined by the data in the file, any type set into the PC before this call is ignored.
1548 
1549   Notes for advanced users:
1550   Most users should not need to know the details of the binary storage
1551   format, since PCLoad() and PCView() completely hide these details.
1552   But for anyone who's interested, the standard binary matrix storage
1553   format is
1554 .vb
1555      has not yet been determined
1556 .ve
1557 
1558 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad()
1559 @*/
1560 PetscErrorCode  PCLoad(PC newdm, PetscViewer viewer)
1561 {
1562   PetscErrorCode ierr;
1563   PetscBool      isbinary;
1564   PetscInt       classid;
1565   char           type[256];
1566 
1567   PetscFunctionBegin;
1568   PetscValidHeaderSpecific(newdm,PC_CLASSID,1);
1569   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1570   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1571   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1572 
1573   ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr);
1574   if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file");
1575   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
1576   ierr = PCSetType(newdm, type);CHKERRQ(ierr);
1577   if (newdm->ops->load) {
1578     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1579   }
1580   PetscFunctionReturn(0);
1581 }
1582 
1583 #include <petscdraw.h>
1584 #if defined(PETSC_HAVE_SAWS)
1585 #include <petscviewersaws.h>
1586 #endif
1587 #undef __FUNCT__
1588 #define __FUNCT__ "PCView"
1589 /*@C
1590    PCView - Prints the PC data structure.
1591 
1592    Collective on PC
1593 
1594    Input Parameters:
1595 +  PC - the PC context
1596 -  viewer - optional visualization context
1597 
1598    Note:
1599    The available visualization contexts include
1600 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1601 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1602          output where only the first processor opens
1603          the file.  All other processors send their
1604          data to the first processor to print.
1605 
1606    The user can open an alternative visualization contexts with
1607    PetscViewerASCIIOpen() (output to a specified file).
1608 
1609    Level: developer
1610 
1611 .keywords: PC, view
1612 
1613 .seealso: KSPView(), PetscViewerASCIIOpen()
1614 @*/
1615 PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1616 {
1617   PCType            cstr;
1618   PetscErrorCode    ierr;
1619   PetscBool         iascii,isstring,isbinary,isdraw;
1620   PetscViewerFormat format;
1621 #if defined(PETSC_HAVE_SAWS)
1622   PetscBool         issaws;
1623 #endif
1624 
1625   PetscFunctionBegin;
1626   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1627   if (!viewer) {
1628     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
1629   }
1630   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1631   PetscCheckSameComm(pc,1,viewer,2);
1632 
1633   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1634   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1635   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1636   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1637 #if defined(PETSC_HAVE_SAWS)
1638   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);CHKERRQ(ierr);
1639 #endif
1640 
1641   if (iascii) {
1642     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1643     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);CHKERRQ(ierr);
1644     if (!pc->setupcalled) {
1645       ierr = PetscViewerASCIIPrintf(viewer,"  PC has not been set up so information may be incomplete\n");CHKERRQ(ierr);
1646     }
1647     if (pc->ops->view) {
1648       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1649       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1650       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1651     }
1652     if (pc->mat) {
1653       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1654       if (pc->pmat == pc->mat) {
1655         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1656         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1657         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1658         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1659       } else {
1660         if (pc->pmat) {
1661           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1662         } else {
1663           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1664         }
1665         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1666         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1667         if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1668         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1669       }
1670       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1671     }
1672   } else if (isstring) {
1673     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1674     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1675     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1676   } else if (isbinary) {
1677     PetscInt    classid = PC_FILE_CLASSID;
1678     MPI_Comm    comm;
1679     PetscMPIInt rank;
1680     char        type[256];
1681 
1682     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1683     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1684     if (!rank) {
1685       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1686       ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr);
1687       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
1688     }
1689     if (pc->ops->view) {
1690       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1691     }
1692   } else if (isdraw) {
1693     PetscDraw draw;
1694     char      str[25];
1695     PetscReal x,y,bottom,h;
1696     PetscInt  n;
1697 
1698     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1699     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1700     if (pc->mat) {
1701       ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr);
1702       ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr);
1703     } else {
1704       ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr);
1705     }
1706     ierr   = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1707     bottom = y - h;
1708     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1709     if (pc->ops->view) {
1710       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1711     }
1712     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1713 #if defined(PETSC_HAVE_SAWS)
1714   } else if (issaws) {
1715     PetscMPIInt rank;
1716 
1717     ierr = PetscObjectName((PetscObject)pc);CHKERRQ(ierr);
1718     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1719     if (!((PetscObject)pc)->amsmem && !rank) {
1720       ierr = PetscObjectViewSAWs((PetscObject)pc,viewer);CHKERRQ(ierr);
1721     }
1722     if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);}
1723     if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1724 #endif
1725   }
1726   PetscFunctionReturn(0);
1727 }
1728 
1729 
1730 #undef __FUNCT__
1731 #define __FUNCT__ "PCSetInitialGuessNonzero"
1732 /*@
1733    PCSetInitialGuessNonzero - Tells the iterative solver that the
1734    initial guess is nonzero; otherwise PC assumes the initial guess
1735    is to be zero (and thus zeros it out before solving).
1736 
1737    Logically Collective on PC
1738 
1739    Input Parameters:
1740 +  pc - iterative context obtained from PCCreate()
1741 -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1742 
1743    Level: Developer
1744 
1745    Notes:
1746     This is a weird function. Since PC's are linear operators on the right hand side they
1747     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1748     PCKSP and PCREDUNDANT  and causes the inner KSP object to use the nonzero
1749     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1750 
1751 
1752 .keywords: PC, set, initial guess, nonzero
1753 
1754 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1755 @*/
1756 PetscErrorCode  PCSetInitialGuessNonzero(PC pc,PetscBool flg)
1757 {
1758   PetscFunctionBegin;
1759   PetscValidLogicalCollectiveBool(pc,flg,2);
1760   pc->nonzero_guess = flg;
1761   PetscFunctionReturn(0);
1762 }
1763 
1764 #undef __FUNCT__
1765 #define __FUNCT__ "PCGetInitialGuessNonzero"
1766 /*@
1767    PCGetInitialGuessNonzero - Determines if the iterative solver assumes that the
1768    initial guess is nonzero; otherwise PC assumes the initial guess
1769    is to be zero (and thus zeros it out before solving).
1770 
1771    Logically Collective on PC
1772 
1773    Input Parameter:
1774 .   pc - iterative context obtained from PCCreate()
1775 
1776    Output Parameter:
1777 .  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1778 
1779    Level: Developer
1780 
1781 .keywords: PC, set, initial guess, nonzero
1782 
1783 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll(), PCSetInitialGuessNonzero()
1784 @*/
1785 PetscErrorCode  PCGetInitialGuessNonzero(PC pc,PetscBool *flg)
1786 {
1787   PetscFunctionBegin;
1788   *flg = pc->nonzero_guess;
1789   PetscFunctionReturn(0);
1790 }
1791 
1792 #undef __FUNCT__
1793 #define __FUNCT__ "PCRegister"
1794 /*@C
1795   PCRegister -  Adds a method to the preconditioner package.
1796 
1797    Not collective
1798 
1799    Input Parameters:
1800 +  name_solver - name of a new user-defined solver
1801 -  routine_create - routine to create method context
1802 
1803    Notes:
1804    PCRegister() may be called multiple times to add several user-defined preconditioners.
1805 
1806    Sample usage:
1807 .vb
1808    PCRegister("my_solver", MySolverCreate);
1809 .ve
1810 
1811    Then, your solver can be chosen with the procedural interface via
1812 $     PCSetType(pc,"my_solver")
1813    or at runtime via the option
1814 $     -pc_type my_solver
1815 
1816    Level: advanced
1817 
1818 .keywords: PC, register
1819 
1820 .seealso: PCRegisterAll(), PCRegisterDestroy()
1821 @*/
1822 PetscErrorCode  PCRegister(const char sname[],PetscErrorCode (*function)(PC))
1823 {
1824   PetscErrorCode ierr;
1825 
1826   PetscFunctionBegin;
1827   ierr = PetscFunctionListAdd(&PCList,sname,function);CHKERRQ(ierr);
1828   PetscFunctionReturn(0);
1829 }
1830 
1831 #undef __FUNCT__
1832 #define __FUNCT__ "PCComputeExplicitOperator"
1833 /*@
1834     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1835 
1836     Collective on PC
1837 
1838     Input Parameter:
1839 .   pc - the preconditioner object
1840 
1841     Output Parameter:
1842 .   mat - the explict preconditioned operator
1843 
1844     Notes:
1845     This computation is done by applying the operators to columns of the
1846     identity matrix.
1847 
1848     Currently, this routine uses a dense matrix format when 1 processor
1849     is used and a sparse format otherwise.  This routine is costly in general,
1850     and is recommended for use only with relatively small systems.
1851 
1852     Level: advanced
1853 
1854 .keywords: PC, compute, explicit, operator
1855 
1856 .seealso: KSPComputeExplicitOperator()
1857 
1858 @*/
1859 PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1860 {
1861   Vec            in,out;
1862   PetscErrorCode ierr;
1863   PetscInt       i,M,m,*rows,start,end;
1864   PetscMPIInt    size;
1865   MPI_Comm       comm;
1866   PetscScalar    *array,one = 1.0;
1867 
1868   PetscFunctionBegin;
1869   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1870   PetscValidPointer(mat,2);
1871 
1872   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1873   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1874 
1875   if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1876   ierr = MatCreateVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1877   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1878   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1879   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1880   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1881   ierr = PetscMalloc1(m+1,&rows);CHKERRQ(ierr);
1882   for (i=0; i<m; i++) rows[i] = start + i;
1883 
1884   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
1885   ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr);
1886   if (size == 1) {
1887     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1888     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
1889   } else {
1890     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1891     ierr = MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);CHKERRQ(ierr);
1892   }
1893   ierr = MatSetOption(*mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
1894 
1895   for (i=0; i<M; i++) {
1896 
1897     ierr = VecSet(in,0.0);CHKERRQ(ierr);
1898     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1899     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1900     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1901 
1902     /* should fix, allowing user to choose side */
1903     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1904 
1905     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1906     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1907     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1908 
1909   }
1910   ierr = PetscFree(rows);CHKERRQ(ierr);
1911   ierr = VecDestroy(&out);CHKERRQ(ierr);
1912   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1913   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1914   PetscFunctionReturn(0);
1915 }
1916 
1917 #undef __FUNCT__
1918 #define __FUNCT__ "PCSetCoordinates"
1919 /*@
1920    PCSetCoordinates - sets the coordinates of all the nodes on the local process
1921 
1922    Collective on PC
1923 
1924    Input Parameters:
1925 +  pc - the solver context
1926 .  dim - the dimension of the coordinates 1, 2, or 3
1927 -  coords - the coordinates
1928 
1929    Level: intermediate
1930 
1931    Notes: coords is an array of the 3D coordinates for the nodes on
1932    the local processor.  So if there are 108 equation on a processor
1933    for a displacement finite element discretization of elasticity (so
1934    that there are 36 = 108/3 nodes) then the array must have 108
1935    double precision values (ie, 3 * 36).  These x y z coordinates
1936    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1937    ... , N-1.z ].
1938 
1939 .seealso: MatSetNearNullSpace
1940 @*/
1941 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1942 {
1943   PetscErrorCode ierr;
1944 
1945   PetscFunctionBegin;
1946   ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr);
1947   PetscFunctionReturn(0);
1948 }
1949