xref: /petsc/src/ksp/pc/interface/precon.c (revision 23ee1639a9fcee8ab09e03ad3e2325b6ad810ca8)
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;
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   if (!pc->setupcalled) {
877     ierr            = PetscInfo(pc,"Setting up PC for first time");CHKERRQ(ierr);
878     pc->flag        = DIFFERENT_NONZERO_PATTERN;
879   } else if (matstate == pc->matstate) {
880     ierr = PetscInfo(pc,"Leaving PC with identical preconditioner since operator is unchanged\n");CHKERRQ(ierr);
881     PetscFunctionReturn(0);
882   } else {
883     PetscObjectState matnonzerostate;
884     ierr = MatGetNonzeroState(pc->pmat,&matnonzerostate);CHKERRQ(ierr);
885     if (matnonzerostate > pc->matnonzerostate) {
886        ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
887        pc->flag            = DIFFERENT_NONZERO_PATTERN;
888        pc->matnonzerostate = matnonzerostate;
889     } else {
890       ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
891       pc->flag            = SAME_NONZERO_PATTERN;
892     }
893   }
894   pc->matstate = matstate;
895 
896   if (!((PetscObject)pc)->type_name) {
897     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
898     ierr = PCSetType(pc,def);CHKERRQ(ierr);
899   }
900 
901   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
902   if (pc->ops->setup) {
903     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
904   }
905   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
906   pc->setupcalled = 1;
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "PCSetUpOnBlocks"
912 /*@
913    PCSetUpOnBlocks - Sets up the preconditioner for each block in
914    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
915    methods.
916 
917    Collective on PC
918 
919    Input Parameters:
920 .  pc - the preconditioner context
921 
922    Level: developer
923 
924 .keywords: PC, setup, blocks
925 
926 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
927 @*/
928 PetscErrorCode  PCSetUpOnBlocks(PC pc)
929 {
930   PetscErrorCode ierr;
931 
932   PetscFunctionBegin;
933   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
934   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
935   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
936   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
937   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "PCSetModifySubMatrices"
943 /*@C
944    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
945    submatrices that arise within certain subdomain-based preconditioners.
946    The basic submatrices are extracted from the preconditioner matrix as
947    usual; the user can then alter these (for example, to set different boundary
948    conditions for each submatrix) before they are used for the local solves.
949 
950    Logically Collective on PC
951 
952    Input Parameters:
953 +  pc - the preconditioner context
954 .  func - routine for modifying the submatrices
955 -  ctx - optional user-defined context (may be null)
956 
957    Calling sequence of func:
958 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
959 
960 .  row - an array of index sets that contain the global row numbers
961          that comprise each local submatrix
962 .  col - an array of index sets that contain the global column numbers
963          that comprise each local submatrix
964 .  submat - array of local submatrices
965 -  ctx - optional user-defined context for private data for the
966          user-defined func routine (may be null)
967 
968    Notes:
969    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
970    KSPSolve().
971 
972    A routine set by PCSetModifySubMatrices() is currently called within
973    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
974    preconditioners.  All other preconditioners ignore this routine.
975 
976    Level: advanced
977 
978 .keywords: PC, set, modify, submatrices
979 
980 .seealso: PCModifySubMatrices(), PCASMGetSubMatrices()
981 @*/
982 PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
983 {
984   PetscFunctionBegin;
985   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
986   pc->modifysubmatrices  = func;
987   pc->modifysubmatricesP = ctx;
988   PetscFunctionReturn(0);
989 }
990 
991 #undef __FUNCT__
992 #define __FUNCT__ "PCModifySubMatrices"
993 /*@C
994    PCModifySubMatrices - Calls an optional user-defined routine within
995    certain preconditioners if one has been set with PCSetModifySubMarices().
996 
997    Collective on PC
998 
999    Input Parameters:
1000 +  pc - the preconditioner context
1001 .  nsub - the number of local submatrices
1002 .  row - an array of index sets that contain the global row numbers
1003          that comprise each local submatrix
1004 .  col - an array of index sets that contain the global column numbers
1005          that comprise each local submatrix
1006 .  submat - array of local submatrices
1007 -  ctx - optional user-defined context for private data for the
1008          user-defined routine (may be null)
1009 
1010    Output Parameter:
1011 .  submat - array of local submatrices (the entries of which may
1012             have been modified)
1013 
1014    Notes:
1015    The user should NOT generally call this routine, as it will
1016    automatically be called within certain preconditioners (currently
1017    block Jacobi, additive Schwarz) if set.
1018 
1019    The basic submatrices are extracted from the preconditioner matrix
1020    as usual; the user can then alter these (for example, to set different
1021    boundary conditions for each submatrix) before they are used for the
1022    local solves.
1023 
1024    Level: developer
1025 
1026 .keywords: PC, modify, submatrices
1027 
1028 .seealso: PCSetModifySubMatrices()
1029 @*/
1030 PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
1031 {
1032   PetscErrorCode ierr;
1033 
1034   PetscFunctionBegin;
1035   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1036   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
1037   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
1038   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
1039   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 #undef __FUNCT__
1044 #define __FUNCT__ "PCSetOperators"
1045 /*@
1046    PCSetOperators - Sets the matrix associated with the linear system and
1047    a (possibly) different one associated with the preconditioner.
1048 
1049    Logically Collective on PC and Mat
1050 
1051    Input Parameters:
1052 +  pc - the preconditioner context
1053 .  Amat - the matrix that defines the linear system
1054 .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
1055 -  flag - flag indicating information about the preconditioner matrix structure
1056    during successive linear solves.  This flag is ignored the first time a
1057    linear system is solved, and thus is irrelevant when solving just one linear
1058    system.
1059 
1060    Notes:
1061    The flag can be used to eliminate unnecessary work in the preconditioner
1062    during the repeated solution of linear systems of the same size.  The
1063    available options are
1064 +    SAME_PRECONDITIONER -
1065        Pmat is identical during successive linear solves.
1066        This option is intended for folks who are using
1067        different Amat and Pmat matrices and wish to reuse the
1068        same preconditioner matrix.  For example, this option
1069        saves work by not recomputing incomplete factorization
1070        for ILU/ICC preconditioners.
1071 .     SAME_NONZERO_PATTERN -
1072        Pmat has the same nonzero structure during
1073        successive linear solves.
1074 -     DIFFERENT_NONZERO_PATTERN -
1075        Pmat does not have the same nonzero structure.
1076 
1077     Passing a NULL for Amat or Pmat removes the matrix that is currently used.
1078 
1079     If you wish to replace either Amat or Pmat but leave the other one untouched then
1080     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
1081     on it and then pass it back in in your call to KSPSetOperators().
1082 
1083    Caution:
1084    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
1085    and does not check the structure of the matrix.  If you erroneously
1086    claim that the structure is the same when it actually is not, the new
1087    preconditioner will not function correctly.  Thus, use this optimization
1088    feature carefully!
1089 
1090    If in doubt about whether your preconditioner matrix has changed
1091    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
1092 
1093    More Notes about Repeated Solution of Linear Systems:
1094    PETSc does NOT reset the matrix entries of either Amat or Pmat
1095    to zero after a linear solve; the user is completely responsible for
1096    matrix assembly.  See the routine MatZeroEntries() if desiring to
1097    zero all elements of a matrix.
1098 
1099    Level: intermediate
1100 
1101 .keywords: PC, set, operators, matrix, linear system
1102 
1103 .seealso: PCGetOperators(), MatZeroEntries()
1104  @*/
1105 PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat)
1106 {
1107   PetscErrorCode   ierr;
1108   PetscInt         m1,n1,m2,n2;
1109 
1110   PetscFunctionBegin;
1111   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1112   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1113   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1114   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1115   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1116   if (pc->setupcalled && Amat && Pmat) {
1117     ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr);
1118     ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr);
1119     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);
1120     ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr);
1121     ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr);
1122     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);
1123   }
1124 
1125   if (Pmat != pc->pmat) {
1126     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1127     pc->matnonzerostate = -1;
1128     pc->matstate        = -1;
1129   }
1130 
1131   /* reference first in case the matrices are the same */
1132   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1133   ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1134   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1135   ierr     = MatDestroy(&pc->pmat);CHKERRQ(ierr);
1136   pc->mat  = Amat;
1137   pc->pmat = Pmat;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "PCSetReusePreconditioner"
1143 /*@
1144    PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
1145 
1146    Logically Collective on PC
1147 
1148    Input Parameters:
1149 +  pc - the preconditioner context
1150 -  flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1151 
1152 .seealso: PCGetOperators(), MatZeroEntries()
1153  @*/
1154 PetscErrorCode  PCSetReusePreconditioner(PC pc,PetscBool flag)
1155 {
1156   PetscFunctionBegin;
1157   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1158   pc->reusepreconditioner = flag;
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 #undef __FUNCT__
1163 #define __FUNCT__ "PCGetOperators"
1164 /*@C
1165    PCGetOperators - Gets the matrix associated with the linear system and
1166    possibly a different one associated with the preconditioner.
1167 
1168    Not collective, though parallel Mats are returned if the PC is parallel
1169 
1170    Input Parameter:
1171 .  pc - the preconditioner context
1172 
1173    Output Parameters:
1174 +  Amat - the matrix defining the linear system
1175 -  Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1176 
1177    Level: intermediate
1178 
1179    Notes: Does not increase the reference count of the matrices, so you should not destroy them
1180 
1181    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1182       are created in PC and returned to the user. In this case, if both operators
1183       mat and pmat are requested, two DIFFERENT operators will be returned. If
1184       only one is requested both operators in the PC will be the same (i.e. as
1185       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1186       The user must set the sizes of the returned matrices and their type etc just
1187       as if the user created them with MatCreate(). For example,
1188 
1189 $         KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1190 $           set size, type, etc of Amat
1191 
1192 $         MatCreate(comm,&mat);
1193 $         KSP/PCSetOperators(ksp/pc,Amat,Amat);
1194 $         PetscObjectDereference((PetscObject)mat);
1195 $           set size, type, etc of Amat
1196 
1197      and
1198 
1199 $         KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1200 $           set size, type, etc of Amat and Pmat
1201 
1202 $         MatCreate(comm,&Amat);
1203 $         MatCreate(comm,&Pmat);
1204 $         KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1205 $         PetscObjectDereference((PetscObject)Amat);
1206 $         PetscObjectDereference((PetscObject)Pmat);
1207 $           set size, type, etc of Amat and Pmat
1208 
1209     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1210     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1211     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1212     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1213     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1214     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1215     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1216     it can be created for you?
1217 
1218 
1219 .keywords: PC, get, operators, matrix, linear system
1220 
1221 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1222 @*/
1223 PetscErrorCode  PCGetOperators(PC pc,Mat *Amat,Mat *Pmat)
1224 {
1225   PetscErrorCode ierr;
1226 
1227   PetscFunctionBegin;
1228   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1229   if (Amat) {
1230     if (!pc->mat) {
1231       if (pc->pmat && !Pmat) {  /* Apmat has been set, but user did not request it, so use for Amat */
1232         pc->mat = pc->pmat;
1233         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1234       } else {                  /* both Amat and Pmat are empty */
1235         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr);
1236         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1237           pc->pmat = pc->mat;
1238           ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1239         }
1240       }
1241     }
1242     *Amat = pc->mat;
1243   }
1244   if (Pmat) {
1245     if (!pc->pmat) {
1246       if (pc->mat && !Amat) {    /* Amat has been set but was not requested, so use for pmat */
1247         pc->pmat = pc->mat;
1248         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1249       } else {
1250         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr);
1251         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1252           pc->mat = pc->pmat;
1253           ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1254         }
1255       }
1256     }
1257     *Pmat = pc->pmat;
1258   }
1259   PetscFunctionReturn(0);
1260 }
1261 
1262 #undef __FUNCT__
1263 #define __FUNCT__ "PCGetOperatorsSet"
1264 /*@C
1265    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1266    possibly a different one associated with the preconditioner have been set in the PC.
1267 
1268    Not collective, though the results on all processes should be the same
1269 
1270    Input Parameter:
1271 .  pc - the preconditioner context
1272 
1273    Output Parameters:
1274 +  mat - the matrix associated with the linear system was set
1275 -  pmat - matrix associated with the preconditioner was set, usually the same
1276 
1277    Level: intermediate
1278 
1279 .keywords: PC, get, operators, matrix, linear system
1280 
1281 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1282 @*/
1283 PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1284 {
1285   PetscFunctionBegin;
1286   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1287   if (mat) *mat = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1288   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1289   PetscFunctionReturn(0);
1290 }
1291 
1292 #undef __FUNCT__
1293 #define __FUNCT__ "PCFactorGetMatrix"
1294 /*@
1295    PCFactorGetMatrix - Gets the factored matrix from the
1296    preconditioner context.  This routine is valid only for the LU,
1297    incomplete LU, Cholesky, and incomplete Cholesky methods.
1298 
1299    Not Collective on PC though Mat is parallel if PC is parallel
1300 
1301    Input Parameters:
1302 .  pc - the preconditioner context
1303 
1304    Output parameters:
1305 .  mat - the factored matrix
1306 
1307    Level: advanced
1308 
1309    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1310 
1311 .keywords: PC, get, factored, matrix
1312 @*/
1313 PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1314 {
1315   PetscErrorCode ierr;
1316 
1317   PetscFunctionBegin;
1318   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1319   PetscValidPointer(mat,2);
1320   if (pc->ops->getfactoredmatrix) {
1321     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1322   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1323   PetscFunctionReturn(0);
1324 }
1325 
1326 #undef __FUNCT__
1327 #define __FUNCT__ "PCSetOptionsPrefix"
1328 /*@C
1329    PCSetOptionsPrefix - Sets the prefix used for searching for all
1330    PC options in the database.
1331 
1332    Logically Collective on PC
1333 
1334    Input Parameters:
1335 +  pc - the preconditioner context
1336 -  prefix - the prefix string to prepend to all PC option requests
1337 
1338    Notes:
1339    A hyphen (-) must NOT be given at the beginning of the prefix name.
1340    The first character of all runtime options is AUTOMATICALLY the
1341    hyphen.
1342 
1343    Level: advanced
1344 
1345 .keywords: PC, set, options, prefix, database
1346 
1347 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1348 @*/
1349 PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1350 {
1351   PetscErrorCode ierr;
1352 
1353   PetscFunctionBegin;
1354   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1355   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1356   PetscFunctionReturn(0);
1357 }
1358 
1359 #undef __FUNCT__
1360 #define __FUNCT__ "PCAppendOptionsPrefix"
1361 /*@C
1362    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1363    PC options in the database.
1364 
1365    Logically Collective on PC
1366 
1367    Input Parameters:
1368 +  pc - the preconditioner context
1369 -  prefix - the prefix string to prepend to all PC option requests
1370 
1371    Notes:
1372    A hyphen (-) must NOT be given at the beginning of the prefix name.
1373    The first character of all runtime options is AUTOMATICALLY the
1374    hyphen.
1375 
1376    Level: advanced
1377 
1378 .keywords: PC, append, options, prefix, database
1379 
1380 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1381 @*/
1382 PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1383 {
1384   PetscErrorCode ierr;
1385 
1386   PetscFunctionBegin;
1387   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1388   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 #undef __FUNCT__
1393 #define __FUNCT__ "PCGetOptionsPrefix"
1394 /*@C
1395    PCGetOptionsPrefix - Gets the prefix used for searching for all
1396    PC options in the database.
1397 
1398    Not Collective
1399 
1400    Input Parameters:
1401 .  pc - the preconditioner context
1402 
1403    Output Parameters:
1404 .  prefix - pointer to the prefix string used, is returned
1405 
1406    Notes: On the fortran side, the user should pass in a string 'prifix' of
1407    sufficient length to hold the prefix.
1408 
1409    Level: advanced
1410 
1411 .keywords: PC, get, options, prefix, database
1412 
1413 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1414 @*/
1415 PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1416 {
1417   PetscErrorCode ierr;
1418 
1419   PetscFunctionBegin;
1420   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1421   PetscValidPointer(prefix,2);
1422   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 #undef __FUNCT__
1427 #define __FUNCT__ "PCPreSolve"
1428 /*@
1429    PCPreSolve - Optional pre-solve phase, intended for any
1430    preconditioner-specific actions that must be performed before
1431    the iterative solve itself.
1432 
1433    Collective on PC
1434 
1435    Input Parameters:
1436 +  pc - the preconditioner context
1437 -  ksp - the Krylov subspace context
1438 
1439    Level: developer
1440 
1441    Sample of Usage:
1442 .vb
1443     PCPreSolve(pc,ksp);
1444     KSPSolve(ksp,b,x);
1445     PCPostSolve(pc,ksp);
1446 .ve
1447 
1448    Notes:
1449    The pre-solve phase is distinct from the PCSetUp() phase.
1450 
1451    KSPSolve() calls this directly, so is rarely called by the user.
1452 
1453 .keywords: PC, pre-solve
1454 
1455 .seealso: PCPostSolve()
1456 @*/
1457 PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1458 {
1459   PetscErrorCode ierr;
1460   Vec            x,rhs;
1461 
1462   PetscFunctionBegin;
1463   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1464   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1465   pc->presolvedone++;
1466   if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice");
1467   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1468   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1469 
1470   if (pc->ops->presolve) {
1471     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1472   }
1473   PetscFunctionReturn(0);
1474 }
1475 
1476 #undef __FUNCT__
1477 #define __FUNCT__ "PCPostSolve"
1478 /*@
1479    PCPostSolve - Optional post-solve phase, intended for any
1480    preconditioner-specific actions that must be performed after
1481    the iterative solve itself.
1482 
1483    Collective on PC
1484 
1485    Input Parameters:
1486 +  pc - the preconditioner context
1487 -  ksp - the Krylov subspace context
1488 
1489    Sample of Usage:
1490 .vb
1491     PCPreSolve(pc,ksp);
1492     KSPSolve(ksp,b,x);
1493     PCPostSolve(pc,ksp);
1494 .ve
1495 
1496    Note:
1497    KSPSolve() calls this routine directly, so it is rarely called by the user.
1498 
1499    Level: developer
1500 
1501 .keywords: PC, post-solve
1502 
1503 .seealso: PCPreSolve(), KSPSolve()
1504 @*/
1505 PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1506 {
1507   PetscErrorCode ierr;
1508   Vec            x,rhs;
1509 
1510   PetscFunctionBegin;
1511   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1512   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1513   pc->presolvedone--;
1514   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1515   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1516   if (pc->ops->postsolve) {
1517     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1518   }
1519   PetscFunctionReturn(0);
1520 }
1521 
1522 #undef __FUNCT__
1523 #define __FUNCT__ "PCLoad"
1524 /*@C
1525   PCLoad - Loads a PC that has been stored in binary  with PCView().
1526 
1527   Collective on PetscViewer
1528 
1529   Input Parameters:
1530 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or
1531            some related function before a call to PCLoad().
1532 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1533 
1534    Level: intermediate
1535 
1536   Notes:
1537    The type is determined by the data in the file, any type set into the PC before this call is ignored.
1538 
1539   Notes for advanced users:
1540   Most users should not need to know the details of the binary storage
1541   format, since PCLoad() and PCView() completely hide these details.
1542   But for anyone who's interested, the standard binary matrix storage
1543   format is
1544 .vb
1545      has not yet been determined
1546 .ve
1547 
1548 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad()
1549 @*/
1550 PetscErrorCode  PCLoad(PC newdm, PetscViewer viewer)
1551 {
1552   PetscErrorCode ierr;
1553   PetscBool      isbinary;
1554   PetscInt       classid;
1555   char           type[256];
1556 
1557   PetscFunctionBegin;
1558   PetscValidHeaderSpecific(newdm,PC_CLASSID,1);
1559   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1560   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1561   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1562 
1563   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
1564   if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file");
1565   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
1566   ierr = PCSetType(newdm, type);CHKERRQ(ierr);
1567   if (newdm->ops->load) {
1568     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1569   }
1570   PetscFunctionReturn(0);
1571 }
1572 
1573 #include <petscdraw.h>
1574 #if defined(PETSC_HAVE_SAWS)
1575 #include <petscviewersaws.h>
1576 #endif
1577 #undef __FUNCT__
1578 #define __FUNCT__ "PCView"
1579 /*@C
1580    PCView - Prints the PC data structure.
1581 
1582    Collective on PC
1583 
1584    Input Parameters:
1585 +  PC - the PC context
1586 -  viewer - optional visualization context
1587 
1588    Note:
1589    The available visualization contexts include
1590 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1591 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1592          output where only the first processor opens
1593          the file.  All other processors send their
1594          data to the first processor to print.
1595 
1596    The user can open an alternative visualization contexts with
1597    PetscViewerASCIIOpen() (output to a specified file).
1598 
1599    Level: developer
1600 
1601 .keywords: PC, view
1602 
1603 .seealso: KSPView(), PetscViewerASCIIOpen()
1604 @*/
1605 PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1606 {
1607   PCType            cstr;
1608   PetscErrorCode    ierr;
1609   PetscBool         iascii,isstring,isbinary,isdraw;
1610   PetscViewerFormat format;
1611 #if defined(PETSC_HAVE_SAWS)
1612   PetscBool         isams;
1613 #endif
1614 
1615   PetscFunctionBegin;
1616   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1617   if (!viewer) {
1618     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
1619   }
1620   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1621   PetscCheckSameComm(pc,1,viewer,2);
1622 
1623   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1624   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1625   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1626   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1627 #if defined(PETSC_HAVE_SAWS)
1628   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&isams);CHKERRQ(ierr);
1629 #endif
1630 
1631   if (iascii) {
1632     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1633     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);CHKERRQ(ierr);
1634     if (!pc->setupcalled) {
1635       ierr = PetscViewerASCIIPrintf(viewer,"  PC has not been set up so information may be incomplete\n");CHKERRQ(ierr);
1636     }
1637     if (pc->ops->view) {
1638       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1639       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1640       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1641     }
1642     if (pc->mat) {
1643       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1644       if (pc->pmat == pc->mat) {
1645         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1646         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1647         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1648         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1649       } else {
1650         if (pc->pmat) {
1651           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1652         } else {
1653           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1654         }
1655         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1656         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1657         if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1658         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1659       }
1660       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1661     }
1662   } else if (isstring) {
1663     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1664     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1665     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1666   } else if (isbinary) {
1667     PetscInt    classid = PC_FILE_CLASSID;
1668     MPI_Comm    comm;
1669     PetscMPIInt rank;
1670     char        type[256];
1671 
1672     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1673     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1674     if (!rank) {
1675       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1676       ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr);
1677       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
1678     }
1679     if (pc->ops->view) {
1680       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1681     }
1682   } else if (isdraw) {
1683     PetscDraw draw;
1684     char      str[25];
1685     PetscReal x,y,bottom,h;
1686     PetscInt  n;
1687 
1688     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1689     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1690     if (pc->mat) {
1691       ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr);
1692       ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr);
1693     } else {
1694       ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr);
1695     }
1696     ierr   = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1697     bottom = y - h;
1698     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1699     if (pc->ops->view) {
1700       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1701     }
1702     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1703 #if defined(PETSC_HAVE_SAWS)
1704   } else if (isams) {
1705     PetscMPIInt rank;
1706 
1707     ierr = PetscObjectName((PetscObject)pc);CHKERRQ(ierr);
1708     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1709     if (!((PetscObject)pc)->amsmem && !rank) {
1710       ierr = PetscObjectViewSAWs((PetscObject)pc,viewer);CHKERRQ(ierr);
1711     }
1712     if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);}
1713     if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1714 #endif
1715   }
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 
1720 #undef __FUNCT__
1721 #define __FUNCT__ "PCSetInitialGuessNonzero"
1722 /*@
1723    PCSetInitialGuessNonzero - Tells the iterative solver that the
1724    initial guess is nonzero; otherwise PC assumes the initial guess
1725    is to be zero (and thus zeros it out before solving).
1726 
1727    Logically Collective on PC
1728 
1729    Input Parameters:
1730 +  pc - iterative context obtained from PCCreate()
1731 -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1732 
1733    Level: Developer
1734 
1735    Notes:
1736     This is a weird function. Since PC's are linear operators on the right hand side they
1737     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1738     PCKSP and PCREDUNDANT  and causes the inner KSP object to use the nonzero
1739     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1740 
1741 
1742 .keywords: PC, set, initial guess, nonzero
1743 
1744 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1745 @*/
1746 PetscErrorCode  PCSetInitialGuessNonzero(PC pc,PetscBool flg)
1747 {
1748   PetscFunctionBegin;
1749   PetscValidLogicalCollectiveBool(pc,flg,2);
1750   pc->nonzero_guess = flg;
1751   PetscFunctionReturn(0);
1752 }
1753 
1754 #undef __FUNCT__
1755 #define __FUNCT__ "PCRegister"
1756 /*@C
1757   PCRegister -  Adds a method to the preconditioner package.
1758 
1759    Not collective
1760 
1761    Input Parameters:
1762 +  name_solver - name of a new user-defined solver
1763 -  routine_create - routine to create method context
1764 
1765    Notes:
1766    PCRegister() may be called multiple times to add several user-defined preconditioners.
1767 
1768    Sample usage:
1769 .vb
1770    PCRegister("my_solver", MySolverCreate);
1771 .ve
1772 
1773    Then, your solver can be chosen with the procedural interface via
1774 $     PCSetType(pc,"my_solver")
1775    or at runtime via the option
1776 $     -pc_type my_solver
1777 
1778    Level: advanced
1779 
1780 .keywords: PC, register
1781 
1782 .seealso: PCRegisterAll(), PCRegisterDestroy()
1783 @*/
1784 PetscErrorCode  PCRegister(const char sname[],PetscErrorCode (*function)(PC))
1785 {
1786   PetscErrorCode ierr;
1787 
1788   PetscFunctionBegin;
1789   ierr = PetscFunctionListAdd(&PCList,sname,function);CHKERRQ(ierr);
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 #undef __FUNCT__
1794 #define __FUNCT__ "PCComputeExplicitOperator"
1795 /*@
1796     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1797 
1798     Collective on PC
1799 
1800     Input Parameter:
1801 .   pc - the preconditioner object
1802 
1803     Output Parameter:
1804 .   mat - the explict preconditioned operator
1805 
1806     Notes:
1807     This computation is done by applying the operators to columns of the
1808     identity matrix.
1809 
1810     Currently, this routine uses a dense matrix format when 1 processor
1811     is used and a sparse format otherwise.  This routine is costly in general,
1812     and is recommended for use only with relatively small systems.
1813 
1814     Level: advanced
1815 
1816 .keywords: PC, compute, explicit, operator
1817 
1818 .seealso: KSPComputeExplicitOperator()
1819 
1820 @*/
1821 PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1822 {
1823   Vec            in,out;
1824   PetscErrorCode ierr;
1825   PetscInt       i,M,m,*rows,start,end;
1826   PetscMPIInt    size;
1827   MPI_Comm       comm;
1828   PetscScalar    *array,one = 1.0;
1829 
1830   PetscFunctionBegin;
1831   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1832   PetscValidPointer(mat,2);
1833 
1834   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1835   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1836 
1837   if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1838   ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1839   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1840   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1841   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1842   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1843   ierr = PetscMalloc1((m+1),&rows);CHKERRQ(ierr);
1844   for (i=0; i<m; i++) rows[i] = start + i;
1845 
1846   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
1847   ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr);
1848   if (size == 1) {
1849     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1850     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
1851   } else {
1852     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1853     ierr = MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);CHKERRQ(ierr);
1854   }
1855 
1856   for (i=0; i<M; i++) {
1857 
1858     ierr = VecSet(in,0.0);CHKERRQ(ierr);
1859     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1860     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1861     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1862 
1863     /* should fix, allowing user to choose side */
1864     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1865 
1866     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1867     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1868     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1869 
1870   }
1871   ierr = PetscFree(rows);CHKERRQ(ierr);
1872   ierr = VecDestroy(&out);CHKERRQ(ierr);
1873   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1874   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1875   PetscFunctionReturn(0);
1876 }
1877 
1878 #undef __FUNCT__
1879 #define __FUNCT__ "PCSetCoordinates"
1880 /*@
1881    PCSetCoordinates - sets the coordinates of all the nodes on the local process
1882 
1883    Collective on PC
1884 
1885    Input Parameters:
1886 +  pc - the solver context
1887 .  dim - the dimension of the coordinates 1, 2, or 3
1888 -  coords - the coordinates
1889 
1890    Level: intermediate
1891 
1892    Notes: coords is an array of the 3D coordinates for the nodes on
1893    the local processor.  So if there are 108 equation on a processor
1894    for a displacement finite element discretization of elasticity (so
1895    that there are 36 = 108/3 nodes) then the array must have 108
1896    double precision values (ie, 3 * 36).  These x y z coordinates
1897    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1898    ... , N-1.z ].
1899 
1900 .seealso: MatSetNearNullSpace
1901 @*/
1902 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1903 {
1904   PetscErrorCode ierr;
1905 
1906   PetscFunctionBegin;
1907   ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr);
1908   PetscFunctionReturn(0);
1909 }
1910