Actual source code: kspimpl.h


  2: #ifndef _KSPIMPL_H
  3:   #define _KSPIMPL_H

  5: #include <petscksp.h>
  6: #include <petscds.h>
  7: #include <petsc/private/petscimpl.h>

  9: /* SUBMANSEC = KSP */

 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);

 18: typedef struct _KSPOps *KSPOps;

 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: };

 41: typedef struct _KSPGuessOps *KSPGuessOps;

 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: };

 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: };

 65: PETSC_EXTERN PetscErrorCode KSPGuessCreate_Fischer(KSPGuess);
 66: PETSC_EXTERN PetscErrorCode KSPGuessCreate_POD(KSPGuess);

 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,
 75:   KSP_SETUP_NEWMATRIX,
 76:   KSP_SETUP_NEWRHS
 77: } KSPSetUpStage;

 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:   KSPGuess  guess;
 90:   PetscBool guess_zero,                                  /* flag for whether initial guess is 0 */
 91:     calc_sings,                                          /* calculate extreme Singular Values */
 92:     calc_ritz,                                           /* calculate (harmonic) Ritz pairs */
 93:     guess_knoll;                                         /* use initial guess of PCApply(ksp->B,b */
 94:   PCSide    pc_side;                                     /* flag for left, right, or symmetric preconditioning */
 95:   PetscInt  normsupporttable[KSP_NORM_MAX][PC_SIDE_MAX]; /* Table of supported norms and pc_side, see KSPSetSupportedNorm() */
 96:   PetscReal rtol,                                        /* relative tolerance */
 97:     abstol,                                              /* absolute tolerance */
 98:     ttol,                                                /* (not set by user)  */
 99:     divtol;                                              /* divergence tolerance */
100:   PetscReal          rnorm0;                             /* initial residual norm (used for divergence testing) */
101:   PetscReal          rnorm;                              /* current residual norm */
102:   KSPConvergedReason reason;
103:   PetscBool          errorifnotconverged; /* create an error if the KSPSolve() does not converge */

105:   Vec        vec_sol, vec_rhs; /* pointer to where user has stashed
106:                                       the solution and rhs, these are
107:                                       never touched by the code, only
108:                                       passed back to the user */
109:   PetscReal *res_hist;         /* If !0 stores residual each at iteration */
110:   PetscReal *res_hist_alloc;   /* If !0 means user did not provide buffer, needs deallocation */
111:   size_t     res_hist_len;     /* current size of residual history array */
112:   size_t     res_hist_max;     /* actual amount of storage in residual history */
113:   PetscBool  res_hist_reset;   /* reset history to length zero for each new solve */
114:   PetscReal *err_hist;         /* If !0 stores error at each iteration */
115:   PetscReal *err_hist_alloc;   /* If !0 means user did not provide buffer, needs deallocation */
116:   size_t     err_hist_len;     /* current size of error history array */
117:   size_t     err_hist_max;     /* actual amount of storage in error history */
118:   PetscBool  err_hist_reset;   /* reset history to length zero for each new solve */

120:   PetscInt  chknorm; /* only compute/check norm if iterations is great than this */
121:   PetscBool lagnorm; /* Lag the residual norm calculation so that it is computed as part of the
122:                                         MPI_Allreduce() for computing the inner products for the next iteration. */

124:   PetscInt nmax; /* maximum number of right-hand sides to be handled simultaneously */

126:   /* --------User (or default) routines (most return -1 on error) --------*/
127:   PetscErrorCode (*monitor[MAXKSPMONITORS])(KSP, PetscInt, PetscReal, void *); /* returns control to user after */
128:   PetscErrorCode (*monitordestroy[MAXKSPMONITORS])(void **);                   /* */
129:   void     *monitorcontext[MAXKSPMONITORS];                                    /* residual calculation, allows user */
130:   PetscInt  numbermonitors;                                                    /* to, for instance, print residual norm, etc. */
131:   PetscBool pauseFinal;                                                        /* Pause all drawing monitor at the final iterate */

