Actual source code: kspimpl.h

  1: #pragma once

  3: #include <petscksp.h>
  4: #include <petscds.h>
  5: #include <petsc/private/petscimpl.h>

  7: /* SUBMANSEC = KSP */

  9: PETSC_EXTERN PetscBool      KSPRegisterAllCalled;
 10: PETSC_EXTERN PetscBool      KSPMonitorRegisterAllCalled;
 11: PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
 12: PETSC_EXTERN PetscErrorCode KSPMonitorRegisterAll(void);
 13: PETSC_EXTERN PetscErrorCode KSPGuessRegisterAll(void);
 14: PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);

 16: typedef struct _KSPOps *KSPOps;

 18: struct _KSPOps {
 19:   PetscErrorCode (*buildsolution)(KSP, Vec, Vec *);      /* Returns a pointer to the solution, or
 20:                                                           calculates the solution in a
 21:                                                           user-provided area. */
 22:   PetscErrorCode (*buildresidual)(KSP, Vec, Vec, Vec *); /* Returns a pointer to the residual, or
 23:                                                           calculates the residual in a
 24:                                                           user-provided area.  */
 25:   PetscErrorCode (*solve)(KSP);                          /* actual solver */
 26:   PetscErrorCode (*matsolve)(KSP, Mat, Mat);             /* multiple dense RHS solver */
 27:   PetscErrorCode (*setup)(KSP);
 28:   PetscErrorCode (*setfromoptions)(KSP, PetscOptionItems *);
 29:   PetscErrorCode (*publishoptions)(KSP);
 30:   PetscErrorCode (*computeextremesingularvalues)(KSP, PetscReal *, PetscReal *);
 31:   PetscErrorCode (*computeeigenvalues)(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);
 32:   PetscErrorCode (*computeritz)(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal *, PetscReal *);
 33:   PetscErrorCode (*destroy)(KSP);
 34:   PetscErrorCode (*view)(KSP, PetscViewer);
 35:   PetscErrorCode (*reset)(KSP);
 36:   PetscErrorCode (*load)(KSP, PetscViewer);
 37: };

 39: typedef struct _KSPGuessOps *KSPGuessOps;

 41: struct _KSPGuessOps {
 42:   PetscErrorCode (*formguess)(KSPGuess, Vec, Vec); /* Form initial guess */
 43:   PetscErrorCode (*update)(KSPGuess, Vec, Vec);    /* Update database */
 44:   PetscErrorCode (*setfromoptions)(KSPGuess);
 45:   PetscErrorCode (*settolerance)(KSPGuess, PetscReal);
 46:   PetscErrorCode (*setup)(KSPGuess);
 47:   PetscErrorCode (*destroy)(KSPGuess);
 48:   PetscErrorCode (*view)(KSPGuess, PetscViewer);
 49:   PetscErrorCode (*reset)(KSPGuess);
 50: };

 52: /*
 53:    Defines the KSPGuess data structure.
 54: */
 55: struct _p_KSPGuess {
 56:   PETSCHEADER(struct _KSPGuessOps);
 57:   KSP              ksp;       /* the parent KSP */
 58:   Mat              A;         /* the current linear operator */
 59:   PetscObjectState omatstate; /* previous linear operator state */
 60:   void            *data;      /* pointer to the specific implementation */
 61: };

 63: PETSC_EXTERN PetscErrorCode KSPGuessCreate_Fischer(KSPGuess);
 64: PETSC_EXTERN PetscErrorCode KSPGuessCreate_POD(KSPGuess);

 66: /*
 67:      Maximum number of monitors you can run with a single KSP
 68: */
 69: #define MAXKSPMONITORS    5
 70: #define MAXKSPREASONVIEWS 5
 71: typedef enum {
 72:   KSP_SETUP_NEW = 0,
 73:   KSP_SETUP_NEWMATRIX,
 74:   KSP_SETUP_NEWRHS
 75: } KSPSetUpStage;

 77: /*
 78:    Defines the KSP data structure.
 79: */
 80: struct _p_KSP {
 81:   PETSCHEADER(struct _KSPOps);
 82:   DM        dm;
 83:   PetscBool dmAuto;   /* DM was created automatically by KSP */
 84:   PetscBool dmActive; /* KSP should use DM for computing operators */
 85:   /*------------------------- User parameters--------------------------*/
 86:   PetscInt  max_it; /* maximum number of iterations */
 87:   PetscInt  min_it; /* minimum number of iterations */
 88:   KSPGuess  guess;
 89:   PetscBool guess_zero,                                  /* flag for whether initial guess is 0 */
 90:     guess_not_read,                                      /* guess is not read, does not need to be zeroed */
 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:   /* Allow declaring convergence when negative curvature is detected */
185:   PetscBool converged_neg_curve;

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

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

194:   PetscInt nestlevel; /* how many levels of nesting does the KSP have */
195: };

