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