133:   PetscErrorCode (*reasonview[MAXKSPREASONVIEWS])(KSP, void *);    /* KSP converged reason view */
134:   PetscErrorCode (*reasonviewdestroy[MAXKSPREASONVIEWS])(void **); /* Optional destroy routine */
135:   void    *reasonviewcontext[MAXKSPREASONVIEWS];                   /* User context */
136:   PetscInt numberreasonviews;                                      /* Number if reason viewers */

138:   PetscErrorCode (*converged)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
139:   PetscErrorCode (*convergeddestroy)(void *);
140:   void *cnvP;

142:   void *user; /* optional user-defined context */

144:   PC pc;

146:   void *data; /* holder for misc stuff associated
147:                                    with a particular iterative solver */

149:   PetscBool         view, viewPre, viewRate, viewMat, viewPMat, viewRhs, viewSol, viewMatExp, viewEV, viewSV, viewEVExp, viewFinalRes, viewPOpExp, viewDScale;
150:   PetscViewer       viewer, viewerPre, viewerRate, viewerMat, viewerPMat, viewerRhs, viewerSol, viewerMatExp, viewerEV, viewerSV, viewerEVExp, viewerFinalRes, viewerPOpExp, viewerDScale;
151:   PetscViewerFormat format, formatPre, formatRate, formatMat, formatPMat, formatRhs, formatSol, formatMatExp, formatEV, formatSV, formatEVExp, formatFinalRes, formatPOpExp, formatDScale;

153:   /* ----------------Default work-area management -------------------- */
154:   PetscInt nwork;
155:   Vec     *work;

157:   KSPSetUpStage setupstage;
158:   PetscBool     setupnewmatrix; /* true if we need to call ksp->ops->setup with KSP_SETUP_NEWMATRIX */

160:   PetscInt its;      /* number of iterations so far computed in THIS linear solve*/
161:   PetscInt totalits; /* number of iterations used by this KSP object since it was created */

163:   PetscBool transpose_solve; /* solve transpose system instead */
164:   struct {
165:     Mat       AT, BT;
166:     PetscBool use_explicittranspose; /* transpose the system explicitly in KSPSolveTranspose */
167:     PetscBool reuse_transpose;       /* reuse the previous transposed system */
168:   } transpose;

170:   KSPNormType normtype; /* type of norm used for convergence tests */

172:   PCSide      pc_side_set;  /* PC type set explicitly by user */
173:   KSPNormType normtype_set; /* Norm type set explicitly by user */

175:   /*   Allow diagonally scaling the matrix before computing the preconditioner or using
176:        the Krylov method. Note this is NOT just Jacobi preconditioning */

178:   PetscBool dscale;     /* diagonal scale system; used with KSPSetDiagonalScale() */
179:   PetscBool dscalefix;  /* unscale system after solve */
180:   PetscBool dscalefix2; /* system has been unscaled */
181:   Vec       diagonal;   /* 1/sqrt(diag of matrix) */
182:   Vec       truediagonal;

184:   PetscInt  setfromoptionscalled;
185:   PetscBool skippcsetfromoptions; /* if set then KSPSetFromOptions() does not call PCSetFromOptions() */

187:   PetscViewer eigviewer; /* Viewer where computed eigenvalues are displayed */

189:   PetscErrorCode (*presolve)(KSP, Vec, Vec, void *);
190:   PetscErrorCode (*postsolve)(KSP, Vec, Vec, void *);
191:   void *prectx, *postctx;
192: };

194: typedef struct { /* dummy data structure used in KSPMonitorDynamicTolerance() */
195:   PetscReal coef;
196:   PetscReal bnrm;
197: } KSPDynTolCtx;

199: typedef struct {
200:   PetscBool initialrtol;    /* default relative residual decrease is computed from initial residual, not rhs */
201:   PetscBool mininitialrtol; /* default relative residual decrease is computed from min of initial residual and rhs */
202:   PetscBool convmaxits;     /* if true, the convergence test returns KSP_CONVERGED_ITS if the maximum number of iterations is reached */
203:   Vec       work;
204: } KSPConvergedDefaultCtx;

206: static inline PetscErrorCode KSPLogResidualHistory(KSP ksp, PetscReal norm)
207: {
208:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
209:   if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len) ksp->res_hist[ksp->res_hist_len++] = norm;
210:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
211:   return 0;
212: }

