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