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