xref: /petsc/src/ksp/pc/interface/precon.c (revision 3338378cc101aa65b47a89f540b98ec4eaa14e3f)
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 = NULL; 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 = NULL;
382   ierr = PCInitializePackage();CHKERRQ(ierr);
383 
384   ierr = PetscHeaderCreate(pc,PC_CLASSID,"PC","Preconditioner","PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
385 
386   pc->mat                  = NULL;
387   pc->pmat                 = NULL;
388   pc->setupcalled          = 0;
389   pc->setfromoptionscalled = 0;
390   pc->data                 = NULL;
391   pc->diagonalscale        = PETSC_FALSE;
392   pc->diagonalscaleleft    = NULL;
393   pc->diagonalscaleright   = NULL;
394 
395   pc->modifysubmatrices  = NULL;
396   pc->modifysubmatricesP = NULL;
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   PetscValidLogicalCollectiveBool(pc,flag,2);
1108   pc->reusepreconditioner = flag;
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 /*@
1113    PCGetReusePreconditioner - Determines if the PC reuses the current preconditioner even if the operator in the preconditioner has changed.
1114 
1115    Not Collective
1116 
1117    Input Parameter:
1118 .  pc - the preconditioner context
1119 
1120    Output Parameter:
1121 .  flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1122 
1123    Level: intermediate
1124 
1125 .seealso: PCGetOperators(), MatZeroEntries(), PCSetReusePreconditioner()
1126  @*/
1127 PetscErrorCode  PCGetReusePreconditioner(PC pc,PetscBool *flag)
1128 {
1129   PetscFunctionBegin;
1130   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1131   PetscValidPointer(flag,2);
1132   *flag = pc->reusepreconditioner;
1133   PetscFunctionReturn(0);
1134 }
1135 
1136 /*@
1137    PCGetOperators - Gets the matrix associated with the linear system and
1138    possibly a different one associated with the preconditioner.
1139 
1140    Not collective, though parallel Mats are returned if the PC is parallel
1141 
1142    Input Parameter:
1143 .  pc - the preconditioner context
1144 
1145    Output Parameters:
1146 +  Amat - the matrix defining the linear system
1147 -  Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1148 
1149    Level: intermediate
1150 
1151    Notes:
1152     Does not increase the reference count of the matrices, so you should not destroy them
1153 
1154    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1155       are created in PC and returned to the user. In this case, if both operators
1156       mat and pmat are requested, two DIFFERENT operators will be returned. If
1157       only one is requested both operators in the PC will be the same (i.e. as
1158       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1159       The user must set the sizes of the returned matrices and their type etc just
1160       as if the user created them with MatCreate(). For example,
1161 
1162 $         KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1163 $           set size, type, etc of Amat
1164 
1165 $         MatCreate(comm,&mat);
1166 $         KSP/PCSetOperators(ksp/pc,Amat,Amat);
1167 $         PetscObjectDereference((PetscObject)mat);
1168 $           set size, type, etc of Amat
1169 
1170      and
1171 
1172 $         KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1173 $           set size, type, etc of Amat and Pmat
1174 
1175 $         MatCreate(comm,&Amat);
1176 $         MatCreate(comm,&Pmat);
1177 $         KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1178 $         PetscObjectDereference((PetscObject)Amat);
1179 $         PetscObjectDereference((PetscObject)Pmat);
1180 $           set size, type, etc of Amat and Pmat
1181 
1182     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1183     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1184     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1185     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1186     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1187     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1188     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1189     it can be created for you?
1190 
1191 
1192 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1193 @*/
1194 PetscErrorCode  PCGetOperators(PC pc,Mat *Amat,Mat *Pmat)
1195 {
1196   PetscErrorCode ierr;
1197 
1198   PetscFunctionBegin;
1199   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1200   if (Amat) {
1201     if (!pc->mat) {
1202       if (pc->pmat && !Pmat) {  /* Apmat has been set, but user did not request it, so use for Amat */
1203         pc->mat = pc->pmat;
1204         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1205       } else {                  /* both Amat and Pmat are empty */
1206         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr);
1207         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1208           pc->pmat = pc->mat;
1209           ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1210         }
1211       }
1212     }
1213     *Amat = pc->mat;
1214   }
1215   if (Pmat) {
1216     if (!pc->pmat) {
1217       if (pc->mat && !Amat) {    /* Amat has been set but was not requested, so use for pmat */
1218         pc->pmat = pc->mat;
1219         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1220       } else {
1221         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr);
1222         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1223           pc->mat = pc->pmat;
1224           ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1225         }
1226       }
1227     }
1228     *Pmat = pc->pmat;
1229   }
1230   PetscFunctionReturn(0);
1231 }
1232 
1233 /*@C
1234    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1235    possibly a different one associated with the preconditioner have been set in the PC.
1236 
1237    Not collective, though the results on all processes should be the same
1238 
1239    Input Parameter:
1240 .  pc - the preconditioner context
1241 
1242    Output Parameters:
1243 +  mat - the matrix associated with the linear system was set
1244 -  pmat - matrix associated with the preconditioner was set, usually the same
1245 
1246    Level: intermediate
1247 
1248 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1249 @*/
1250 PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1251 {
1252   PetscFunctionBegin;
1253   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1254   if (mat) *mat = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1255   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1256   PetscFunctionReturn(0);
1257 }
1258 
1259 /*@
1260    PCFactorGetMatrix - Gets the factored matrix from the
1261    preconditioner context.  This routine is valid only for the LU,
1262    incomplete LU, Cholesky, and incomplete Cholesky methods.
1263 
1264    Not Collective on PC though Mat is parallel if PC is parallel
1265 
1266    Input Parameters:
1267 .  pc - the preconditioner context
1268 
1269    Output parameters:
1270 .  mat - the factored matrix
1271 
1272    Level: advanced
1273 
1274    Notes:
1275     Does not increase the reference count for the matrix so DO NOT destroy it
1276 
1277 @*/
1278 PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1279 {
1280   PetscErrorCode ierr;
1281 
1282   PetscFunctionBegin;
1283   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1284   PetscValidPointer(mat,2);
1285   if (pc->ops->getfactoredmatrix) {
1286     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1287   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1288   PetscFunctionReturn(0);
1289 }
1290 
1291 /*@C
1292    PCSetOptionsPrefix - Sets the prefix used for searching for all
1293    PC options in the database.
1294 
1295    Logically Collective on PC
1296 
1297    Input Parameters:
1298 +  pc - the preconditioner context
1299 -  prefix - the prefix string to prepend to all PC option requests
1300 
1301    Notes:
1302    A hyphen (-) must NOT be given at the beginning of the prefix name.
1303    The first character of all runtime options is AUTOMATICALLY the
1304    hyphen.
1305 
1306    Level: advanced
1307 
1308 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1309 @*/
1310 PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1311 {
1312   PetscErrorCode ierr;
1313 
1314   PetscFunctionBegin;
1315   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1316   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 /*@C
1321    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1322    PC options in the database.
1323 
1324    Logically Collective on PC
1325 
1326    Input Parameters:
1327 +  pc - the preconditioner context
1328 -  prefix - the prefix string to prepend to all PC option requests
1329 
1330    Notes:
1331    A hyphen (-) must NOT be given at the beginning of the prefix name.
1332    The first character of all runtime options is AUTOMATICALLY the
1333    hyphen.
1334 
1335    Level: advanced
1336 
1337 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1338 @*/
1339 PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1340 {
1341   PetscErrorCode ierr;
1342 
1343   PetscFunctionBegin;
1344   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1345   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1346   PetscFunctionReturn(0);
1347 }
1348 
1349 /*@C
1350    PCGetOptionsPrefix - Gets the prefix used for searching for all
1351    PC options in the database.
1352 
1353    Not Collective
1354 
1355    Input Parameters:
1356 .  pc - the preconditioner context
1357 
1358    Output Parameters:
1359 .  prefix - pointer to the prefix string used, is returned
1360 
1361    Notes:
1362     On the fortran side, the user should pass in a string 'prifix' of
1363    sufficient length to hold the prefix.
1364 
1365    Level: advanced
1366 
1367 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1368 @*/
1369 PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1370 {
1371   PetscErrorCode ierr;
1372 
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1375   PetscValidPointer(prefix,2);
1376   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 /*
1381    Indicates the right hand side will be changed by KSPSolve(), this occurs for a few
1382   preconditioners including BDDC and Eisentat that transform the equations before applying
1383   the Krylov methods
1384 */
1385 PETSC_INTERN PetscErrorCode  PCPreSolveChangeRHS(PC pc,PetscBool *change)
1386 {
1387   PetscErrorCode ierr;
1388 
1389   PetscFunctionBegin;
1390   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1391   PetscValidPointer(change,2);
1392   *change = PETSC_FALSE;
1393   ierr = PetscTryMethod(pc,"PCPreSolveChangeRHS_C",(PC,PetscBool*),(pc,change));CHKERRQ(ierr);
1394   PetscFunctionReturn(0);
1395 }
1396 
1397 /*@
1398    PCPreSolve - Optional pre-solve phase, intended for any
1399    preconditioner-specific actions that must be performed before
1400    the iterative solve itself.
1401 
1402    Collective on PC
1403 
1404    Input Parameters:
1405 +  pc - the preconditioner context
1406 -  ksp - the Krylov subspace context
1407 
1408    Level: developer
1409 
1410    Sample of Usage:
1411 .vb
1412     PCPreSolve(pc,ksp);
1413     KSPSolve(ksp,b,x);
1414     PCPostSolve(pc,ksp);
1415 .ve
1416 
1417    Notes:
1418    The pre-solve phase is distinct from the PCSetUp() phase.
1419 
1420    KSPSolve() calls this directly, so is rarely called by the user.
1421 
1422 .seealso: PCPostSolve()
1423 @*/
1424 PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1425 {
1426   PetscErrorCode ierr;
1427   Vec            x,rhs;
1428 
1429   PetscFunctionBegin;
1430   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1431   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1432   pc->presolvedone++;
1433   if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice");
1434   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1435   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1436 
1437   if (pc->ops->presolve) {
1438     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1439   }
1440   PetscFunctionReturn(0);
1441 }
1442 
1443 /*@
1444    PCPostSolve - Optional post-solve phase, intended for any
1445    preconditioner-specific actions that must be performed after
1446    the iterative solve itself.
1447 
1448    Collective on PC
1449 
1450    Input Parameters:
1451 +  pc - the preconditioner context
1452 -  ksp - the Krylov subspace context
1453 
1454    Sample of Usage:
1455 .vb
1456     PCPreSolve(pc,ksp);
1457     KSPSolve(ksp,b,x);
1458     PCPostSolve(pc,ksp);
1459 .ve
1460 
1461    Note:
1462    KSPSolve() calls this routine directly, so it is rarely called by the user.
1463 
1464    Level: developer
1465 
1466 .seealso: PCPreSolve(), KSPSolve()
1467 @*/
1468 PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1469 {
1470   PetscErrorCode ierr;
1471   Vec            x,rhs;
1472 
1473   PetscFunctionBegin;
1474   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1475   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1476   pc->presolvedone--;
1477   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1478   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1479   if (pc->ops->postsolve) {
1480     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1481   }
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 /*@C
1486   PCLoad - Loads a PC that has been stored in binary  with PCView().
1487 
1488   Collective on PetscViewer
1489 
1490   Input Parameters:
1491 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or
1492            some related function before a call to PCLoad().
1493 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1494 
1495    Level: intermediate
1496 
1497   Notes:
1498    The type is determined by the data in the file, any type set into the PC before this call is ignored.
1499 
1500   Notes for advanced users:
1501   Most users should not need to know the details of the binary storage
1502   format, since PCLoad() and PCView() completely hide these details.
1503   But for anyone who's interested, the standard binary matrix storage
1504   format is
1505 .vb
1506      has not yet been determined
1507 .ve
1508 
1509 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad()
1510 @*/
1511 PetscErrorCode  PCLoad(PC newdm, PetscViewer viewer)
1512 {
1513   PetscErrorCode ierr;
1514   PetscBool      isbinary;
1515   PetscInt       classid;
1516   char           type[256];
1517 
1518   PetscFunctionBegin;
1519   PetscValidHeaderSpecific(newdm,PC_CLASSID,1);
1520   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1521   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1522   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1523 
1524   ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr);
1525   if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file");
1526   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
1527   ierr = PCSetType(newdm, type);CHKERRQ(ierr);
1528   if (newdm->ops->load) {
1529     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1530   }
1531   PetscFunctionReturn(0);
1532 }
1533 
1534 #include <petscdraw.h>
1535 #if defined(PETSC_HAVE_SAWS)
1536 #include <petscviewersaws.h>
1537 #endif
1538 
1539 /*@C
1540    PCViewFromOptions - View from Options
1541 
1542    Collective on PC
1543 
1544    Input Parameters:
1545 +  A - the PC context
1546 .  obj - Optional object
1547 -  name - command line option
1548 
1549    Level: intermediate
1550 .seealso:  PC, PCView, PetscObjectViewFromOptions(), PCCreate()
1551 @*/
1552 PetscErrorCode  PCViewFromOptions(PC A,PetscObject obj,const char name[])
1553 {
1554   PetscErrorCode ierr;
1555 
1556   PetscFunctionBegin;
1557   PetscValidHeaderSpecific(A,PC_CLASSID,1);
1558   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
1559   PetscFunctionReturn(0);
1560 }
1561 
1562 /*@C
1563    PCView - Prints the PC data structure.
1564 
1565    Collective on PC
1566 
1567    Input Parameters:
1568 +  PC - the PC context
1569 -  viewer - optional visualization context
1570 
1571    Note:
1572    The available visualization contexts include
1573 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1574 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1575          output where only the first processor opens
1576          the file.  All other processors send their
1577          data to the first processor to print.
1578 
1579    The user can open an alternative visualization contexts with
1580    PetscViewerASCIIOpen() (output to a specified file).
1581 
1582    Level: developer
1583 
1584 .seealso: KSPView(), PetscViewerASCIIOpen()
1585 @*/
1586 PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1587 {
1588   PCType         cstr;
1589   PetscErrorCode ierr;
1590   PetscBool      iascii,isstring,isbinary,isdraw;
1591 #if defined(PETSC_HAVE_SAWS)
1592   PetscBool      issaws;
1593 #endif
1594 
1595   PetscFunctionBegin;
1596   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1597   if (!viewer) {
1598     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
1599   }
1600   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1601   PetscCheckSameComm(pc,1,viewer,2);
1602 
1603   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1604   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1605   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1606   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1607 #if defined(PETSC_HAVE_SAWS)
1608   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);CHKERRQ(ierr);
1609 #endif
1610 
1611   if (iascii) {
1612     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);CHKERRQ(ierr);
1613     if (!pc->setupcalled) {
1614       ierr = PetscViewerASCIIPrintf(viewer,"  PC has not been set up so information may be incomplete\n");CHKERRQ(ierr);
1615     }
1616     if (pc->ops->view) {
1617       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1618       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1619       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1620     }
1621     if (pc->mat) {
1622       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1623       if (pc->pmat == pc->mat) {
1624         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1625         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1626         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1627         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1628       } else {
1629         if (pc->pmat) {
1630           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1631         } else {
1632           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1633         }
1634         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1635         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1636         if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1637         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1638       }
1639       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1640     }
1641   } else if (isstring) {
1642     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1643     ierr = PetscViewerStringSPrintf(viewer," PCType: %-7.7s",cstr);CHKERRQ(ierr);
1644     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1645     if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);}
1646     if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1647   } else if (isbinary) {
1648     PetscInt    classid = PC_FILE_CLASSID;
1649     MPI_Comm    comm;
1650     PetscMPIInt rank;
1651     char        type[256];
1652 
1653     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1654     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1655     if (!rank) {
1656       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
1657       ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr);
1658       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
1659     }
1660     if (pc->ops->view) {
1661       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1662     }
1663   } else if (isdraw) {
1664     PetscDraw draw;
1665     char      str[25];
1666     PetscReal x,y,bottom,h;
1667     PetscInt  n;
1668 
1669     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1670     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1671     if (pc->mat) {
1672       ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr);
1673       ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr);
1674     } else {
1675       ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr);
1676     }
1677     ierr   = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1678     bottom = y - h;
1679     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1680     if (pc->ops->view) {
1681       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1682     }
1683     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1684 #if defined(PETSC_HAVE_SAWS)
1685   } else if (issaws) {
1686     PetscMPIInt rank;
1687 
1688     ierr = PetscObjectName((PetscObject)pc);CHKERRQ(ierr);
1689     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1690     if (!((PetscObject)pc)->amsmem && !rank) {
1691       ierr = PetscObjectViewSAWs((PetscObject)pc,viewer);CHKERRQ(ierr);
1692     }
1693     if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);}
1694     if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1695 #endif
1696   }
1697   PetscFunctionReturn(0);
1698 }
1699 
1700 /*@C
1701   PCRegister -  Adds a method to the preconditioner package.
1702 
1703    Not collective
1704 
1705    Input Parameters:
1706 +  name_solver - name of a new user-defined solver
1707 -  routine_create - routine to create method context
1708 
1709    Notes:
1710    PCRegister() may be called multiple times to add several user-defined preconditioners.
1711 
1712    Sample usage:
1713 .vb
1714    PCRegister("my_solver", MySolverCreate);
1715 .ve
1716 
1717    Then, your solver can be chosen with the procedural interface via
1718 $     PCSetType(pc,"my_solver")
1719    or at runtime via the option
1720 $     -pc_type my_solver
1721 
1722    Level: advanced
1723 
1724 .seealso: PCRegisterAll()
1725 @*/
1726 PetscErrorCode  PCRegister(const char sname[],PetscErrorCode (*function)(PC))
1727 {
1728   PetscErrorCode ierr;
1729 
1730   PetscFunctionBegin;
1731   ierr = PCInitializePackage();CHKERRQ(ierr);
1732   ierr = PetscFunctionListAdd(&PCList,sname,function);CHKERRQ(ierr);
1733   PetscFunctionReturn(0);
1734 }
1735 
1736 static PetscErrorCode MatMult_PC(Mat A,Vec X,Vec Y)
1737 {
1738   PC             pc;
1739   PetscErrorCode ierr;
1740 
1741   PetscFunctionBegin;
1742   ierr = MatShellGetContext(A,&pc);CHKERRQ(ierr);
1743   ierr = PCApply(pc,X,Y);CHKERRQ(ierr);
1744   PetscFunctionReturn(0);
1745 }
1746 
1747 /*@
1748     PCComputeOperator - Computes the explicit preconditioned operator.
1749 
1750     Collective on PC
1751 
1752     Input Parameter:
1753 +   pc - the preconditioner object
1754 -   mattype - the matrix type to be used for the operator
1755 
1756     Output Parameter:
1757 .   mat - the explict preconditioned operator
1758 
1759     Notes:
1760     This computation is done by applying the operators to columns of the identity matrix.
1761     This routine is costly in general, and is recommended for use only with relatively small systems.
1762     Currently, this routine uses a dense matrix format when mattype == NULL
1763 
1764     Level: advanced
1765 
1766 .seealso: KSPComputeOperator(), MatType
1767 
1768 @*/
1769 PetscErrorCode  PCComputeOperator(PC pc,MatType mattype,Mat *mat)
1770 {
1771   PetscErrorCode ierr;
1772   PetscInt       N,M,m,n;
1773   Mat            A,Apc;
1774 
1775   PetscFunctionBegin;
1776   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1777   PetscValidPointer(mat,3);
1778   ierr = PCGetOperators(pc,&A,NULL);CHKERRQ(ierr);
1779   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1780   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1781   ierr = MatCreateShell(PetscObjectComm((PetscObject)pc),m,n,M,N,pc,&Apc);CHKERRQ(ierr);
1782   ierr = MatShellSetOperation(Apc,MATOP_MULT,(void (*)(void))MatMult_PC);CHKERRQ(ierr);
1783   ierr = MatComputeOperator(Apc,mattype,mat);CHKERRQ(ierr);
1784   ierr = MatDestroy(&Apc);CHKERRQ(ierr);
1785   PetscFunctionReturn(0);
1786 }
1787 
1788 /*@
1789    PCSetCoordinates - sets the coordinates of all the nodes on the local process
1790 
1791    Collective on PC
1792 
1793    Input Parameters:
1794 +  pc - the solver context
1795 .  dim - the dimension of the coordinates 1, 2, or 3
1796 .  nloc - the blocked size of the coordinates array
1797 -  coords - the coordinates array
1798 
1799    Level: intermediate
1800 
1801    Notes:
1802    coords is an array of the dim coordinates for the nodes on
1803    the local processor, of size dim*nloc.
1804    If there are 108 equation on a processor
1805    for a displacement finite element discretization of elasticity (so
1806    that there are nloc = 36 = 108/3 nodes) then the array must have 108
1807    double precision values (ie, 3 * 36).  These x y z coordinates
1808    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1809    ... , N-1.z ].
1810 
1811 .seealso: MatSetNearNullSpace()
1812 @*/
1813 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1814 {
1815   PetscErrorCode ierr;
1816 
1817   PetscFunctionBegin;
1818   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1819   PetscValidLogicalCollectiveInt(pc,dim,2);
1820   ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr);
1821   PetscFunctionReturn(0);
1822 }
1823 
1824 /*@
1825    PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1826 
1827    Logically Collective on PC
1828 
1829    Input Parameters:
1830 +  pc - the precondition context
1831 
1832    Output Parameter:
1833 -  num_levels - the number of levels
1834 .  interpolations - the interpolation matrices (size of num_levels-1)
1835 
1836    Level: advanced
1837 
1838 .keywords: MG, GAMG, BoomerAMG, multigrid, interpolation, level
1839 
1840 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetInterpolation(), PCGetCoarseOperators()
1841 @*/
1842 PetscErrorCode PCGetInterpolations(PC pc,PetscInt *num_levels,Mat *interpolations[])
1843 {
1844   PetscErrorCode ierr;
1845 
1846   PetscFunctionBegin;
1847   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1848   PetscValidPointer(num_levels,2);
1849   PetscValidPointer(interpolations,3);
1850   ierr = PetscUseMethod(pc,"PCGetInterpolations_C",(PC,PetscInt*,Mat*[]),(pc,num_levels,interpolations));CHKERRQ(ierr);
1851   PetscFunctionReturn(0);
1852 }
1853 
1854 /*@
1855    PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
1856 
1857    Logically Collective on PC
1858 
1859    Input Parameters:
1860 +  pc - the precondition context
1861 
1862    Output Parameter:
1863 -  num_levels - the number of levels
1864 .  coarseOperators - the coarse operator matrices (size of num_levels-1)
1865 
1866    Level: advanced
1867 
1868 .keywords: MG, GAMG, BoomerAMG, get, multigrid, interpolation, level
1869 
1870 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale(), PCMGGetInterpolation(), PCGetInterpolations()
1871 @*/
1872 PetscErrorCode PCGetCoarseOperators(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1873 {
1874   PetscErrorCode ierr;
1875 
1876   PetscFunctionBegin;
1877   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1878   PetscValidPointer(num_levels,2);
1879   PetscValidPointer(coarseOperators,3);
1880   ierr = PetscUseMethod(pc,"PCGetCoarseOperators_C",(PC,PetscInt*,Mat*[]),(pc,num_levels,coarseOperators));CHKERRQ(ierr);
1881   PetscFunctionReturn(0);
1882 }
1883