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