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