197: typedef struct { /* dummy data structure used in KSPMonitorDynamicTolerance() */
198:   PetscReal coef;
199:   PetscReal bnrm;
200: } KSPDynTolCtx;

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

209: static inline PetscErrorCode KSPLogResidualHistory(KSP ksp, PetscReal norm)
210: {
211:   PetscFunctionBegin;
212:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
213:   if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len) ksp->res_hist[ksp->res_hist_len++] = norm;
214:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: static inline PetscErrorCode KSPLogErrorHistory(KSP ksp)
219: {
220:   DM dm;

222:   PetscFunctionBegin;
223:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
224:   PetscCall(KSPGetDM(ksp, &dm));
225:   if (dm && ksp->err_hist && ksp->err_hist_max > ksp->err_hist_len) {
226:     PetscSimplePointFn *exactSol;
227:     void               *exactCtx;
228:     PetscDS             ds;
229:     Vec                 u;
230:     PetscReal           error;
231:     PetscInt            Nf;

233:     PetscCall(KSPBuildSolution(ksp, NULL, &u));
234:     /* TODO Was needed to correct for Newton solution, but I just need to set a solution */
235:     //PetscCall(VecScale(u, -1.0));
236:     /* TODO Case when I have a solution */
237:     if (0) {
238:       PetscCall(DMGetDS(dm, &ds));
239:       PetscCall(PetscDSGetNumFields(ds, &Nf));
240:       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle number of fields %" PetscInt_FMT " > 1 right now", Nf);
241:       PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol, &exactCtx));
242:       PetscCall(DMComputeL2FieldDiff(dm, 0.0, &exactSol, &exactCtx, u, &error));
243:     } else {
244:       /* The null solution A 0 = 0 */
245:       PetscCall(VecNorm(u, NORM_2, &error));
246:     }
247:     ksp->err_hist[ksp->err_hist_len++] = error;
248:   }
249:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

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

262: static inline PetscErrorCode KSPSetNoisy_Private(Vec v)
263: {
264:   PetscScalar *a;
265:   PetscInt     n, istart;

267:   PetscFunctionBegin;
268:   PetscCall(VecGetOwnershipRange(v, &istart, NULL));
269:   PetscCall(VecGetLocalSize(v, &n));
270:   PetscCall(VecGetArrayWrite(v, &a));
271:   for (PetscInt i = 0; i < n; ++i) a[i] = KSPNoisyHash_Private(i + istart);
272:   PetscCall(VecRestoreArrayWrite(v, &a));
273:   PetscFunctionReturn(PETSC_SUCCESS);
274: }

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

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

280: typedef struct _p_DMKSP  *DMKSP;
281: typedef struct _DMKSPOps *DMKSPOps;
282: struct _DMKSPOps {
283:   KSPComputeOperatorsFn    *computeoperators;
284:   KSPComputeRHSFn          *computerhs;
285:   KSPComputeInitialGuessFn *computeinitialguess;
286:   PetscErrorCode (*destroy)(DMKSP *);
287:   PetscErrorCode (*duplicate)(DMKSP, DMKSP);
288: };

