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