214: static inline PetscErrorCode KSPLogErrorHistory(KSP ksp)
215: {
216:   DM dm;

218:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
219:   KSPGetDM(ksp, &dm);
220:   if (dm && ksp->err_hist && ksp->err_hist_max > ksp->err_hist_len) {
221:     PetscSimplePointFunc exactSol;
222:     void                *exactCtx;
223:     PetscDS              ds;
224:     Vec                  u;
225:     PetscReal            error;
226:     PetscInt             Nf;

228:     KSPBuildSolution(ksp, NULL, &u);
229:     /* TODO Was needed to correct for Newton solution, but I just need to set a solution */
230:     //VecScale(u, -1.0);
231:     /* TODO Case when I have a solution */
232:     if (0) {
233:       DMGetDS(dm, &ds);
234:       PetscDSGetNumFields(ds, &Nf);
236:       PetscDSGetExactSolution(ds, 0, &exactSol, &exactCtx);
237:       DMComputeL2FieldDiff(dm, 0.0, &exactSol, &exactCtx, u, &error);
238:     } else {
239:       /* The null solution A 0 = 0 */
240:       VecNorm(u, NORM_2, &error);
241:     }
242:     ksp->err_hist[ksp->err_hist_len++] = error;
243:   }
244:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
245:   return 0;
246: }

248: static inline PetscScalar KSPNoisyHash_Private(PetscInt xx)
249: {
250:   unsigned int x = (unsigned int)xx;
251:   x              = ((x >> 16) ^ x) * 0x45d9f3b;
252:   x              = ((x >> 16) ^ x) * 0x45d9f3b;
253:   x              = ((x >> 16) ^ x);
254:   return (PetscScalar)((PetscInt64)x - 2147483648) * 5.e-10; /* center around zero, scaled about -1. to 1.*/
255: }

257: static inline PetscErrorCode KSPSetNoisy_Private(Vec v)
258: {
259:   PetscScalar *a;
260:   PetscInt     n, istart;

262:   VecGetOwnershipRange(v, &istart, NULL);
263:   VecGetLocalSize(v, &n);
264:   VecGetArrayWrite(v, &a);
265:   for (PetscInt i = 0; i < n; ++i) a[i] = KSPNoisyHash_Private(i + istart);
266:   VecRestoreArrayWrite(v, &a);
267:   return 0;
268: }

270: PETSC_INTERN PetscErrorCode KSPSetUpNorms_Private(KSP, PetscBool, KSPNormType *, PCSide *);

272: PETSC_INTERN PetscErrorCode KSPPlotEigenContours_Private(KSP, PetscInt, const PetscReal *, const PetscReal *);

274: typedef struct _p_DMKSP  *DMKSP;
275: typedef struct _DMKSPOps *DMKSPOps;
276: struct _DMKSPOps {
277:   PetscErrorCode (*computeoperators)(KSP, Mat, Mat, void *);
278:   PetscErrorCode (*computerhs)(KSP, Vec, void *);
279:   PetscErrorCode (*computeinitialguess)(KSP, Vec, void *);
280:   PetscErrorCode (*destroy)(DMKSP *);
281:   PetscErrorCode (*duplicate)(DMKSP, DMKSP);
282: };

284: struct _p_DMKSP {
285:   PETSCHEADER(struct _DMKSPOps);
286:   void *operatorsctx;
287:   void *rhsctx;
288:   void *initialguessctx;
289:   void *data;

291:   /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
292:    * copy-on-write. When DMGetDMKSPWrite() sees a request using a different DM, it makes a copy. Thus, if a user
293:    * only interacts directly with one level, e.g., using KSPSetComputeOperators(), then coarse levels are constructed by
294:    * PCMG, then the user changes the routine with another call to KSPSetComputeOperators(), it automatically propagates
295:    * to all the levels. If instead, they get out a specific level and set the routines on that level, subsequent changes
296:    * to the original level will no longer propagate to that level.
297:    */
298:   DM originaldm;

300:   void (*fortran_func_pointers[3])(void); /* Store our own function pointers so they are associated with the DMKSP instead of the DM */
301: };
302: PETSC_EXTERN PetscErrorCode DMGetDMKSP(DM, DMKSP *);
303: PETSC_EXTERN PetscErrorCode DMGetDMKSPWrite(DM, DMKSP *);
304: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);

