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