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