306: /*
307:        These allow the various Krylov methods to apply to either the linear system or its transpose.
308: */
309: static inline PetscErrorCode KSP_RemoveNullSpace(KSP ksp, Vec y)
310: {
311:   if (ksp->pc_side == PC_LEFT) {
312:     Mat          A;
313:     MatNullSpace nullsp;

315:     PCGetOperators(ksp->pc, &A, NULL);
316:     MatGetNullSpace(A, &nullsp);
317:     if (nullsp) MatNullSpaceRemove(nullsp, y);
318:   }
319:   return 0;
320: }

322: static inline PetscErrorCode KSP_RemoveNullSpaceTranspose(KSP ksp, Vec y)
323: {
324:   if (ksp->pc_side == PC_LEFT) {
325:     Mat          A;
326:     MatNullSpace nullsp;

328:     PCGetOperators(ksp->pc, &A, NULL);
329:     MatGetTransposeNullSpace(A, &nullsp);
330:     if (nullsp) MatNullSpaceRemove(nullsp, y);
331:   }
332:   return 0;
333: }

335: static inline PetscErrorCode KSP_MatMult(KSP ksp, Mat A, Vec x, Vec y)
336: {
337:   if (ksp->transpose_solve) MatMultTranspose(A, x, y);
338:   else MatMult(A, x, y);
339:   return 0;
340: }

342: static inline PetscErrorCode KSP_MatMultTranspose(KSP ksp, Mat A, Vec x, Vec y)
343: {
344:   if (ksp->transpose_solve) MatMult(A, x, y);
345:   else MatMultTranspose(A, x, y);
346:   return 0;
347: }

349: static inline PetscErrorCode KSP_MatMultHermitianTranspose(KSP ksp, Mat A, Vec x, Vec y)
350: {
351:   if (!ksp->transpose_solve) MatMultHermitianTranspose(A, x, y);
352:   else {
353:     Vec w;

355:     VecDuplicate(x, &w);
356:     VecCopy(x, w);
357:     VecConjugate(w);
358:     MatMult(A, w, y);
359:     VecDestroy(&w);
360:     VecConjugate(y);
361:   }
362:   return 0;
363: }

365: static inline PetscErrorCode KSP_PCApply(KSP ksp, Vec x, Vec y)
366: {
367:   if (ksp->transpose_solve) {
368:     PCApplyTranspose(ksp->pc, x, y);
369:     KSP_RemoveNullSpaceTranspose(ksp, y);
370:   } else {
371:     PCApply(ksp->pc, x, y);
372:     KSP_RemoveNullSpace(ksp, y);
373:   }
374:   return 0;
375: }

377: static inline PetscErrorCode KSP_PCApplyTranspose(KSP ksp, Vec x, Vec y)
378: {
379:   if (ksp->transpose_solve) {
380:     PCApply(ksp->pc, x, y);
381:     KSP_RemoveNullSpace(ksp, y);
382:   } else {
383:     PCApplyTranspose(ksp->pc, x, y);
384:     KSP_RemoveNullSpaceTranspose(ksp, y);
385:   }
386:   return 0;
387: }

389: static inline PetscErrorCode KSP_PCApplyHermitianTranspose(KSP ksp, Vec x, Vec y)
390: {
391:   VecConjugate(x);
392:   KSP_PCApplyTranspose(ksp, x, y);
393:   VecConjugate(x);
394:   VecConjugate(y);
395:   return 0;
396: }

398: static inline PetscErrorCode KSP_PCApplyBAorAB(KSP ksp, Vec x, Vec y, Vec w)
399: {
400:   if (ksp->transpose_solve) {
401:     PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w);
402:     KSP_RemoveNullSpaceTranspose(ksp, y);
403:   } else {
404:     PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w);
405:     KSP_RemoveNullSpace(ksp, y);
406:   }
407:   return 0;
408: }

410: static inline PetscErrorCode KSP_PCApplyBAorABTranspose(KSP ksp, Vec x, Vec y, Vec w)
411: {
412:   if (ksp->transpose_solve) PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w);
413:   else PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w);
414:   return 0;
415: }