290: /*S
291:    DMKSP - Object held by a `DM` that contains all the callback functions and their contexts needed by a `KSP`

293:    Level: developer

295:    Notes:
296:    Users provides callback functions and their contexts to `KSP` using, for example, `KSPSetComputeRHS()`. These values are stored
297:    in a `DMKSP` that is contained in the `DM` associated with the `KSP`. If no `DM` was provided by
298:    the user with `KSPSetDM()` it is automatically created by `KSPGetDM()` with `DMShellCreate()`.

300:    Users very rarely need to worked directly with the `DMKSP` object, rather they work with the `KSP` and the `DM` they created

302:    Multiple `DM` can share a single `DMKSP`, often each `DM` is associated with
303:    a grid refinement level. `DMGetDMKSP()` returns the `DMKSP` associated with a `DM`. `DMGetDMKSPWrite()` returns a unique
304:    `DMKSP` that is only associated with the current `DM`, making a copy of the shared `DMKSP` if needed (copy-on-write).

306:    Developer Notes:
307:    It is rather subtle why `DMKSP`, `DMSNES`, and `DMTS` are needed instead of simply storing the user callback functions and contexts in `DM` or `KSP`, `SNES`, or `TS`.
308:    It is to support composable solvers such as geometric multigrid. We want, by default, the same callback functions and contexts for all the levels in the computation,
309:    but we need to also support different callbacks and contexts on each level. The copy-on-write approach of `DMGetDMKSPWrite()` makes this possible.

311:    The `originaldm` inside the `DMKSP` is NOT reference counted (to prevent a reference count loop between a `DM` and a `DMKSP`).
312:    The `DM` on which this context was first created is cached here to implement one-way
313:    copy-on-write. When `DMGetDMKSPWrite()` sees a request using a different `DM`, it makes a copy of the `TSDM`. Thus, if a user
314:    only interacts directly with one level, e.g., using `TSSetIFunction()`, then coarse levels of a multilevel item
315:    integrator are built, then the user changes the routine with another call to `TSSetIFunction()`, it automatically
316:    propagates to all the levels. If instead, they get out a specific level and set the function on that level,
317:    subsequent changes to the original level will no longer propagate to that level.

319: .seealso: [](ch_ts), `KSP`, `KSPCreate()`, `DM`, `DMGetDMKSPWrite()`, `DMGetDMKSP()`,  `DMSNES`, `DMTS`, `DMKSPSetComputeOperators()`, `DMKSPGetComputeOperators()`,
320:           `DMKSPSetComputeRHS()`, `DMKSPSetComputeInitialGuess()`
321: S*/
322: struct _p_DMKSP {
323:   PETSCHEADER(struct _DMKSPOps);
324:   void *operatorsctx;
325:   void *rhsctx;
326:   void *initialguessctx;
327:   void *data;

329:   /* See developer note for `DMKSP` above */
330:   DM originaldm;

332:   void (*fortran_func_pointers[3])(void); /* Store our own function pointers so they are associated with the DMKSP instead of the DM */
333: };
334: PETSC_EXTERN PetscErrorCode DMGetDMKSP(DM, DMKSP *);
335: PETSC_EXTERN PetscErrorCode DMGetDMKSPWrite(DM, DMKSP *);
336: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM);

338: /*
339:        These allow the various Krylov methods to apply to either the linear system or its transpose.
340: */
341: static inline PetscErrorCode KSP_RemoveNullSpace(KSP ksp, Vec y)
342: {
343:   PetscFunctionBegin;
344:   if (ksp->pc_side == PC_LEFT) {
345:     Mat          A;
346:     MatNullSpace nullsp;

348:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
349:     PetscCall(MatGetNullSpace(A, &nullsp));
350:     if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, y));
351:   }
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: static inline PetscErrorCode KSP_RemoveNullSpaceTranspose(KSP ksp, Vec y)
356: {
357:   PetscFunctionBegin;
358:   if (ksp->pc_side == PC_LEFT) {
359:     Mat          A;
360:     MatNullSpace nullsp;

362:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
363:     PetscCall(MatGetTransposeNullSpace(A, &nullsp));
364:     if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, y));
365:   }
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: static inline PetscErrorCode KSP_MatMult(KSP ksp, Mat A, Vec x, Vec y)
370: {
371:   PetscFunctionBegin;
372:   if (ksp->transpose_solve) PetscCall(MatMultTranspose(A, x, y));
373:   else PetscCall(MatMult(A, x, y));
374:   PetscFunctionReturn(PETSC_SUCCESS);
375: }

