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