xref: /petsc/include/petsc/private/kspimpl.h (revision d68d54c7c77e67f355f59cef0b2b3547c8e78a09)
1 
2 #ifndef _KSPIMPL_H
3   #define _KSPIMPL_H
4 
5   #include <petscksp.h>
6   #include <petscds.h>
7   #include <petsc/private/petscimpl.h>
8 
9 /* SUBMANSEC = KSP */
10 
11 PETSC_EXTERN PetscBool      KSPRegisterAllCalled;
12 PETSC_EXTERN PetscBool      KSPMonitorRegisterAllCalled;
13 PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
14 PETSC_EXTERN PetscErrorCode KSPMonitorRegisterAll(void);
15 PETSC_EXTERN PetscErrorCode KSPGuessRegisterAll(void);
16 PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);
17 
18 typedef struct _KSPOps *KSPOps;
19 
20 struct _KSPOps {
21   PetscErrorCode (*buildsolution)(KSP, Vec, Vec *);      /* Returns a pointer to the solution, or
22                                                           calculates the solution in a
23                                                           user-provided area. */
24   PetscErrorCode (*buildresidual)(KSP, Vec, Vec, Vec *); /* Returns a pointer to the residual, or
25                                                           calculates the residual in a
26                                                           user-provided area.  */
27   PetscErrorCode (*solve)(KSP);                          /* actual solver */
28   PetscErrorCode (*matsolve)(KSP, Mat, Mat);             /* multiple dense RHS solver */
29   PetscErrorCode (*setup)(KSP);
30   PetscErrorCode (*setfromoptions)(KSP, PetscOptionItems *);
31   PetscErrorCode (*publishoptions)(KSP);
32   PetscErrorCode (*computeextremesingularvalues)(KSP, PetscReal *, PetscReal *);
33   PetscErrorCode (*computeeigenvalues)(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);
34   PetscErrorCode (*computeritz)(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal *, PetscReal *);
35   PetscErrorCode (*destroy)(KSP);
36   PetscErrorCode (*view)(KSP, PetscViewer);
37   PetscErrorCode (*reset)(KSP);
38   PetscErrorCode (*load)(KSP, PetscViewer);
39 };
40 
41 typedef struct _KSPGuessOps *KSPGuessOps;
42 
43 struct _KSPGuessOps {
44   PetscErrorCode (*formguess)(KSPGuess, Vec, Vec); /* Form initial guess */
45   PetscErrorCode (*update)(KSPGuess, Vec, Vec);    /* Update database */
46   PetscErrorCode (*setfromoptions)(KSPGuess);
47   PetscErrorCode (*settolerance)(KSPGuess, PetscReal);
48   PetscErrorCode (*setup)(KSPGuess);
49   PetscErrorCode (*destroy)(KSPGuess);
50   PetscErrorCode (*view)(KSPGuess, PetscViewer);
51   PetscErrorCode (*reset)(KSPGuess);
52 };
53 
54 /*
55    Defines the KSPGuess data structure.
56 */
57 struct _p_KSPGuess {
58   PETSCHEADER(struct _KSPGuessOps);
59   KSP              ksp;       /* the parent KSP */
60   Mat              A;         /* the current linear operator */
61   PetscObjectState omatstate; /* previous linear operator state */
62   void            *data;      /* pointer to the specific implementation */
63 };
64 
65 PETSC_EXTERN PetscErrorCode KSPGuessCreate_Fischer(KSPGuess);
66 PETSC_EXTERN PetscErrorCode KSPGuessCreate_POD(KSPGuess);
67 
68   /*
69      Maximum number of monitors you can run with a single KSP
70 */
71   #define MAXKSPMONITORS    5
72   #define MAXKSPREASONVIEWS 5
73 typedef enum {
74   KSP_SETUP_NEW = 0,
75   KSP_SETUP_NEWMATRIX,
76   KSP_SETUP_NEWRHS
77 } KSPSetUpStage;
78 
79 /*
80    Defines the KSP data structure.
81 */
82 struct _p_KSP {
83   PETSCHEADER(struct _KSPOps);
84   DM        dm;
85   PetscBool dmAuto;   /* DM was created automatically by KSP */
86   PetscBool dmActive; /* KSP should use DM for computing operators */
87   /*------------------------- User parameters--------------------------*/
88   PetscInt  max_it; /* maximum number of iterations */
89   PetscInt  min_it; /* minimum number of iterations */
90   KSPGuess  guess;
91   PetscBool guess_zero,                                  /* flag for whether initial guess is 0 */
92     guess_not_read,                                      /* guess is not read, does not need to be zeroed */
93     calc_sings,                                          /* calculate extreme Singular Values */
94     calc_ritz,                                           /* calculate (harmonic) Ritz pairs */
95     guess_knoll;                                         /* use initial guess of PCApply(ksp->B,b */
96   PCSide    pc_side;                                     /* flag for left, right, or symmetric preconditioning */
97   PetscInt  normsupporttable[KSP_NORM_MAX][PC_SIDE_MAX]; /* Table of supported norms and pc_side, see KSPSetSupportedNorm() */
98   PetscReal rtol,                                        /* relative tolerance */
99     abstol,                                              /* absolute tolerance */
100     ttol,                                                /* (not set by user)  */
101     divtol;                                              /* divergence tolerance */
102   PetscReal          rnorm0;                             /* initial residual norm (used for divergence testing) */
103   PetscReal          rnorm;                              /* current residual norm */
104   KSPConvergedReason reason;
105   PetscBool          errorifnotconverged; /* create an error if the KSPSolve() does not converge */
106 
107   Vec        vec_sol, vec_rhs; /* pointer to where user has stashed
108                                       the solution and rhs, these are
109                                       never touched by the code, only
110                                       passed back to the user */
111   PetscReal *res_hist;         /* If !0 stores residual each at iteration */
112   PetscReal *res_hist_alloc;   /* If !0 means user did not provide buffer, needs deallocation */
113   size_t     res_hist_len;     /* current size of residual history array */
114   size_t     res_hist_max;     /* actual amount of storage in residual history */
115   PetscBool  res_hist_reset;   /* reset history to length zero for each new solve */
116   PetscReal *err_hist;         /* If !0 stores error at each iteration */
117   PetscReal *err_hist_alloc;   /* If !0 means user did not provide buffer, needs deallocation */
118   size_t     err_hist_len;     /* current size of error history array */
119   size_t     err_hist_max;     /* actual amount of storage in error history */
120   PetscBool  err_hist_reset;   /* reset history to length zero for each new solve */
121 
122   PetscInt  chknorm; /* only compute/check norm if iterations is great than this */
123   PetscBool lagnorm; /* Lag the residual norm calculation so that it is computed as part of the
124                                         MPI_Allreduce() for computing the inner products for the next iteration. */
125 
126   PetscInt nmax; /* maximum number of right-hand sides to be handled simultaneously */
127 
128   /* --------User (or default) routines (most return -1 on error) --------*/
129   PetscErrorCode (*monitor[MAXKSPMONITORS])(KSP, PetscInt, PetscReal, void *); /* returns control to user after */
130   PetscErrorCode (*monitordestroy[MAXKSPMONITORS])(void **);                   /* */
131   void     *monitorcontext[MAXKSPMONITORS];                                    /* residual calculation, allows user */
132   PetscInt  numbermonitors;                                                    /* to, for instance, print residual norm, etc. */
133   PetscBool pauseFinal;                                                        /* Pause all drawing monitor at the final iterate */
134 
135   PetscErrorCode (*reasonview[MAXKSPREASONVIEWS])(KSP, void *);    /* KSP converged reason view */
136   PetscErrorCode (*reasonviewdestroy[MAXKSPREASONVIEWS])(void **); /* Optional destroy routine */
137   void    *reasonviewcontext[MAXKSPREASONVIEWS];                   /* User context */
138   PetscInt numberreasonviews;                                      /* Number if reason viewers */
139 
140   PetscErrorCode (*converged)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
141   PetscErrorCode (*convergeddestroy)(void *);
142   void *cnvP;
143 
144   void *user; /* optional user-defined context */
145 
146   PC pc;
147 
148   void *data; /* holder for misc stuff associated
149                                    with a particular iterative solver */
150 
151   PetscBool         view, viewPre, viewRate, viewMat, viewPMat, viewRhs, viewSol, viewMatExp, viewEV, viewSV, viewEVExp, viewFinalRes, viewPOpExp, viewDScale;
152   PetscViewer       viewer, viewerPre, viewerRate, viewerMat, viewerPMat, viewerRhs, viewerSol, viewerMatExp, viewerEV, viewerSV, viewerEVExp, viewerFinalRes, viewerPOpExp, viewerDScale;
153   PetscViewerFormat format, formatPre, formatRate, formatMat, formatPMat, formatRhs, formatSol, formatMatExp, formatEV, formatSV, formatEVExp, formatFinalRes, formatPOpExp, formatDScale;
154 
155   /* ----------------Default work-area management -------------------- */
156   PetscInt nwork;
157   Vec     *work;
158 
159   KSPSetUpStage setupstage;
160   PetscBool     setupnewmatrix; /* true if we need to call ksp->ops->setup with KSP_SETUP_NEWMATRIX */
161 
162   PetscInt its;      /* number of iterations so far computed in THIS linear solve*/
163   PetscInt totalits; /* number of iterations used by this KSP object since it was created */
164 
165   PetscBool transpose_solve; /* solve transpose system instead */
166   struct {
167     Mat       AT, BT;
168     PetscBool use_explicittranspose; /* transpose the system explicitly in KSPSolveTranspose */
169     PetscBool reuse_transpose;       /* reuse the previous transposed system */
170   } transpose;
171 
172   KSPNormType normtype; /* type of norm used for convergence tests */
173 
174   PCSide      pc_side_set;  /* PC type set explicitly by user */
175   KSPNormType normtype_set; /* Norm type set explicitly by user */
176 
177   /*   Allow diagonally scaling the matrix before computing the preconditioner or using
178        the Krylov method. Note this is NOT just Jacobi preconditioning */
179 
180   PetscBool dscale;     /* diagonal scale system; used with KSPSetDiagonalScale() */
181   PetscBool dscalefix;  /* unscale system after solve */
182   PetscBool dscalefix2; /* system has been unscaled */
183   Vec       diagonal;   /* 1/sqrt(diag of matrix) */
184   Vec       truediagonal;
185 
186   /* Allow declaring convergence when negative curvature is detected */
187   PetscBool converged_neg_curve;
188 
189   PetscInt  setfromoptionscalled;
190   PetscBool skippcsetfromoptions; /* if set then KSPSetFromOptions() does not call PCSetFromOptions() */
191 
192   PetscErrorCode (*presolve)(KSP, Vec, Vec, void *);
193   PetscErrorCode (*postsolve)(KSP, Vec, Vec, void *);
194   void *prectx, *postctx;
195 
196   PetscInt nestlevel; /* how many levels of nesting does the KSP have */
197 };
198 
199 typedef struct { /* dummy data structure used in KSPMonitorDynamicTolerance() */
200   PetscReal coef;
201   PetscReal bnrm;
202 } KSPDynTolCtx;
203 
204 typedef struct {
205   PetscBool initialrtol;    /* default relative residual decrease is computed from initial residual, not rhs */
206   PetscBool mininitialrtol; /* default relative residual decrease is computed from min of initial residual and rhs */
207   PetscBool convmaxits;     /* if true, the convergence test returns KSP_CONVERGED_ITS if the maximum number of iterations is reached */
208   Vec       work;
209 } KSPConvergedDefaultCtx;
210 
211 static inline PetscErrorCode KSPLogResidualHistory(KSP ksp, PetscReal norm)
212 {
213   PetscFunctionBegin;
214   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
215   if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len) ksp->res_hist[ksp->res_hist_len++] = norm;
216   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
217   PetscFunctionReturn(PETSC_SUCCESS);
218 }
219 
220 static inline PetscErrorCode KSPLogErrorHistory(KSP ksp)
221 {
222   DM dm;
223 
224   PetscFunctionBegin;
225   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
226   PetscCall(KSPGetDM(ksp, &dm));
227   if (dm && ksp->err_hist && ksp->err_hist_max > ksp->err_hist_len) {
228     PetscSimplePointFunc exactSol;
229     void                *exactCtx;
230     PetscDS              ds;
231     Vec                  u;
232     PetscReal            error;
233     PetscInt             Nf;
234 
235     PetscCall(KSPBuildSolution(ksp, NULL, &u));
236     /* TODO Was needed to correct for Newton solution, but I just need to set a solution */
237     //PetscCall(VecScale(u, -1.0));
238     /* TODO Case when I have a solution */
239     if (0) {
240       PetscCall(DMGetDS(dm, &ds));
241       PetscCall(PetscDSGetNumFields(ds, &Nf));
242       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle number of fields %" PetscInt_FMT " > 1 right now", Nf);
243       PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol, &exactCtx));
244       PetscCall(DMComputeL2FieldDiff(dm, 0.0, &exactSol, &exactCtx, u, &error));
245     } else {
246       /* The null solution A 0 = 0 */
247       PetscCall(VecNorm(u, NORM_2, &error));
248     }
249     ksp->err_hist[ksp->err_hist_len++] = error;
250   }
251   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
252   PetscFunctionReturn(PETSC_SUCCESS);
253 }
254 
255 static inline PetscScalar KSPNoisyHash_Private(PetscInt xx)
256 {
257   unsigned int x = (unsigned int)xx;
258   x              = ((x >> 16) ^ x) * 0x45d9f3b;
259   x              = ((x >> 16) ^ x) * 0x45d9f3b;
260   x              = ((x >> 16) ^ x);
261   return (PetscScalar)((PetscInt64)x - 2147483648) * 5.e-10; /* center around zero, scaled about -1. to 1.*/
262 }
263 
264 static inline PetscErrorCode KSPSetNoisy_Private(Vec v)
265 {
266   PetscScalar *a;
267   PetscInt     n, istart;
268 
269   PetscFunctionBegin;
270   PetscCall(VecGetOwnershipRange(v, &istart, NULL));
271   PetscCall(VecGetLocalSize(v, &n));
272   PetscCall(VecGetArrayWrite(v, &a));
273   for (PetscInt i = 0; i < n; ++i) a[i] = KSPNoisyHash_Private(i + istart);
274   PetscCall(VecRestoreArrayWrite(v, &a));
275   PetscFunctionReturn(PETSC_SUCCESS);
276 }
277 
278 PETSC_INTERN PetscErrorCode KSPSetUpNorms_Private(KSP, PetscBool, KSPNormType *, PCSide *);
279 
280 PETSC_INTERN PetscErrorCode KSPPlotEigenContours_Private(KSP, PetscInt, const PetscReal *, const PetscReal *);
281 
282 typedef struct _p_DMKSP  *DMKSP;
283 typedef struct _DMKSPOps *DMKSPOps;
284 struct _DMKSPOps {
285   PetscErrorCode (*computeoperators)(KSP, Mat, Mat, void *);
286   PetscErrorCode (*computerhs)(KSP, Vec, void *);
287   PetscErrorCode (*computeinitialguess)(KSP, Vec, void *);
288   PetscErrorCode (*destroy)(DMKSP *);
289   PetscErrorCode (*duplicate)(DMKSP, DMKSP);
290 };
291 
292 struct _p_DMKSP {
293   PETSCHEADER(struct _DMKSPOps);
294   void *operatorsctx;
295   void *rhsctx;
296   void *initialguessctx;
297   void *data;
298 
299   /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
300    * copy-on-write. When DMGetDMKSPWrite() sees a request using a different DM, it makes a copy. Thus, if a user
301    * only interacts directly with one level, e.g., using KSPSetComputeOperators(), then coarse levels are constructed by
302    * PCMG, then the user changes the routine with another call to KSPSetComputeOperators(), it automatically propagates
303    * to all the levels. If instead, they get out a specific level and set the routines on that level, subsequent changes
304    * to the original level will no longer propagate to that level.
305    */
306   DM originaldm;
307 
308   void (*fortran_func_pointers[3])(void); /* Store our own function pointers so they are associated with the DMKSP instead of the DM */
309 };
310 PETSC_EXTERN PetscErrorCode DMGetDMKSP(DM, DMKSP *);
311 PETSC_EXTERN PetscErrorCode DMGetDMKSPWrite(DM, DMKSP *);
312 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);
313 
314 /*
315        These allow the various Krylov methods to apply to either the linear system or its transpose.
316 */
317 static inline PetscErrorCode KSP_RemoveNullSpace(KSP ksp, Vec y)
318 {
319   PetscFunctionBegin;
320   if (ksp->pc_side == PC_LEFT) {
321     Mat          A;
322     MatNullSpace nullsp;
323 
324     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
325     PetscCall(MatGetNullSpace(A, &nullsp));
326     if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, y));
327   }
328   PetscFunctionReturn(PETSC_SUCCESS);
329 }
330 
331 static inline PetscErrorCode KSP_RemoveNullSpaceTranspose(KSP ksp, Vec y)
332 {
333   PetscFunctionBegin;
334   if (ksp->pc_side == PC_LEFT) {
335     Mat          A;
336     MatNullSpace nullsp;
337 
338     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
339     PetscCall(MatGetTransposeNullSpace(A, &nullsp));
340     if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, y));
341   }
342   PetscFunctionReturn(PETSC_SUCCESS);
343 }
344 
345 static inline PetscErrorCode KSP_MatMult(KSP ksp, Mat A, Vec x, Vec y)
346 {
347   PetscFunctionBegin;
348   if (ksp->transpose_solve) PetscCall(MatMultTranspose(A, x, y));
349   else PetscCall(MatMult(A, x, y));
350   PetscFunctionReturn(PETSC_SUCCESS);
351 }
352 
353 static inline PetscErrorCode KSP_MatMultTranspose(KSP ksp, Mat A, Vec x, Vec y)
354 {
355   PetscFunctionBegin;
356   if (ksp->transpose_solve) PetscCall(MatMult(A, x, y));
357   else PetscCall(MatMultTranspose(A, x, y));
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 static inline PetscErrorCode KSP_MatMultHermitianTranspose(KSP ksp, Mat A, Vec x, Vec y)
362 {
363   PetscFunctionBegin;
364   if (!ksp->transpose_solve) PetscCall(MatMultHermitianTranspose(A, x, y));
365   else {
366     Vec w;
367 
368     PetscCall(VecDuplicate(x, &w));
369     PetscCall(VecCopy(x, w));
370     PetscCall(VecConjugate(w));
371     PetscCall(MatMult(A, w, y));
372     PetscCall(VecDestroy(&w));
373     PetscCall(VecConjugate(y));
374   }
375   PetscFunctionReturn(PETSC_SUCCESS);
376 }
377 
378 static inline PetscErrorCode KSP_PCApply(KSP ksp, Vec x, Vec y)
379 {
380   PetscFunctionBegin;
381   if (ksp->transpose_solve) {
382     PetscCall(PCApplyTranspose(ksp->pc, x, y));
383     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
384   } else {
385     PetscCall(PCApply(ksp->pc, x, y));
386     PetscCall(KSP_RemoveNullSpace(ksp, y));
387   }
388   PetscFunctionReturn(PETSC_SUCCESS);
389 }
390 
391 static inline PetscErrorCode KSP_PCApplyTranspose(KSP ksp, Vec x, Vec y)
392 {
393   PetscFunctionBegin;
394   if (ksp->transpose_solve) {
395     PetscCall(PCApply(ksp->pc, x, y));
396     PetscCall(KSP_RemoveNullSpace(ksp, y));
397   } else {
398     PetscCall(PCApplyTranspose(ksp->pc, x, y));
399     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
400   }
401   PetscFunctionReturn(PETSC_SUCCESS);
402 }
403 
404 static inline PetscErrorCode KSP_PCApplyHermitianTranspose(KSP ksp, Vec x, Vec y)
405 {
406   PetscFunctionBegin;
407   PetscCall(VecConjugate(x));
408   PetscCall(KSP_PCApplyTranspose(ksp, x, y));
409   PetscCall(VecConjugate(x));
410   PetscCall(VecConjugate(y));
411   PetscFunctionReturn(PETSC_SUCCESS);
412 }
413 
414 static inline PetscErrorCode KSP_PCMatApply(KSP ksp, Mat X, Mat Y)
415 {
416   PetscFunctionBegin;
417   if (ksp->transpose_solve) {
418     PetscBool flg;
419     PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp->pc, &flg, PCNONE, PCICC, PCCHOLESKY, ""));
420     PetscCheck(flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "PCMatApplyTranspose() not yet implemented for nonsymmetric PC");
421   }
422   PetscCall(PCMatApply(ksp->pc, X, Y));
423   PetscFunctionReturn(PETSC_SUCCESS);
424 }
425 
426 static inline PetscErrorCode KSP_PCMatApplyTranspose(KSP ksp, Mat X, Mat Y)
427 {
428   PetscFunctionBegin;
429   if (!ksp->transpose_solve) {
430     PetscBool flg;
431     PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp->pc, &flg, PCNONE, PCICC, PCCHOLESKY, ""));
432     PetscCheck(flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "PCMatApplyTranspose() not yet implemented for nonsymmetric PC");
433   }
434   PetscCall(PCMatApply(ksp->pc, X, Y));
435   PetscFunctionReturn(PETSC_SUCCESS);
436 }
437 
438 static inline PetscErrorCode KSP_PCApplyBAorAB(KSP ksp, Vec x, Vec y, Vec w)
439 {
440   PetscFunctionBegin;
441   if (ksp->transpose_solve) {
442     PetscCall(PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w));
443     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
444   } else {
445     PetscCall(PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w));
446     PetscCall(KSP_RemoveNullSpace(ksp, y));
447   }
448   PetscFunctionReturn(PETSC_SUCCESS);
449 }
450 
451 static inline PetscErrorCode KSP_PCApplyBAorABTranspose(KSP ksp, Vec x, Vec y, Vec w)
452 {
453   PetscFunctionBegin;
454   if (ksp->transpose_solve) PetscCall(PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w));
455   else PetscCall(PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w));
456   PetscFunctionReturn(PETSC_SUCCESS);
457 }
458 
459 PETSC_EXTERN PetscLogEvent KSP_GMRESOrthogonalization;
460 PETSC_EXTERN PetscLogEvent KSP_SetUp;
461 PETSC_EXTERN PetscLogEvent KSP_Solve;
462 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_0;
463 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_1;
464 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_2;
465 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_3;
466 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_4;
467 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_S;
468 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_L;
469 PETSC_EXTERN PetscLogEvent KSP_Solve_FS_U;
470 PETSC_EXTERN PetscLogEvent KSP_SolveTranspose;
471 PETSC_EXTERN PetscLogEvent KSP_MatSolve;
472 PETSC_EXTERN PetscLogEvent KSP_MatSolveTranspose;
473 
474 PETSC_INTERN PetscErrorCode MatGetSchurComplement_Basic(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
475 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *);
476 
477   /*MC
478    KSPCheckDot - Checks if the result of a dot product used by the corresponding `KSP` contains Inf or NaN. These indicate that the previous
479       application of the preconditioner generated an error. Sets a `KSPConvergedReason` and returns if the `PC` set a `PCFailedReason`.
480 
481    Collective
482 
483    Input Parameter:
484 .  ksp - the linear solver `KSP` context.
485 
486    Output Parameter:
487 .  beta - the result of the inner product
488 
489    Level: developer
490 
491    Developer Notes:
492    Used to manage returning from `KSP` solvers whose preconditioners have failed, possibly only a subset of MPI ranks, in some way
493 
494    It uses the fact that `KSP` piggy-backs the collectivity of certain error conditions on the results of norms and inner products.
495 
496 .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckNorm()`, `KSPCheckSolve()`
497 M*/
498   #define KSPCheckDot(ksp, beta) \
499     do { \
500       if (PetscIsInfOrNanScalar(beta)) { \
501         PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged due to Nan or Inf inner product"); \
502         { \
503           PCFailedReason pcreason; \
504           PetscCall(PCReduceFailedReason(ksp->pc)); \
505           PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason)); \
506           if (pcreason) { \
507             ksp->reason = KSP_DIVERGED_PC_FAILED; \
508             PetscCall(VecSetInf(ksp->vec_sol)); \
509           } else { \
510             ksp->reason = KSP_DIVERGED_NANORINF; \
511           } \
512           PetscFunctionReturn(PETSC_SUCCESS); \
513         } \
514       } \
515     } while (0)
516 
517   /*MC
518    KSPCheckNorm - Checks if the result of a norm used by the corresponding `KSP` contains `inf` or `NaN`. These indicate that the previous
519       application of the preconditioner generated an error. Sets a `KSPConvergedReason` and returns if the `PC` set a `PCFailedReason`.
520 
521    Collective
522 
523    Input Parameter:
524 .  ksp - the linear solver `KSP` context.
525 
526    Output Parameter:
527 .  beta - the result of the norm
528 
529    Level: developer
530 
531    Developer Notes:
532    Used to manage returning from `KSP` solvers whose preconditioners have failed, possibly only a subset of MPI ranks, in some way.
533 
534    It uses the fact that `KSP` piggy-backs the collectivity of certain error conditions on the results of norms and inner products.
535 
536 .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckDot()`, `KSPCheckSolve()`
537 M*/
538   #define KSPCheckNorm(ksp, beta) \
539     do { \
540       if (PetscIsInfOrNanReal(beta)) { \
541         PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged due to Nan or Inf norm"); \
542         { \
543           PCFailedReason pcreason; \
544           PetscCall(PCReduceFailedReason(ksp->pc)); \
545           PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason)); \
546           if (pcreason) { \
547             ksp->reason = KSP_DIVERGED_PC_FAILED; \
548             PetscCall(VecSetInf(ksp->vec_sol)); \
549             ksp->rnorm = beta; \
550           } else { \
551             ksp->reason = KSP_DIVERGED_NANORINF; \
552             ksp->rnorm  = beta; \
553           } \
554           PetscFunctionReturn(PETSC_SUCCESS); \
555         } \
556       } \
557     } while (0)
558 
559 #endif
560 
561 PETSC_INTERN PetscErrorCode KSPMonitorMakeKey_Internal(const char[], PetscViewerType, PetscViewerFormat, char[]);
562 PETSC_INTERN PetscErrorCode KSPMonitorRange_Private(KSP, PetscInt, PetscReal *);
563