377: static inline PetscErrorCode KSP_MatMultTranspose(KSP ksp, Mat A, Vec x, Vec y)
378: {
379:   PetscFunctionBegin;
380:   if (ksp->transpose_solve) PetscCall(MatMult(A, x, y));
381:   else PetscCall(MatMultTranspose(A, x, y));
382:   PetscFunctionReturn(PETSC_SUCCESS);
383: }

385: static inline PetscErrorCode KSP_MatMultHermitianTranspose(KSP ksp, Mat A, Vec x, Vec y)
386: {
387:   PetscFunctionBegin;
388:   if (!ksp->transpose_solve) PetscCall(MatMultHermitianTranspose(A, x, y));
389:   else {
390:     Vec w;

392:     PetscCall(VecDuplicate(x, &w));
393:     PetscCall(VecCopy(x, w));
394:     PetscCall(VecConjugate(w));
395:     PetscCall(MatMult(A, w, y));
396:     PetscCall(VecDestroy(&w));
397:     PetscCall(VecConjugate(y));
398:   }
399:   PetscFunctionReturn(PETSC_SUCCESS);
400: }

402: static inline PetscErrorCode KSP_PCApply(KSP ksp, Vec x, Vec y)
403: {
404:   PetscFunctionBegin;
405:   if (ksp->transpose_solve) {
406:     PetscCall(PCApplyTranspose(ksp->pc, x, y));
407:     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
408:   } else {
409:     PetscCall(PCApply(ksp->pc, x, y));
410:     PetscCall(KSP_RemoveNullSpace(ksp, y));
411:   }
412:   PetscFunctionReturn(PETSC_SUCCESS);
413: }

415: static inline PetscErrorCode KSP_PCApplyTranspose(KSP ksp, Vec x, Vec y)
416: {
417:   PetscFunctionBegin;
418:   if (ksp->transpose_solve) {
419:     PetscCall(PCApply(ksp->pc, x, y));
420:     PetscCall(KSP_RemoveNullSpace(ksp, y));
421:   } else {
422:     PetscCall(PCApplyTranspose(ksp->pc, x, y));
423:     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
424:   }
425:   PetscFunctionReturn(PETSC_SUCCESS);
426: }

428: static inline PetscErrorCode KSP_PCApplyHermitianTranspose(KSP ksp, Vec x, Vec y)
429: {
430:   PetscFunctionBegin;
431:   PetscCall(VecConjugate(x));
432:   PetscCall(KSP_PCApplyTranspose(ksp, x, y));
433:   PetscCall(VecConjugate(x));
434:   PetscCall(VecConjugate(y));
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: static inline PetscErrorCode KSP_PCMatApply(KSP ksp, Mat X, Mat Y)
439: {
440:   PetscFunctionBegin;
441:   if (ksp->transpose_solve) {
442:     PetscBool flg;
443:     PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp->pc, &flg, PCNONE, PCICC, PCCHOLESKY, ""));
444:     PetscCheck(flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "PCMatApplyTranspose() not yet implemented for nonsymmetric PC");
445:   }
446:   PetscCall(PCMatApply(ksp->pc, X, Y));
447:   PetscFunctionReturn(PETSC_SUCCESS);
448: }

450: static inline PetscErrorCode KSP_PCMatApplyTranspose(KSP ksp, Mat X, Mat Y)
451: {
452:   PetscFunctionBegin;
453:   if (!ksp->transpose_solve) {
454:     PetscBool flg;
455:     PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp->pc, &flg, PCNONE, PCICC, PCCHOLESKY, ""));
456:     PetscCheck(flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "PCMatApplyTranspose() not yet implemented for nonsymmetric PC");
457:   }
458:   PetscCall(PCMatApply(ksp->pc, X, Y));
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }

462: static inline PetscErrorCode KSP_PCApplyBAorAB(KSP ksp, Vec x, Vec y, Vec w)
463: {
464:   PetscFunctionBegin;
465:   if (ksp->transpose_solve) {
466:     PetscCall(PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w));
467:     PetscCall(KSP_RemoveNullSpaceTranspose(ksp, y));
468:   } else {
469:     PetscCall(PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w));
470:     PetscCall(KSP_RemoveNullSpace(ksp, y));
471:   }
472:   PetscFunctionReturn(PETSC_SUCCESS);
473: }

