xref: /petsc/src/ksp/pc/interface/precon.c (revision b90c6cbe8a28dbbfeb8633979a88a1769a41a59e)
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",(void (**)(void))&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 AMS then destroy it */
120   ierr = PetscObjectDepublish((*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 #if 0
289 #undef __FUNCT__
290 #define __FUNCT__ "PCPublish_Petsc"
291 static PetscErrorCode PCPublish_Petsc(PetscObject obj)
292 {
293   PetscFunctionBegin;
294   PetscFunctionReturn(0);
295 }
296 #endif
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "PCSetUseAmat"
300 /*@
301    PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
302    operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
303    TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
304 
305    Logically Collective on PC
306 
307    Input Parameters:
308 +  pc - the preconditioner context
309 -  flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
310 
311    Options Database Key:
312 .  -pc_use_amat <true,false>
313 
314    Notes:
315    For the common case in which the linear system matrix and the matrix used to construct the
316    preconditioner are identical, this routine is does nothing.
317 
318    Level: intermediate
319 
320 .seealso: PCGetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
321 @*/
322 PetscErrorCode  PCSetUseAmat(PC pc,PetscBool flg)
323 {
324   PetscFunctionBegin;
325   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
326   pc->useAmat = flg;
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "PCGetUseAmat"
332 /*@
333    PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
334    operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
335    TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
336 
337    Logically Collective on PC
338 
339    Input Parameter:
340 .  pc - the preconditioner context
341 
342    Output Parameter:
343 .  flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
344 
345    Notes:
346    For the common case in which the linear system matrix and the matrix used to construct the
347    preconditioner are identical, this routine is does nothing.
348 
349    Level: intermediate
350 
351 .seealso: PCSetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
352 @*/
353 PetscErrorCode  PCGetUseAmat(PC pc,PetscBool *flg)
354 {
355   PetscFunctionBegin;
356   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
357   *flg = pc->useAmat;
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "PCCreate"
363 /*@
364    PCCreate - Creates a preconditioner context.
365 
366    Collective on MPI_Comm
367 
368    Input Parameter:
369 .  comm - MPI communicator
370 
371    Output Parameter:
372 .  pc - location to put the preconditioner context
373 
374    Notes:
375    The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
376    in parallel. For dense matrices it is always PCNONE.
377 
378    Level: developer
379 
380 .keywords: PC, create, context
381 
382 .seealso: PCSetUp(), PCApply(), PCDestroy()
383 @*/
384 PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
385 {
386   PC             pc;
387   PetscErrorCode ierr;
388 
389   PetscFunctionBegin;
390   PetscValidPointer(newpc,1);
391   *newpc = 0;
392 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
393   ierr = PCInitializePackage(NULL);CHKERRQ(ierr);
394 #endif
395 
396   ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,"PC","Preconditioner","PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
397 
398   pc->mat                  = 0;
399   pc->pmat                 = 0;
400   pc->setupcalled          = 0;
401   pc->setfromoptionscalled = 0;
402   pc->data                 = 0;
403   pc->diagonalscale        = PETSC_FALSE;
404   pc->diagonalscaleleft    = 0;
405   pc->diagonalscaleright   = 0;
406   pc->reuse                = 0;
407 
408   pc->modifysubmatrices  = 0;
409   pc->modifysubmatricesP = 0;
410 
411   *newpc = pc;
412   PetscFunctionReturn(0);
413 
414 }
415 
416 /* -------------------------------------------------------------------------------*/
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "PCApply"
420 /*@
421    PCApply - Applies the preconditioner to a vector.
422 
423    Collective on PC and Vec
424 
425    Input Parameters:
426 +  pc - the preconditioner context
427 -  x - input vector
428 
429    Output Parameter:
430 .  y - output vector
431 
432    Level: developer
433 
434 .keywords: PC, apply
435 
436 .seealso: PCApplyTranspose(), PCApplyBAorAB()
437 @*/
438 PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
439 {
440   PetscErrorCode ierr;
441 
442   PetscFunctionBegin;
443   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
444   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
445   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
446   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
447   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
448   if (pc->setupcalled < 2) {
449     ierr = PCSetUp(pc);CHKERRQ(ierr);
450   }
451   if (!pc->ops->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply");
452   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
453   ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr);
454   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
455   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
456   PetscFunctionReturn(0);
457 }
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "PCApplySymmetricLeft"
461 /*@
462    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
463 
464    Collective on PC and Vec
465 
466    Input Parameters:
467 +  pc - the preconditioner context
468 -  x - input vector
469 
470    Output Parameter:
471 .  y - output vector
472 
473    Notes:
474    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
475 
476    Level: developer
477 
478 .keywords: PC, apply, symmetric, left
479 
480 .seealso: PCApply(), PCApplySymmetricRight()
481 @*/
482 PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
483 {
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
488   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
489   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
490   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
491   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
492   if (pc->setupcalled < 2) {
493     ierr = PCSetUp(pc);CHKERRQ(ierr);
494   }
495   if (!pc->ops->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
496   ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
497   ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr);
498   ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
499   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "PCApplySymmetricRight"
505 /*@
506    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
507 
508    Collective on PC and Vec
509 
510    Input Parameters:
511 +  pc - the preconditioner context
512 -  x - input vector
513 
514    Output Parameter:
515 .  y - output vector
516 
517    Level: developer
518 
519    Notes:
520    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
521 
522 .keywords: PC, apply, symmetric, right
523 
524 .seealso: PCApply(), PCApplySymmetricLeft()
525 @*/
526 PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
527 {
528   PetscErrorCode ierr;
529 
530   PetscFunctionBegin;
531   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
532   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
533   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
534   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
535   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
536   if (pc->setupcalled < 2) {
537     ierr = PCSetUp(pc);CHKERRQ(ierr);
538   }
539   if (!pc->ops->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
540   ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
541   ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr);
542   ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
543   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
544   PetscFunctionReturn(0);
545 }
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "PCApplyTranspose"
549 /*@
550    PCApplyTranspose - Applies the transpose of preconditioner to a vector.
551 
552    Collective on PC and Vec
553 
554    Input Parameters:
555 +  pc - the preconditioner context
556 -  x - input vector
557 
558    Output Parameter:
559 .  y - output vector
560 
561    Notes: For complex numbers this applies the non-Hermitian transpose.
562 
563    Developer Notes: We need to implement a PCApplyHermitianTranspose()
564 
565    Level: developer
566 
567 .keywords: PC, apply, transpose
568 
569 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
570 @*/
571 PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
572 {
573   PetscErrorCode ierr;
574 
575   PetscFunctionBegin;
576   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
577   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
578   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
579   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
580   ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);
581   if (pc->setupcalled < 2) {
582     ierr = PCSetUp(pc);CHKERRQ(ierr);
583   }
584   if (!pc->ops->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply transpose");
585   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
586   ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr);
587   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
588   ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 #undef __FUNCT__
593 #define __FUNCT__ "PCApplyTransposeExists"
594 /*@
595    PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
596 
597    Collective on PC and Vec
598 
599    Input Parameters:
600 .  pc - the preconditioner context
601 
602    Output Parameter:
603 .  flg - PETSC_TRUE if a transpose operation is defined
604 
605    Level: developer
606 
607 .keywords: PC, apply, transpose
608 
609 .seealso: PCApplyTranspose()
610 @*/
611 PetscErrorCode  PCApplyTransposeExists(PC pc,PetscBool  *flg)
612 {
613   PetscFunctionBegin;
614   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
615   PetscValidPointer(flg,2);
616   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
617   else *flg = PETSC_FALSE;
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "PCApplyBAorAB"
623 /*@
624    PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
625 
626    Collective on PC and Vec
627 
628    Input Parameters:
629 +  pc - the preconditioner context
630 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
631 .  x - input vector
632 -  work - work vector
633 
634    Output Parameter:
635 .  y - output vector
636 
637    Level: developer
638 
639    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
640    specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling.
641 
642 .keywords: PC, apply, operator
643 
644 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
645 @*/
646 PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
647 {
648   PetscErrorCode ierr;
649 
650   PetscFunctionBegin;
651   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
652   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
653   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
654   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
655   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
656   ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);
657   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");
658   if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
659 
660   if (pc->setupcalled < 2) {
661     ierr = PCSetUp(pc);CHKERRQ(ierr);
662   }
663 
664   if (pc->diagonalscale) {
665     if (pc->ops->applyBA) {
666       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
667       ierr = VecDuplicate(x,&work2);CHKERRQ(ierr);
668       ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr);
669       ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr);
670       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
671       ierr = VecDestroy(&work2);CHKERRQ(ierr);
672     } else if (side == PC_RIGHT) {
673       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
674       ierr = PCApply(pc,y,work);CHKERRQ(ierr);
675       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
676       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
677     } else if (side == PC_LEFT) {
678       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
679       ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr);
680       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
681       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
682     } else if (side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
683   } else {
684     if (pc->ops->applyBA) {
685       ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr);
686     } else if (side == PC_RIGHT) {
687       ierr = PCApply(pc,x,work);CHKERRQ(ierr);
688       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
689     } else if (side == PC_LEFT) {
690       ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr);
691       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
692     } else if (side == PC_SYMMETRIC) {
693       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
694       ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr);
695       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
696       ierr = VecCopy(y,work);CHKERRQ(ierr);
697       ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr);
698     }
699   }
700   ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "PCApplyBAorABTranspose"
706 /*@
707    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
708    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
709    NOT tr(B*A) = tr(A)*tr(B).
710 
711    Collective on PC and Vec
712 
713    Input Parameters:
714 +  pc - the preconditioner context
715 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
716 .  x - input vector
717 -  work - work vector
718 
719    Output Parameter:
720 .  y - output vector
721 
722 
723    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
724       defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
725 
726     Level: developer
727 
728 .keywords: PC, apply, operator, transpose
729 
730 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
731 @*/
732 PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
733 {
734   PetscErrorCode ierr;
735 
736   PetscFunctionBegin;
737   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
738   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
739   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
740   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
741   if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
742   ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);
743   if (pc->ops->applyBAtranspose) {
744     ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr);
745     ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
746     PetscFunctionReturn(0);
747   }
748   if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
749 
750   if (pc->setupcalled < 2) {
751     ierr = PCSetUp(pc);CHKERRQ(ierr);
752   }
753 
754   if (side == PC_RIGHT) {
755     ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr);
756     ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr);
757   } else if (side == PC_LEFT) {
758     ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr);
759     ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr);
760   }
761   /* add support for PC_SYMMETRIC */
762   ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);
763   PetscFunctionReturn(0);
764 }
765 
766 /* -------------------------------------------------------------------------------*/
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "PCApplyRichardsonExists"
770 /*@
771    PCApplyRichardsonExists - Determines whether a particular preconditioner has a
772    built-in fast application of Richardson's method.
773 
774    Not Collective
775 
776    Input Parameter:
777 .  pc - the preconditioner
778 
779    Output Parameter:
780 .  exists - PETSC_TRUE or PETSC_FALSE
781 
782    Level: developer
783 
784 .keywords: PC, apply, Richardson, exists
785 
786 .seealso: PCApplyRichardson()
787 @*/
788 PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscBool  *exists)
789 {
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
792   PetscValidIntPointer(exists,2);
793   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
794   else *exists = PETSC_FALSE;
795   PetscFunctionReturn(0);
796 }
797 
798 #undef __FUNCT__
799 #define __FUNCT__ "PCApplyRichardson"
800 /*@
801    PCApplyRichardson - Applies several steps of Richardson iteration with
802    the particular preconditioner. This routine is usually used by the
803    Krylov solvers and not the application code directly.
804 
805    Collective on PC
806 
807    Input Parameters:
808 +  pc  - the preconditioner context
809 .  b   - the right hand side
810 .  w   - one work vector
811 .  rtol - relative decrease in residual norm convergence criteria
812 .  abstol - absolute residual norm convergence criteria
813 .  dtol - divergence residual norm increase criteria
814 .  its - the number of iterations to apply.
815 -  guesszero - if the input x contains nonzero initial guess
816 
817    Output Parameter:
818 +  outits - number of iterations actually used (for SOR this always equals its)
819 .  reason - the reason the apply terminated
820 -  y - the solution (also contains initial guess if guesszero is PETSC_FALSE
821 
822    Notes:
823    Most preconditioners do not support this function. Use the command
824    PCApplyRichardsonExists() to determine if one does.
825 
826    Except for the multigrid PC this routine ignores the convergence tolerances
827    and always runs for the number of iterations
828 
829    Level: developer
830 
831 .keywords: PC, apply, Richardson
832 
833 .seealso: PCApplyRichardsonExists()
834 @*/
835 PetscErrorCode  PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
836 {
837   PetscErrorCode ierr;
838 
839   PetscFunctionBegin;
840   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
841   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
842   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
843   PetscValidHeaderSpecific(w,VEC_CLASSID,4);
844   if (b == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"b and y must be different vectors");
845   if (pc->setupcalled < 2) {
846     ierr = PCSetUp(pc);CHKERRQ(ierr);
847   }
848   if (!pc->ops->applyrichardson) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply richardson");
849   ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr);
850   PetscFunctionReturn(0);
851 }
852 
853 /*
854       a setupcall of 0 indicates never setup,
855                      1 needs to be resetup,
856                      2 does not need any changes.
857 */
858 #undef __FUNCT__
859 #define __FUNCT__ "PCSetUp"
860 /*@
861    PCSetUp - Prepares for the use of a preconditioner.
862 
863    Collective on PC
864 
865    Input Parameter:
866 .  pc - the preconditioner context
867 
868    Level: developer
869 
870 .keywords: PC, setup
871 
872 .seealso: PCCreate(), PCApply(), PCDestroy()
873 @*/
874 PetscErrorCode  PCSetUp(PC pc)
875 {
876   PetscErrorCode ierr;
877   const char     *def;
878 
879   PetscFunctionBegin;
880   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
881   if (!pc->mat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");
882 
883   if (pc->setupcalled > 1) {
884     ierr = PetscInfo(pc,"Setting PC with identical preconditioner\n");CHKERRQ(ierr);
885     PetscFunctionReturn(0);
886   } else if (!pc->setupcalled) {
887     ierr = PetscInfo(pc,"Setting up new PC\n");CHKERRQ(ierr);
888   } else if (pc->flag == SAME_NONZERO_PATTERN) {
889     ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
890   } else {
891     ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
892   }
893 
894   if (!((PetscObject)pc)->type_name) {
895     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
896     ierr = PCSetType(pc,def);CHKERRQ(ierr);
897   }
898 
899   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
900   if (pc->ops->setup) {
901     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
902   }
903   pc->setupcalled = 2;
904 
905   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
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 -  flag - flag indicating information about the preconditioner matrix structure
1055    during successive linear solves.  This flag is ignored the first time a
1056    linear system is solved, and thus is irrelevant when solving just one linear
1057    system.
1058 
1059    Notes:
1060    The flag can be used to eliminate unnecessary work in the preconditioner
1061    during the repeated solution of linear systems of the same size.  The
1062    available options are
1063 +    SAME_PRECONDITIONER -
1064        Pmat is identical during successive linear solves.
1065        This option is intended for folks who are using
1066        different Amat and Pmat matrices and wish to reuse the
1067        same preconditioner matrix.  For example, this option
1068        saves work by not recomputing incomplete factorization
1069        for ILU/ICC preconditioners.
1070 .     SAME_NONZERO_PATTERN -
1071        Pmat has the same nonzero structure during
1072        successive linear solves.
1073 -     DIFFERENT_NONZERO_PATTERN -
1074        Pmat does not have the same nonzero structure.
1075 
1076     Passing a NULL for Amat or Pmat removes the matrix that is currently used.
1077 
1078     If you wish to replace either Amat or Pmat but leave the other one untouched then
1079     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
1080     on it and then pass it back in in your call to KSPSetOperators().
1081 
1082    Caution:
1083    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
1084    and does not check the structure of the matrix.  If you erroneously
1085    claim that the structure is the same when it actually is not, the new
1086    preconditioner will not function correctly.  Thus, use this optimization
1087    feature carefully!
1088 
1089    If in doubt about whether your preconditioner matrix has changed
1090    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
1091 
1092    More Notes about Repeated Solution of Linear Systems:
1093    PETSc does NOT reset the matrix entries of either Amat or Pmat
1094    to zero after a linear solve; the user is completely responsible for
1095    matrix assembly.  See the routine MatZeroEntries() if desiring to
1096    zero all elements of a matrix.
1097 
1098    Level: intermediate
1099 
1100 .keywords: PC, set, operators, matrix, linear system
1101 
1102 .seealso: PCGetOperators(), MatZeroEntries()
1103  @*/
1104 PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
1105 {
1106   PetscErrorCode ierr;
1107   PetscInt       m1,n1,m2,n2;
1108 
1109   PetscFunctionBegin;
1110   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1111   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1112   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1113   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1114   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1115   if (pc->setupcalled && Amat && Pmat) {
1116     ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr);
1117     ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr);
1118     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);
1119     ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr);
1120     ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr);
1121     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);
1122   }
1123 
1124   /* reference first in case the matrices are the same */
1125   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1126   ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1127   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1128   ierr     = MatDestroy(&pc->pmat);CHKERRQ(ierr);
1129   pc->mat  = Amat;
1130   pc->pmat = Pmat;
1131 
1132   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1133     pc->setupcalled = 1;
1134   }
1135   pc->flag = flag;
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "PCGetOperators"
1141 /*@C
1142    PCGetOperators - Gets the matrix associated with the linear system and
1143    possibly a different one associated with the preconditioner.
1144 
1145    Not collective, though parallel Mats are returned if the PC is parallel
1146 
1147    Input Parameter:
1148 .  pc - the preconditioner context
1149 
1150    Output Parameters:
1151 +  Amat - the matrix defining the linear system
1152 .  Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1153 -  flag - flag indicating information about the preconditioner
1154           matrix structure.  See PCSetOperators() for details.
1155 
1156    Level: intermediate
1157 
1158    Notes: Does not increase the reference count of the matrices, so you should not destroy them
1159 
1160    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1161       are created in PC and returned to the user. In this case, if both operators
1162       mat and pmat are requested, two DIFFERENT operators will be returned. If
1163       only one is requested both operators in the PC will be the same (i.e. as
1164       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1165       The user must set the sizes of the returned matrices and their type etc just
1166       as if the user created them with MatCreate(). For example,
1167 
1168 $         KSP/PCGetOperators(ksp/pc,&Amat,NULL,NULL); is equivalent to
1169 $           set size, type, etc of Amat
1170 
1171 $         MatCreate(comm,&mat);
1172 $         KSP/PCSetOperators(ksp/pc,Amat,Amat,SAME_NONZERO_PATTERN);
1173 $         PetscObjectDereference((PetscObject)mat);
1174 $           set size, type, etc of Amat
1175 
1176      and
1177 
1178 $         KSP/PCGetOperators(ksp/pc,&Amat,&Pmat,NULL); is equivalent to
1179 $           set size, type, etc of Amat and Pmat
1180 
1181 $         MatCreate(comm,&Amat);
1182 $         MatCreate(comm,&Pmat);
1183 $         KSP/PCSetOperators(ksp/pc,Amat,Pmat,SAME_NONZERO_PATTERN);
1184 $         PetscObjectDereference((PetscObject)Amat);
1185 $         PetscObjectDereference((PetscObject)Pmat);
1186 $           set size, type, etc of Amat and Pmat
1187 
1188     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1189     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1190     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1191     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1192     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1193     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1194     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1195     it can be created for you?
1196 
1197 
1198 .keywords: PC, get, operators, matrix, linear system
1199 
1200 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1201 @*/
1202 PetscErrorCode  PCGetOperators(PC pc,Mat *Amat,Mat *Pmat,MatStructure *flag)
1203 {
1204   PetscErrorCode ierr;
1205 
1206   PetscFunctionBegin;
1207   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1208   if (Amat) {
1209     if (!pc->mat) {
1210       if (pc->pmat && !Pmat) {  /* Apmat has been set, but user did not request it, so use for Amat */
1211         pc->mat = pc->pmat;
1212         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1213       } else {                  /* both Amat and Pmat are empty */
1214         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr);
1215         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1216           pc->pmat = pc->mat;
1217           ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1218         }
1219       }
1220     }
1221     *Amat = pc->mat;
1222   }
1223   if (Pmat) {
1224     if (!pc->pmat) {
1225       if (pc->mat && !Amat) {    /* Amat has been set but was not requested, so use for pmat */
1226         pc->pmat = pc->mat;
1227         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1228       } else {
1229         ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr);
1230         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1231           pc->mat = pc->pmat;
1232           ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1233         }
1234       }
1235     }
1236     *Pmat = pc->pmat;
1237   }
1238   if (flag) *flag = pc->flag;
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNCT__
1243 #define __FUNCT__ "PCGetOperatorsSet"
1244 /*@C
1245    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1246    possibly a different one associated with the preconditioner have been set in the PC.
1247 
1248    Not collective, though the results on all processes should be the same
1249 
1250    Input Parameter:
1251 .  pc - the preconditioner context
1252 
1253    Output Parameters:
1254 +  mat - the matrix associated with the linear system was set
1255 -  pmat - matrix associated with the preconditioner was set, usually the same
1256 
1257    Level: intermediate
1258 
1259 .keywords: PC, get, operators, matrix, linear system
1260 
1261 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1262 @*/
1263 PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1264 {
1265   PetscFunctionBegin;
1266   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1267   if (mat) *mat = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1268   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNCT__
1273 #define __FUNCT__ "PCFactorGetMatrix"
1274 /*@
1275    PCFactorGetMatrix - Gets the factored matrix from the
1276    preconditioner context.  This routine is valid only for the LU,
1277    incomplete LU, Cholesky, and incomplete Cholesky methods.
1278 
1279    Not Collective on PC though Mat is parallel if PC is parallel
1280 
1281    Input Parameters:
1282 .  pc - the preconditioner context
1283 
1284    Output parameters:
1285 .  mat - the factored matrix
1286 
1287    Level: advanced
1288 
1289    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1290 
1291 .keywords: PC, get, factored, matrix
1292 @*/
1293 PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1294 {
1295   PetscErrorCode ierr;
1296 
1297   PetscFunctionBegin;
1298   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1299   PetscValidPointer(mat,2);
1300   if (pc->ops->getfactoredmatrix) {
1301     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1302   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 #undef __FUNCT__
1307 #define __FUNCT__ "PCSetOptionsPrefix"
1308 /*@C
1309    PCSetOptionsPrefix - Sets the prefix used for searching for all
1310    PC options in the database.
1311 
1312    Logically Collective on PC
1313 
1314    Input Parameters:
1315 +  pc - the preconditioner context
1316 -  prefix - the prefix string to prepend to all PC option requests
1317 
1318    Notes:
1319    A hyphen (-) must NOT be given at the beginning of the prefix name.
1320    The first character of all runtime options is AUTOMATICALLY the
1321    hyphen.
1322 
1323    Level: advanced
1324 
1325 .keywords: PC, set, options, prefix, database
1326 
1327 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1328 @*/
1329 PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1330 {
1331   PetscErrorCode ierr;
1332 
1333   PetscFunctionBegin;
1334   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1335   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 #undef __FUNCT__
1340 #define __FUNCT__ "PCAppendOptionsPrefix"
1341 /*@C
1342    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1343    PC options in the database.
1344 
1345    Logically Collective on PC
1346 
1347    Input Parameters:
1348 +  pc - the preconditioner context
1349 -  prefix - the prefix string to prepend to all PC option requests
1350 
1351    Notes:
1352    A hyphen (-) must NOT be given at the beginning of the prefix name.
1353    The first character of all runtime options is AUTOMATICALLY the
1354    hyphen.
1355 
1356    Level: advanced
1357 
1358 .keywords: PC, append, options, prefix, database
1359 
1360 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1361 @*/
1362 PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1363 {
1364   PetscErrorCode ierr;
1365 
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1368   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 #undef __FUNCT__
1373 #define __FUNCT__ "PCGetOptionsPrefix"
1374 /*@C
1375    PCGetOptionsPrefix - Gets the prefix used for searching for all
1376    PC options in the database.
1377 
1378    Not Collective
1379 
1380    Input Parameters:
1381 .  pc - the preconditioner context
1382 
1383    Output Parameters:
1384 .  prefix - pointer to the prefix string used, is returned
1385 
1386    Notes: On the fortran side, the user should pass in a string 'prifix' of
1387    sufficient length to hold the prefix.
1388 
1389    Level: advanced
1390 
1391 .keywords: PC, get, options, prefix, database
1392 
1393 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1394 @*/
1395 PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1396 {
1397   PetscErrorCode ierr;
1398 
1399   PetscFunctionBegin;
1400   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1401   PetscValidPointer(prefix,2);
1402   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "PCPreSolve"
1408 /*@
1409    PCPreSolve - Optional pre-solve phase, intended for any
1410    preconditioner-specific actions that must be performed before
1411    the iterative solve itself.
1412 
1413    Collective on PC
1414 
1415    Input Parameters:
1416 +  pc - the preconditioner context
1417 -  ksp - the Krylov subspace context
1418 
1419    Level: developer
1420 
1421    Sample of Usage:
1422 .vb
1423     PCPreSolve(pc,ksp);
1424     KSPSolve(ksp,b,x);
1425     PCPostSolve(pc,ksp);
1426 .ve
1427 
1428    Notes:
1429    The pre-solve phase is distinct from the PCSetUp() phase.
1430 
1431    KSPSolve() calls this directly, so is rarely called by the user.
1432 
1433 .keywords: PC, pre-solve
1434 
1435 .seealso: PCPostSolve()
1436 @*/
1437 PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1438 {
1439   PetscErrorCode ierr;
1440   Vec            x,rhs;
1441 
1442   PetscFunctionBegin;
1443   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1444   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1445   pc->presolvedone++;
1446   if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice");
1447   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1448   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1449 
1450   if (pc->ops->presolve) {
1451     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1452   }
1453   PetscFunctionReturn(0);
1454 }
1455 
1456 #undef __FUNCT__
1457 #define __FUNCT__ "PCPostSolve"
1458 /*@
1459    PCPostSolve - Optional post-solve phase, intended for any
1460    preconditioner-specific actions that must be performed after
1461    the iterative solve itself.
1462 
1463    Collective on PC
1464 
1465    Input Parameters:
1466 +  pc - the preconditioner context
1467 -  ksp - the Krylov subspace context
1468 
1469    Sample of Usage:
1470 .vb
1471     PCPreSolve(pc,ksp);
1472     KSPSolve(ksp,b,x);
1473     PCPostSolve(pc,ksp);
1474 .ve
1475 
1476    Note:
1477    KSPSolve() calls this routine directly, so it is rarely called by the user.
1478 
1479    Level: developer
1480 
1481 .keywords: PC, post-solve
1482 
1483 .seealso: PCPreSolve(), KSPSolve()
1484 @*/
1485 PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1486 {
1487   PetscErrorCode ierr;
1488   Vec            x,rhs;
1489 
1490   PetscFunctionBegin;
1491   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1492   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1493   pc->presolvedone--;
1494   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1495   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1496   if (pc->ops->postsolve) {
1497     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1498   }
1499   PetscFunctionReturn(0);
1500 }
1501 
1502 #undef __FUNCT__
1503 #define __FUNCT__ "PCLoad"
1504 /*@C
1505   PCLoad - Loads a PC that has been stored in binary  with PCView().
1506 
1507   Collective on PetscViewer
1508 
1509   Input Parameters:
1510 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or
1511            some related function before a call to PCLoad().
1512 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1513 
1514    Level: intermediate
1515 
1516   Notes:
1517    The type is determined by the data in the file, any type set into the PC before this call is ignored.
1518 
1519   Notes for advanced users:
1520   Most users should not need to know the details of the binary storage
1521   format, since PCLoad() and PCView() completely hide these details.
1522   But for anyone who's interested, the standard binary matrix storage
1523   format is
1524 .vb
1525      has not yet been determined
1526 .ve
1527 
1528 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad()
1529 @*/
1530 PetscErrorCode  PCLoad(PC newdm, PetscViewer viewer)
1531 {
1532   PetscErrorCode ierr;
1533   PetscBool      isbinary;
1534   PetscInt       classid;
1535   char           type[256];
1536 
1537   PetscFunctionBegin;
1538   PetscValidHeaderSpecific(newdm,PC_CLASSID,1);
1539   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1540   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1541   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1542 
1543   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
1544   if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file");
1545   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
1546   ierr = PCSetType(newdm, type);CHKERRQ(ierr);
1547   if (newdm->ops->load) {
1548     ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1549   }
1550   PetscFunctionReturn(0);
1551 }
1552 
1553 #include <petscdraw.h>
1554 #undef __FUNCT__
1555 #define __FUNCT__ "PCView"
1556 /*@C
1557    PCView - Prints the PC data structure.
1558 
1559    Collective on PC
1560 
1561    Input Parameters:
1562 +  PC - the PC context
1563 -  viewer - optional visualization context
1564 
1565    Note:
1566    The available visualization contexts include
1567 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1568 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1569          output where only the first processor opens
1570          the file.  All other processors send their
1571          data to the first processor to print.
1572 
1573    The user can open an alternative visualization contexts with
1574    PetscViewerASCIIOpen() (output to a specified file).
1575 
1576    Level: developer
1577 
1578 .keywords: PC, view
1579 
1580 .seealso: KSPView(), PetscViewerASCIIOpen()
1581 @*/
1582 PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1583 {
1584   PCType            cstr;
1585   PetscErrorCode    ierr;
1586   PetscBool         iascii,isstring,isbinary,isdraw;
1587   PetscViewerFormat format;
1588 
1589   PetscFunctionBegin;
1590   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1591   if (!viewer) {
1592     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
1593   }
1594   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1595   PetscCheckSameComm(pc,1,viewer,2);
1596 
1597   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1598   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1599   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1600   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1601 
1602   if (iascii) {
1603     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1604     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer,"PC Object");CHKERRQ(ierr);
1605     if (pc->ops->view) {
1606       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1607       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1608       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1609     }
1610     if (pc->mat) {
1611       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1612       if (pc->pmat == pc->mat) {
1613         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1614         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1615         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1616         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1617       } else {
1618         if (pc->pmat) {
1619           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1620         } else {
1621           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1622         }
1623         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1624         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1625         if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1626         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1627       }
1628       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1629     }
1630   } else if (isstring) {
1631     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1632     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1633     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1634   } else if (isbinary) {
1635     PetscInt    classid = PC_FILE_CLASSID;
1636     MPI_Comm    comm;
1637     PetscMPIInt rank;
1638     char        type[256];
1639 
1640     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1641     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1642     if (!rank) {
1643       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1644       ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr);
1645       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
1646     }
1647     if (pc->ops->view) {
1648       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1649     }
1650   } else if (isdraw) {
1651     PetscDraw draw;
1652     char      str[25];
1653     PetscReal x,y,bottom,h;
1654     PetscInt  n;
1655 
1656     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1657     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1658     if (pc->mat) {
1659       ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr);
1660       ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr);
1661     } else {
1662       ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr);
1663     }
1664     ierr   = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1665     bottom = y - h;
1666     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1667     if (pc->ops->view) {
1668       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1669     }
1670     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1671   }
1672   PetscFunctionReturn(0);
1673 }
1674 
1675 
1676 #undef __FUNCT__
1677 #define __FUNCT__ "PCSetInitialGuessNonzero"
1678 /*@
1679    PCSetInitialGuessNonzero - Tells the iterative solver that the
1680    initial guess is nonzero; otherwise PC assumes the initial guess
1681    is to be zero (and thus zeros it out before solving).
1682 
1683    Logically Collective on PC
1684 
1685    Input Parameters:
1686 +  pc - iterative context obtained from PCCreate()
1687 -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1688 
1689    Level: Developer
1690 
1691    Notes:
1692     This is a weird function. Since PC's are linear operators on the right hand side they
1693     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1694     PCKSP, PCREDUNDANT and PCHMPI and causes the inner KSP object to use the nonzero
1695     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1696 
1697 
1698 .keywords: PC, set, initial guess, nonzero
1699 
1700 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1701 @*/
1702 PetscErrorCode  PCSetInitialGuessNonzero(PC pc,PetscBool flg)
1703 {
1704   PetscFunctionBegin;
1705   PetscValidLogicalCollectiveBool(pc,flg,2);
1706   pc->nonzero_guess = flg;
1707   PetscFunctionReturn(0);
1708 }
1709 
1710 #undef __FUNCT__
1711 #define __FUNCT__ "PCRegister"
1712 /*@C
1713   PCRegister - See PCRegisterDynamic()
1714 
1715   Level: advanced
1716 @*/
1717 PetscErrorCode  PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1718 {
1719   PetscErrorCode ierr;
1720   char           fullname[PETSC_MAX_PATH_LEN];
1721 
1722   PetscFunctionBegin;
1723   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
1724   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&PCList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1725   PetscFunctionReturn(0);
1726 }
1727 
1728 #undef __FUNCT__
1729 #define __FUNCT__ "PCComputeExplicitOperator"
1730 /*@
1731     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1732 
1733     Collective on PC
1734 
1735     Input Parameter:
1736 .   pc - the preconditioner object
1737 
1738     Output Parameter:
1739 .   mat - the explict preconditioned operator
1740 
1741     Notes:
1742     This computation is done by applying the operators to columns of the
1743     identity matrix.
1744 
1745     Currently, this routine uses a dense matrix format when 1 processor
1746     is used and a sparse format otherwise.  This routine is costly in general,
1747     and is recommended for use only with relatively small systems.
1748 
1749     Level: advanced
1750 
1751 .keywords: PC, compute, explicit, operator
1752 
1753 .seealso: KSPComputeExplicitOperator()
1754 
1755 @*/
1756 PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1757 {
1758   Vec            in,out;
1759   PetscErrorCode ierr;
1760   PetscInt       i,M,m,*rows,start,end;
1761   PetscMPIInt    size;
1762   MPI_Comm       comm;
1763   PetscScalar    *array,one = 1.0;
1764 
1765   PetscFunctionBegin;
1766   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1767   PetscValidPointer(mat,2);
1768 
1769   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1770   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1771 
1772   if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1773   ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1774   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1775   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1776   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1777   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1778   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr);
1779   for (i=0; i<m; i++) rows[i] = start + i;
1780 
1781   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
1782   ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr);
1783   if (size == 1) {
1784     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1785     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
1786   } else {
1787     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1788     ierr = MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);CHKERRQ(ierr);
1789   }
1790 
1791   for (i=0; i<M; i++) {
1792 
1793     ierr = VecSet(in,0.0);CHKERRQ(ierr);
1794     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1795     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1796     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1797 
1798     /* should fix, allowing user to choose side */
1799     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1800 
1801     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1802     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1803     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1804 
1805   }
1806   ierr = PetscFree(rows);CHKERRQ(ierr);
1807   ierr = VecDestroy(&out);CHKERRQ(ierr);
1808   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1809   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1810   PetscFunctionReturn(0);
1811 }
1812 
1813 #undef __FUNCT__
1814 #define __FUNCT__ "PCSetCoordinates"
1815 /*@
1816    PCSetCoordinates - sets the coordinates of all the nodes on the local process
1817 
1818    Collective on PC
1819 
1820    Input Parameters:
1821 +  pc - the solver context
1822 .  dim - the dimension of the coordinates 1, 2, or 3
1823 -  coords - the coordinates
1824 
1825    Level: intermediate
1826 
1827    Notes: coords is an array of the 3D coordinates for the nodes on
1828    the local processor.  So if there are 108 equation on a processor
1829    for a displacement finite element discretization of elasticity (so
1830    that there are 36 = 108/3 nodes) then the array must have 108
1831    double precision values (ie, 3 * 36).  These x y z coordinates
1832    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1833    ... , N-1.z ].
1834 
1835 .seealso: MatSetNearNullSpace
1836 @*/
1837 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1838 {
1839   PetscErrorCode ierr;
1840 
1841   PetscFunctionBegin;
1842   ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr);
1843   PetscFunctionReturn(0);
1844 }
1845