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