417: PETSC_EXTERN PetscLogEvent KSP_GMRESOrthogonalization;
418: PETSC_EXTERN PetscLogEvent KSP_SetUp;
419: PETSC_EXTERN PetscLogEvent KSP_Solve;
420: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_0;
421: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_1;
422: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_2;
423: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_3;
424: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_4;
425: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_S;
426: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_L;
427: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_U;
428: PETSC_EXTERN PetscLogEvent KSP_SolveTranspose;
429: PETSC_EXTERN PetscLogEvent KSP_MatSolve;

431: PETSC_INTERN PetscErrorCode MatGetSchurComplement_Basic(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *);
432: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *);

434:   /*MC
435:    KSPCheckDot - Checks if the result of a dot product used by the corresponding `KSP` contains Inf or NaN. These indicate that the previous
436:       application of the preconditioner generated an error. Sets a `KSPConvergedReason` and returns if the `PC` set a `PCFailedReason`.

438:    Collective

440:    Input Parameter:
441: .  ksp - the linear solver `KSP` context.

443:    Output Parameter:
444: .  beta - the result of the inner product

446:    Level: developer

448:    Developer Notes:
449:    Used to manage returning from `KSP` solvers whose preconditioners have failed, possibly only a subset of MPI ranks, in some way

451:    It uses the fact that `KSP` piggy-backs the collectivity of certain error conditions on the results of norms and inner products.

453: .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckNorm()`, `KSPCheckSolve()`
454: M*/
455:   #define KSPCheckDot(ksp, beta) \
456:     do { \
457:       if (PetscIsInfOrNanScalar(beta)) { \
459:         { \
460:           PCFailedReason pcreason; \
461:           PetscInt       sendbuf, recvbuf; \
462:           PCGetFailedReasonRank(ksp->pc, &pcreason); \
463:           sendbuf = (PetscInt)pcreason; \
464:           MPI_Allreduce(&sendbuf, &recvbuf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ksp)); \
465:           if (recvbuf) { \
466:             PCSetFailedReason(ksp->pc, (PCFailedReason)recvbuf); \
467:             ksp->reason = KSP_DIVERGED_PC_FAILED; \
468:             VecSetInf(ksp->vec_sol); \
469:           } else { \
470:             ksp->reason = KSP_DIVERGED_NANORINF; \
471:           } \
472:           return 0; \
473:         } \
474:       } \
475:     } while (0)

477:   /*MC
478:    KSPCheckNorm - Checks if the result of a norm 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`.

481:    Collective

483:    Input Parameter:
484: .  ksp - the linear solver `KSP` context.

486:    Output Parameter:
487: .  beta - the result of the norm

489:    Level: developer

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.

494:    It uses the fact that `KSP` piggy-backs the collectivity of certain error conditions on the results of norms and inner products.

496: .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckDot()`, `KSPCheckSolve()`
497: M*/
498:   #define KSPCheckNorm(ksp, beta) \
499:     do { \
500:       if (PetscIsInfOrNanReal(beta)) { \
502:         { \
503:           PCFailedReason pcreason; \
504:           PetscInt       sendbuf, recvbuf; \
505:           PCGetFailedReasonRank(ksp->pc, &pcreason); \
506:           sendbuf = (PetscInt)pcreason; \
507:           MPI_Allreduce(&sendbuf, &recvbuf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ksp)); \
508:           if (recvbuf) { \
509:             PCSetFailedReason(ksp->pc, (PCFailedReason)recvbuf); \
510:             ksp->reason = KSP_DIVERGED_PC_FAILED; \
511:             VecSetInf(ksp->vec_sol); \
512:             ksp->rnorm = beta; \
513:           } else { \
514:             PCSetFailedReason(ksp->pc, PC_NOERROR); \
515:             ksp->reason = KSP_DIVERGED_NANORINF; \
516:             ksp->rnorm  = beta; \
517:           } \
518:           return 0; \
519:         } \
520:       } \
521:     } while (0)

523: #endif

525: PETSC_INTERN PetscErrorCode KSPMonitorMakeKey_Internal(const char[], PetscViewerType, PetscViewerFormat, char[]);
526: PETSC_INTERN PetscErrorCode KSPMonitorRange_Private(KSP, PetscInt, PetscReal *);