475: static inline PetscErrorCode KSP_PCApplyBAorABTranspose(KSP ksp, Vec x, Vec y, Vec w)
476: {
477:   PetscFunctionBegin;
478:   if (ksp->transpose_solve) PetscCall(PCApplyBAorAB(ksp->pc, ksp->pc_side, x, y, w));
479:   else PetscCall(PCApplyBAorABTranspose(ksp->pc, ksp->pc_side, x, y, w));
480:   PetscFunctionReturn(PETSC_SUCCESS);
481: }

483: PETSC_EXTERN PetscLogEvent KSP_GMRESOrthogonalization;
484: PETSC_EXTERN PetscLogEvent KSP_SetUp;
485: PETSC_EXTERN PetscLogEvent KSP_Solve;
486: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_0;
487: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_1;
488: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_2;
489: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_3;
490: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_4;
491: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_S;
492: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_L;
493: PETSC_EXTERN PetscLogEvent KSP_Solve_FS_U;
494: PETSC_EXTERN PetscLogEvent KSP_SolveTranspose;
495: PETSC_EXTERN PetscLogEvent KSP_MatSolve;
496: PETSC_EXTERN PetscLogEvent KSP_MatSolveTranspose;

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

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

505:    Collective

507:    Input Parameter:
508: .  ksp - the linear solver `KSP` context.

510:    Output Parameter:
511: .  beta - the result of the inner product

513:    Level: developer

515:    Developer Notes:
516:    Used to manage returning from `KSP` solvers collectively whose preconditioners have failed, possibly only a subset of MPI processes, in some way

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

520: .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckNorm()`, `KSPCheckSolve()`,
521:           `KSPSetErrorIfNotConverged()`
522: M*/
523: #define KSPCheckDot(ksp, beta) \
524:   do { \
525:     if (PetscIsInfOrNanScalar(beta)) { \
526:       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged due to Nan or Inf inner product"); \
527:       { \
528:         PCFailedReason pcreason; \
529:         PetscCall(PCReduceFailedReason(ksp->pc)); \
530:         PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason)); \
531:         if (pcreason) { \
532:           ksp->reason = KSP_DIVERGED_PC_FAILED; \
533:           PetscCall(VecSetInf(ksp->vec_sol)); \
534:         } else { \
535:           ksp->reason = KSP_DIVERGED_NANORINF; \
536:         } \
537:         PetscFunctionReturn(PETSC_SUCCESS); \
538:       } \
539:     } \
540:   } while (0)

542: /*MC
543:    KSPCheckNorm - Checks if the result of a norm used by the corresponding `KSP` contains `inf` or `NaN`. These indicate that the previous
544:       application of the preconditioner generated an error. Sets a `KSPConvergedReason` and returns if the `PC` set a `PCFailedReason`.

546:    Collective

548:    Input Parameter:
549: .  ksp - the linear solver `KSP` context.

551:    Output Parameter:
552: .  beta - the result of the norm

554:    Level: developer

556:    Developer Notes:
557:    Used to manage returning from `KSP` solvers collectively whose preconditioners have failed, possibly only a subset of MPI processes, in some way.

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

561: .seealso: `PCFailedReason`, `KSPConvergedReason`, `PCGetFailedReasonRank()`, `KSP`, `KSPCreate()`, `KSPSetType()`, `KSP`, `KSPCheckDot()`, `KSPCheckSolve()`,
562:           `KSPSetErrorIfNotConverged()`
563: M*/
564: #define KSPCheckNorm(ksp, beta) \
565:   do { \
566:     if (PetscIsInfOrNanReal(beta)) { \
567:       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged due to Nan or Inf norm"); \
568:       { \
569:         PCFailedReason pcreason; \
570:         PetscCall(PCReduceFailedReason(ksp->pc)); \
571:         PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason)); \
572:         if (pcreason) { \
573:           ksp->reason = KSP_DIVERGED_PC_FAILED; \
574:           PetscCall(VecSetInf(ksp->vec_sol)); \
575:           ksp->rnorm = beta; \
576:         } else { \
577:           ksp->reason = KSP_DIVERGED_NANORINF; \
578:           ksp->rnorm  = beta; \
579:         } \
580:         PetscFunctionReturn(PETSC_SUCCESS); \
581:       } \
582:     } \
583:   } while (0)

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