xref: /petsc/src/ksp/pc/interface/precon.c (revision fbf9dbe564678ed6eff1806adbc4c4f01b9743f4)
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 .  pc - 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 Note:
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) {
998     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since operator is unchanged\n"));
999     PetscFunctionReturn(PETSC_SUCCESS);
1000   } else {
1001     if (matnonzerostate > pc->matnonzerostate) {
1002       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
1003       pc->flag = DIFFERENT_NONZERO_PATTERN;
1004     } else {
1005       PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
1006       pc->flag = SAME_NONZERO_PATTERN;
1007     }
1008   }
1009   pc->matstate        = matstate;
1010   pc->matnonzerostate = matnonzerostate;
1011 
1012   if (!((PetscObject)pc)->type_name) {
1013     PetscCall(PCGetDefaultType_Private(pc, &def));
1014     PetscCall(PCSetType(pc, def));
1015   }
1016 
1017   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
1018   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
1019   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
1020   if (pc->ops->setup) {
1021     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
1022     PetscCall(KSPInitializePackage());
1023     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
1024     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1025     PetscUseTypeMethod(pc, setup);
1026     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
1027     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
1028   }
1029   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1030   if (!pc->setupcalled) pc->setupcalled = 1;
1031   PetscFunctionReturn(PETSC_SUCCESS);
1032 }
1033 
1034 /*@
1035    PCSetUpOnBlocks - Sets up the preconditioner for each block in
1036    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
1037    methods.
1038 
1039    Collective
1040 
1041    Input Parameter:
1042 .  pc - the preconditioner context
1043 
1044    Level: developer
1045 
1046    Note:
1047    For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1048    called on the outer `PC`, this routine ensures it is called.
1049 
1050 .seealso: `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCSetUp()`
1051 @*/
1052 PetscErrorCode PCSetUpOnBlocks(PC pc)
1053 {
1054   PetscFunctionBegin;
1055   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1056   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
1057   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1058   PetscUseTypeMethod(pc, setuponblocks);
1059   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
1060   PetscFunctionReturn(PETSC_SUCCESS);
1061 }
1062 
1063 /*@C
1064    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
1065    submatrices that arise within certain subdomain-based preconditioners.
1066    The basic submatrices are extracted from the preconditioner matrix as
1067    usual; the user can then alter these (for example, to set different boundary
1068    conditions for each submatrix) before they are used for the local solves.
1069 
1070    Logically Collective
1071 
1072    Input Parameters:
1073 +  pc - the preconditioner context
1074 .  func - routine for modifying the submatrices
1075 -  ctx - optional user-defined context (may be null)
1076 
1077    Calling sequence of `func`:
1078 $  PetscErrorCode func(PC pc, PetscInt nsub, IS *row, IS *col, Mat *submat, void *ctx);
1079 +  pc - the preconditioner context
1080 .  row - an array of index sets that contain the global row numbers
1081          that comprise each local submatrix
1082 .  col - an array of index sets that contain the global column numbers
1083          that comprise each local submatrix
1084 .  submat - array of local submatrices
1085 -  ctx - optional user-defined context for private data for the
1086          user-defined func routine (may be null)
1087 
1088    Level: advanced
1089 
1090    Notes:
1091    `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1092    `KSPSolve()`.
1093 
1094    A routine set by `PCSetModifySubMatrices()` is currently called within
1095    the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
1096    preconditioners.  All other preconditioners ignore this routine.
1097 
1098 .seealso: `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
1099 @*/
1100 PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC, PetscInt, const IS[], const IS[], Mat[], void *), void *ctx)
1101 {
1102   PetscFunctionBegin;
1103   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1104   pc->modifysubmatrices  = func;
1105   pc->modifysubmatricesP = ctx;
1106   PetscFunctionReturn(PETSC_SUCCESS);
1107 }
1108 
1109 /*@C
1110    PCModifySubMatrices - Calls an optional user-defined routine within
1111    certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
1112 
1113    Collective
1114 
1115    Input Parameters:
1116 +  pc - the preconditioner context
1117 .  nsub - the number of local submatrices
1118 .  row - an array of index sets that contain the global row numbers
1119          that comprise each local submatrix
1120 .  col - an array of index sets that contain the global column numbers
1121          that comprise each local submatrix
1122 .  submat - array of local submatrices
1123 -  ctx - optional user-defined context for private data for the
1124          user-defined routine (may be null)
1125 
1126    Output Parameter:
1127 .  submat - array of local submatrices (the entries of which may
1128             have been modified)
1129 
1130    Level: developer
1131 
1132    Notes:
1133    The user should NOT generally call this routine, as it will
1134    automatically be called within certain preconditioners (currently
1135    block Jacobi, additive Schwarz) if set.
1136 
1137    The basic submatrices are extracted from the preconditioner matrix
1138    as usual; the user can then alter these (for example, to set different
1139    boundary conditions for each submatrix) before they are used for the
1140    local solves.
1141 
1142 .seealso: `PC`, `PCSetModifySubMatrices()`
1143 @*/
1144 PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1145 {
1146   PetscFunctionBegin;
1147   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1148   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
1149   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
1150   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
1151   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
1152   PetscFunctionReturn(PETSC_SUCCESS);
1153 }
1154 
1155 /*@
1156    PCSetOperators - Sets the matrix associated with the linear system and
1157    a (possibly) different one associated with the preconditioner.
1158 
1159    Logically Collective
1160 
1161    Input Parameters:
1162 +  pc - the preconditioner context
1163 .  Amat - the matrix that defines the linear system
1164 -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
1165 
1166    Level: intermediate
1167 
1168    Notes:
1169     Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
1170 
1171     If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1172     first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1173     on it and then pass it back in in your call to `KSPSetOperators()`.
1174 
1175    More Notes about Repeated Solution of Linear Systems:
1176    PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
1177    to zero after a linear solve; the user is completely responsible for
1178    matrix assembly.  See the routine `MatZeroEntries()` if desiring to
1179    zero all elements of a matrix.
1180 
1181 .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`
1182  @*/
1183 PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1184 {
1185   PetscInt m1, n1, m2, n2;
1186 
1187   PetscFunctionBegin;
1188   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1189   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1190   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
1191   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
1192   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
1193   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
1194     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
1195     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
1196     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);
1197     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
1198     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
1199     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);
1200   }
1201 
1202   if (Pmat != pc->pmat) {
1203     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1204     pc->matnonzerostate = -1;
1205     pc->matstate        = -1;
1206   }
1207 
1208   /* reference first in case the matrices are the same */
1209   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
1210   PetscCall(MatDestroy(&pc->mat));
1211   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
1212   PetscCall(MatDestroy(&pc->pmat));
1213   pc->mat  = Amat;
1214   pc->pmat = Pmat;
1215   PetscFunctionReturn(PETSC_SUCCESS);
1216 }
1217 
1218 /*@
1219    PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
1220 
1221    Logically Collective
1222 
1223    Input Parameters:
1224 +  pc - the preconditioner context
1225 -  flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1226 
1227     Level: intermediate
1228 
1229    Note:
1230    Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1231    prevents this.
1232 
1233 .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
1234  @*/
1235 PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1236 {
1237   PetscFunctionBegin;
1238   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1239   PetscValidLogicalCollectiveBool(pc, flag, 2);
1240   pc->reusepreconditioner = flag;
1241   PetscFunctionReturn(PETSC_SUCCESS);
1242 }
1243 
1244 /*@
1245    PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1246 
1247    Not Collective
1248 
1249    Input Parameter:
1250 .  pc - the preconditioner context
1251 
1252    Output Parameter:
1253 .  flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1254 
1255    Level: intermediate
1256 
1257 .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1258  @*/
1259 PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1260 {
1261   PetscFunctionBegin;
1262   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1263   PetscValidBoolPointer(flag, 2);
1264   *flag = pc->reusepreconditioner;
1265   PetscFunctionReturn(PETSC_SUCCESS);
1266 }
1267 
1268 /*@
1269    PCGetOperators - Gets the matrix associated with the linear system and
1270    possibly a different one associated with the preconditioner.
1271 
1272    Not Collective, though parallel mats are returned if pc is parallel
1273 
1274    Input Parameter:
1275 .  pc - the preconditioner context
1276 
1277    Output Parameters:
1278 +  Amat - the matrix defining the linear system
1279 -  Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1280 
1281    Level: intermediate
1282 
1283    Note:
1284     Does not increase the reference count of the matrices, so you should not destroy them
1285 
1286    Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1287       are created in `PC` and returned to the user. In this case, if both operators
1288       mat and pmat are requested, two DIFFERENT operators will be returned. If
1289       only one is requested both operators in the PC will be the same (i.e. as
1290       if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
1291       The user must set the sizes of the returned matrices and their type etc just
1292       as if the user created them with `MatCreate()`. For example,
1293 
1294 .vb
1295          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1296            set size, type, etc of Amat
1297 
1298          MatCreate(comm,&mat);
1299          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1300          PetscObjectDereference((PetscObject)mat);
1301            set size, type, etc of Amat
1302 .ve
1303 
1304      and
1305 
1306 .vb
1307          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1308            set size, type, etc of Amat and Pmat
1309 
1310          MatCreate(comm,&Amat);
1311          MatCreate(comm,&Pmat);
1312          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1313          PetscObjectDereference((PetscObject)Amat);
1314          PetscObjectDereference((PetscObject)Pmat);
1315            set size, type, etc of Amat and Pmat
1316 .ve
1317 
1318     The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1319     of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1320     managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1321     at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1322     the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1323     you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1324     Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
1325     it can be created for you?
1326 
1327 .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
1328 @*/
1329 PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1330 {
1331   PetscFunctionBegin;
1332   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1333   if (Amat) {
1334     if (!pc->mat) {
1335       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
1336         pc->mat = pc->pmat;
1337         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1338       } else { /* both Amat and Pmat are empty */
1339         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1340         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1341           pc->pmat = pc->mat;
1342           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1343         }
1344       }
1345     }
1346     *Amat = pc->mat;
1347   }
1348   if (Pmat) {
1349     if (!pc->pmat) {
1350       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
1351         pc->pmat = pc->mat;
1352         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1353       } else {
1354         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1355         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1356           pc->mat = pc->pmat;
1357           PetscCall(PetscObjectReference((PetscObject)pc->mat));
1358         }
1359       }
1360     }
1361     *Pmat = pc->pmat;
1362   }
1363   PetscFunctionReturn(PETSC_SUCCESS);
1364 }
1365 
1366 /*@C
1367    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1368    possibly a different one associated with the preconditioner have been set in the `PC`.
1369 
1370    Not Collective, though the results on all processes should be the same
1371 
1372    Input Parameter:
1373 .  pc - the preconditioner context
1374 
1375    Output Parameters:
1376 +  mat - the matrix associated with the linear system was set
1377 -  pmat - matrix associated with the preconditioner was set, usually the same
1378 
1379    Level: intermediate
1380 
1381 .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1382 @*/
1383 PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1384 {
1385   PetscFunctionBegin;
1386   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1387   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1388   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1389   PetscFunctionReturn(PETSC_SUCCESS);
1390 }
1391 
1392 /*@
1393    PCFactorGetMatrix - Gets the factored matrix from the
1394    preconditioner context.  This routine is valid only for the `PCLU`,
1395    `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
1396 
1397    Not Collective though mat is parallel if pc is parallel
1398 
1399    Input Parameter:
1400 .  pc - the preconditioner context
1401 
1402    Output parameters:
1403 .  mat - the factored matrix
1404 
1405    Level: advanced
1406 
1407    Note:
1408     Does not increase the reference count for the matrix so DO NOT destroy it
1409 
1410 .seealso: `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
1411 @*/
1412 PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1413 {
1414   PetscFunctionBegin;
1415   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1416   PetscValidPointer(mat, 2);
1417   PetscCall(PCFactorSetUpMatSolverType(pc));
1418   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
1419   PetscFunctionReturn(PETSC_SUCCESS);
1420 }
1421 
1422 /*@C
1423    PCSetOptionsPrefix - Sets the prefix used for searching for all
1424    `PC` options in the database.
1425 
1426    Logically Collective
1427 
1428    Input Parameters:
1429 +  pc - the preconditioner context
1430 -  prefix - the prefix string to prepend to all `PC` option requests
1431 
1432    Note:
1433    A hyphen (-) must NOT be given at the beginning of the prefix name.
1434    The first character of all runtime options is AUTOMATICALLY the
1435    hyphen.
1436 
1437    Level: advanced
1438 
1439 .seealso: `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
1440 @*/
1441 PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1442 {
1443   PetscFunctionBegin;
1444   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1445   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
1446   PetscFunctionReturn(PETSC_SUCCESS);
1447 }
1448 
1449 /*@C
1450    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1451    `PC` options in the database.
1452 
1453    Logically Collective
1454 
1455    Input Parameters:
1456 +  pc - the preconditioner context
1457 -  prefix - the prefix string to prepend to all `PC` option requests
1458 
1459    Note:
1460    A hyphen (-) must NOT be given at the beginning of the prefix name.
1461    The first character of all runtime options is AUTOMATICALLY the
1462    hyphen.
1463 
1464    Level: advanced
1465 
1466 .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
1467 @*/
1468 PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1469 {
1470   PetscFunctionBegin;
1471   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1472   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
1473   PetscFunctionReturn(PETSC_SUCCESS);
1474 }
1475 
1476 /*@C
1477    PCGetOptionsPrefix - Gets the prefix used for searching for all
1478    PC options in the database.
1479 
1480    Not Collective
1481 
1482    Input Parameter:
1483 .  pc - the preconditioner context
1484 
1485    Output Parameter:
1486 .  prefix - pointer to the prefix string used, is returned
1487 
1488    Fortran Note:
1489    The user should pass in a string 'prefix' of
1490    sufficient length to hold the prefix.
1491 
1492    Level: advanced
1493 
1494 .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
1495 @*/
1496 PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1497 {
1498   PetscFunctionBegin;
1499   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1500   PetscValidPointer(prefix, 2);
1501   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
1502   PetscFunctionReturn(PETSC_SUCCESS);
1503 }
1504 
1505 /*
1506    Indicates the right hand side will be changed by KSPSolve(), this occurs for a few
1507   preconditioners including BDDC and Eisentat that transform the equations before applying
1508   the Krylov methods
1509 */
1510 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1511 {
1512   PetscFunctionBegin;
1513   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1514   PetscValidBoolPointer(change, 2);
1515   *change = PETSC_FALSE;
1516   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
1517   PetscFunctionReturn(PETSC_SUCCESS);
1518 }
1519 
1520 /*@
1521    PCPreSolve - Optional pre-solve phase, intended for any
1522    preconditioner-specific actions that must be performed before
1523    the iterative solve itself.
1524 
1525    Collective
1526 
1527    Input Parameters:
1528 +  pc - the preconditioner context
1529 -  ksp - the Krylov subspace context
1530 
1531    Level: developer
1532 
1533    Sample of Usage:
1534 .vb
1535     PCPreSolve(pc,ksp);
1536     KSPSolve(ksp,b,x);
1537     PCPostSolve(pc,ksp);
1538 .ve
1539 
1540    Notes:
1541    The pre-solve phase is distinct from the `PCSetUp()` phase.
1542 
1543    `KSPSolve()` calls this directly, so is rarely called by the user.
1544 
1545 .seealso: `PC`, `PCPostSolve()`
1546 @*/
1547 PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1548 {
1549   Vec x, rhs;
1550 
1551   PetscFunctionBegin;
1552   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1553   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1554   pc->presolvedone++;
1555   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
1556   PetscCall(KSPGetSolution(ksp, &x));
1557   PetscCall(KSPGetRhs(ksp, &rhs));
1558 
1559   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
1560   else if (pc->presolve) PetscCall((pc->presolve)(pc, ksp));
1561   PetscFunctionReturn(PETSC_SUCCESS);
1562 }
1563 
1564 /*@C
1565    PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1566    preconditioner-specific actions that must be performed before
1567    the iterative solve itself.
1568 
1569    Logically Collective
1570 
1571    Input Parameters:
1572 +   pc - the preconditioner object
1573 -   presolve - the function to call before the solve
1574 
1575    Calling sequence of `presolve`:
1576 $  PetscErrorCode presolve(PC pc, KSP ksp)
1577 +  pc - the `PC` context
1578 -  ksp - the `KSP` context
1579 
1580    Level: developer
1581 
1582 .seealso: `PC`, `PCSetUp()`, `PCPreSolve()`
1583 @*/
1584 PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC, KSP))
1585 {
1586   PetscFunctionBegin;
1587   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1588   pc->presolve = presolve;
1589   PetscFunctionReturn(PETSC_SUCCESS);
1590 }
1591 
1592 /*@
1593    PCPostSolve - Optional post-solve phase, intended for any
1594    preconditioner-specific actions that must be performed after
1595    the iterative solve itself.
1596 
1597    Collective
1598 
1599    Input Parameters:
1600 +  pc - the preconditioner context
1601 -  ksp - the Krylov subspace context
1602 
1603    Sample of Usage:
1604 .vb
1605     PCPreSolve(pc,ksp);
1606     KSPSolve(ksp,b,x);
1607     PCPostSolve(pc,ksp);
1608 .ve
1609 
1610    Note:
1611    `KSPSolve()` calls this routine directly, so it is rarely called by the user.
1612 
1613    Level: developer
1614 
1615 .seealso: `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
1616 @*/
1617 PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1618 {
1619   Vec x, rhs;
1620 
1621   PetscFunctionBegin;
1622   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1623   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1624   pc->presolvedone--;
1625   PetscCall(KSPGetSolution(ksp, &x));
1626   PetscCall(KSPGetRhs(ksp, &rhs));
1627   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
1628   PetscFunctionReturn(PETSC_SUCCESS);
1629 }
1630 
1631 /*@C
1632   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
1633 
1634   Collective
1635 
1636   Input Parameters:
1637 + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1638            some related function before a call to `PCLoad()`.
1639 - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
1640 
1641    Level: intermediate
1642 
1643   Note:
1644    The type is determined by the data in the file, any type set into the PC before this call is ignored.
1645 
1646 .seealso: `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
1647 @*/
1648 PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1649 {
1650   PetscBool isbinary;
1651   PetscInt  classid;
1652   char      type[256];
1653 
1654   PetscFunctionBegin;
1655   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
1656   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1657   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1658   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1659 
1660   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
1661   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
1662   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
1663   PetscCall(PCSetType(newdm, type));
1664   PetscTryTypeMethod(newdm, load, viewer);
1665   PetscFunctionReturn(PETSC_SUCCESS);
1666 }
1667 
1668 #include <petscdraw.h>
1669 #if defined(PETSC_HAVE_SAWS)
1670   #include <petscviewersaws.h>
1671 #endif
1672 
1673 /*@C
1674    PCViewFromOptions - View from the `PC` based on options in the database
1675 
1676    Collective
1677 
1678    Input Parameters:
1679 +  A - the `PC` context
1680 .  obj - Optional object that provides the options prefix
1681 -  name - command line option
1682 
1683    Level: intermediate
1684 
1685 .seealso: `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1686 @*/
1687 PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1688 {
1689   PetscFunctionBegin;
1690   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
1691   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1692   PetscFunctionReturn(PETSC_SUCCESS);
1693 }
1694 
1695 /*@C
1696    PCView - Prints information about the `PC`
1697 
1698    Collective
1699 
1700    Input Parameters:
1701 +  PC - the `PC` context
1702 -  viewer - optional visualization context
1703 
1704    Level: developer
1705 
1706    Notes:
1707    The available visualization contexts include
1708 +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1709 -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
1710          output where only the first processor opens
1711          the file.  All other processors send their
1712          data to the first processor to print.
1713 
1714    The user can open an alternative visualization contexts with
1715    `PetscViewerASCIIOpen()` (output to a specified file).
1716 
1717 .seealso: `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
1718 @*/
1719 PetscErrorCode PCView(PC pc, PetscViewer viewer)
1720 {
1721   PCType    cstr;
1722   PetscBool iascii, isstring, isbinary, isdraw;
1723 #if defined(PETSC_HAVE_SAWS)
1724   PetscBool issaws;
1725 #endif
1726 
1727   PetscFunctionBegin;
1728   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1729   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
1730   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1731   PetscCheckSameComm(pc, 1, viewer, 2);
1732 
1733   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1734   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1735   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1736   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1737 #if defined(PETSC_HAVE_SAWS)
1738   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
1739 #endif
1740 
1741   if (iascii) {
1742     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
1743     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
1744     PetscCall(PetscViewerASCIIPushTab(viewer));
1745     PetscTryTypeMethod(pc, view, viewer);
1746     PetscCall(PetscViewerASCIIPopTab(viewer));
1747     if (pc->mat) {
1748       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1749       if (pc->pmat == pc->mat) {
1750         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
1751         PetscCall(PetscViewerASCIIPushTab(viewer));
1752         PetscCall(MatView(pc->mat, viewer));
1753         PetscCall(PetscViewerASCIIPopTab(viewer));
1754       } else {
1755         if (pc->pmat) {
1756           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
1757         } else {
1758           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
1759         }
1760         PetscCall(PetscViewerASCIIPushTab(viewer));
1761         PetscCall(MatView(pc->mat, viewer));
1762         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
1763         PetscCall(PetscViewerASCIIPopTab(viewer));
1764       }
1765       PetscCall(PetscViewerPopFormat(viewer));
1766     }
1767   } else if (isstring) {
1768     PetscCall(PCGetType(pc, &cstr));
1769     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1770     PetscTryTypeMethod(pc, view, viewer);
1771     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1772     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1773   } else if (isbinary) {
1774     PetscInt    classid = PC_FILE_CLASSID;
1775     MPI_Comm    comm;
1776     PetscMPIInt rank;
1777     char        type[256];
1778 
1779     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1780     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1781     if (rank == 0) {
1782       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
1783       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
1784       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
1785     }
1786     PetscTryTypeMethod(pc, view, viewer);
1787   } else if (isdraw) {
1788     PetscDraw draw;
1789     char      str[25];
1790     PetscReal x, y, bottom, h;
1791     PetscInt  n;
1792 
1793     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1794     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
1795     if (pc->mat) {
1796       PetscCall(MatGetSize(pc->mat, &n, NULL));
1797       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
1798     } else {
1799       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
1800     }
1801     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
1802     bottom = y - h;
1803     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1804     PetscTryTypeMethod(pc, view, viewer);
1805     PetscCall(PetscDrawPopCurrentPoint(draw));
1806 #if defined(PETSC_HAVE_SAWS)
1807   } else if (issaws) {
1808     PetscMPIInt rank;
1809 
1810     PetscCall(PetscObjectName((PetscObject)pc));
1811     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1812     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
1813     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1814     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1815 #endif
1816   }
1817   PetscFunctionReturn(PETSC_SUCCESS);
1818 }
1819 
1820 /*@C
1821   PCRegister -  Adds a method to the preconditioner package.
1822 
1823    Not collective
1824 
1825    Input Parameters:
1826 +  sname - name of a new user-defined solver
1827 -  function - routine to create method context
1828 
1829    Sample usage:
1830 .vb
1831    PCRegister("my_solver", MySolverCreate);
1832 .ve
1833 
1834    Then, your solver can be chosen with the procedural interface via
1835 $     PCSetType(pc, "my_solver")
1836    or at runtime via the option
1837 $     -pc_type my_solver
1838 
1839    Level: advanced
1840 
1841    Note:
1842    `PCRegister()` may be called multiple times to add several user-defined preconditioners.
1843 
1844 .seealso: `PC`, `PCType`, `PCRegisterAll()`
1845 @*/
1846 PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1847 {
1848   PetscFunctionBegin;
1849   PetscCall(PCInitializePackage());
1850   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
1851   PetscFunctionReturn(PETSC_SUCCESS);
1852 }
1853 
1854 static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1855 {
1856   PC pc;
1857 
1858   PetscFunctionBegin;
1859   PetscCall(MatShellGetContext(A, &pc));
1860   PetscCall(PCApply(pc, X, Y));
1861   PetscFunctionReturn(PETSC_SUCCESS);
1862 }
1863 
1864 /*@
1865     PCComputeOperator - Computes the explicit preconditioned operator.
1866 
1867     Collective
1868 
1869     Input Parameters:
1870 +   pc - the preconditioner object
1871 -   mattype - the matrix type to be used for the operator
1872 
1873     Output Parameter:
1874 .   mat - the explicit preconditioned operator
1875 
1876     Level: advanced
1877 
1878     Note:
1879     This computation is done by applying the operators to columns of the identity matrix.
1880     This routine is costly in general, and is recommended for use only with relatively small systems.
1881     Currently, this routine uses a dense matrix format when mattype == NULL
1882 
1883 .seealso: `PC`, `KSPComputeOperator()`, `MatType`
1884 @*/
1885 PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1886 {
1887   PetscInt N, M, m, n;
1888   Mat      A, Apc;
1889 
1890   PetscFunctionBegin;
1891   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1892   PetscValidPointer(mat, 3);
1893   PetscCall(PCGetOperators(pc, &A, NULL));
1894   PetscCall(MatGetLocalSize(A, &m, &n));
1895   PetscCall(MatGetSize(A, &M, &N));
1896   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
1897   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
1898   PetscCall(MatComputeOperator(Apc, mattype, mat));
1899   PetscCall(MatDestroy(&Apc));
1900   PetscFunctionReturn(PETSC_SUCCESS);
1901 }
1902 
1903 /*@
1904    PCSetCoordinates - sets the coordinates of all the nodes on the local process
1905 
1906    Collective
1907 
1908    Input Parameters:
1909 +  pc - the solver context
1910 .  dim - the dimension of the coordinates 1, 2, or 3
1911 .  nloc - the blocked size of the coordinates array
1912 -  coords - the coordinates array
1913 
1914    Level: intermediate
1915 
1916    Note:
1917    `coords` is an array of the dim coordinates for the nodes on
1918    the local processor, of size `dim`*`nloc`.
1919    If there are 108 equation on a processor
1920    for a displacement finite element discretization of elasticity (so
1921    that there are nloc = 36 = 108/3 nodes) then the array must have 108
1922    double precision values (ie, 3 * 36).  These x y z coordinates
1923    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1924    ... , N-1.z ].
1925 
1926 .seealso: `PC`, `MatSetNearNullSpace()`
1927 @*/
1928 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1929 {
1930   PetscFunctionBegin;
1931   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1932   PetscValidLogicalCollectiveInt(pc, dim, 2);
1933   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
1934   PetscFunctionReturn(PETSC_SUCCESS);
1935 }
1936 
1937 /*@
1938    PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1939 
1940    Logically Collective
1941 
1942    Input Parameter:
1943 .  pc - the precondition context
1944 
1945    Output Parameters:
1946 +  num_levels - the number of levels
1947 -  interpolations - the interpolation matrices (size of num_levels-1)
1948 
1949    Level: advanced
1950 
1951    Developer Note:
1952    Why is this here instead of in `PCMG` etc?
1953 
1954 .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
1955 @*/
1956 PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
1957 {
1958   PetscFunctionBegin;
1959   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1960   PetscValidIntPointer(num_levels, 2);
1961   PetscValidPointer(interpolations, 3);
1962   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
1963   PetscFunctionReturn(PETSC_SUCCESS);
1964 }
1965 
1966 /*@
1967    PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
1968 
1969    Logically Collective
1970 
1971    Input Parameter:
1972 .  pc - the precondition context
1973 
1974    Output Parameters:
1975 +  num_levels - the number of levels
1976 -  coarseOperators - the coarse operator matrices (size of num_levels-1)
1977 
1978    Level: advanced
1979 
1980    Developer Note:
1981    Why is this here instead of in `PCMG` etc?
1982 
1983 .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
1984 @*/
1985 PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
1986 {
1987   PetscFunctionBegin;
1988   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1989   PetscValidIntPointer(num_levels, 2);
1990   PetscValidPointer(coarseOperators, 3);
1991   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
1992   PetscFunctionReturn(PETSC_SUCCESS);
1993 }
1994