Actual source code: itcreate.c

  1: /*
  2:      The basic KSP routines, Create, View etc. are here.
  3: */
  4: #include <petsc/private/kspimpl.h>

  6: /* Logging support */
  7: PetscClassId  KSP_CLASSID;
  8: PetscClassId  DMKSP_CLASSID;
  9: PetscClassId  KSPGUESS_CLASSID;
 10: PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose, KSP_MatSolve, KSP_MatSolveTranspose;

 12: /*
 13:    Contains the list of registered KSP routines
 14: */
 15: PetscFunctionList KSPList              = NULL;
 16: PetscBool         KSPRegisterAllCalled = PETSC_FALSE;

 18: /*
 19:    Contains the list of registered KSP monitors
 20: */
 21: PetscFunctionList KSPMonitorList              = NULL;
 22: PetscFunctionList KSPMonitorCreateList        = NULL;
 23: PetscFunctionList KSPMonitorDestroyList       = NULL;
 24: PetscBool         KSPMonitorRegisterAllCalled = PETSC_FALSE;

 26: /*@
 27:   KSPLoad - Loads a `KSP` that has been stored in a `PETSCVIEWERBINARY`  with `KSPView()`.

 29:   Collective

 31:   Input Parameters:
 32: + newdm  - the newly loaded `KSP`, this needs to have been created with `KSPCreate()` or
 33:            some related function before a call to `KSPLoad()`.
 34: - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`

 36:   Level: intermediate

 38:   Note:
 39:   The type is determined by the data in the file, any type set into the `KSP` before this call is ignored.

 41: .seealso: [](ch_ksp), `KSP`, `PetscViewerBinaryOpen()`, `KSPView()`, `MatLoad()`, `VecLoad()`
 42: @*/
 43: PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer)
 44: {
 45:   PetscBool isbinary;
 46:   PetscInt  classid;
 47:   char      type[256];
 48:   PC        pc;

 50:   PetscFunctionBegin;
 53:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
 54:   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");

 56:   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
 57:   PetscCheck(classid == KSP_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not KSP next in file");
 58:   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
 59:   PetscCall(KSPSetType(newdm, type));
 60:   PetscTryTypeMethod(newdm, load, viewer);
 61:   PetscCall(KSPGetPC(newdm, &pc));
 62:   PetscCall(PCLoad(pc, viewer));
 63:   PetscFunctionReturn(PETSC_SUCCESS);
 64: }

 66: #include <petscdraw.h>
 67: #if defined(PETSC_HAVE_SAWS)
 68: #include <petscviewersaws.h>
 69: #endif
 70: /*@
 71:   KSPView - Prints the various parameters currently set in the `KSP` object. For example, the convergence tolerances and `KSPType`.
 72:   Also views the `PC` and `Mat` contained by the `KSP` with `PCView()` and `MatView()`.

 74:   Collective

 76:   Input Parameters:
 77: + ksp    - the Krylov space context
 78: - viewer - visualization context

 80:   Options Database Key:
 81: . -ksp_view - print the `KSP` data structure at the end of each `KSPSolve()` call

 83:   Level: beginner

 85:   Notes:
 86:   The available visualization contexts include
 87: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
 88: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
 89:   output where only the first processor opens
 90:   the file.  All other processors send their
 91:   data to the first processor to print.

 93:   The available formats include
 94: +     `PETSC_VIEWER_DEFAULT` - standard output (default)
 95: -     `PETSC_VIEWER_ASCII_INFO_DETAIL` - more verbose output for `PCBJACOBI` and `PCASM`

 97:   The user can open an alternative visualization context with
 98:   `PetscViewerASCIIOpen()` - output to a specified file.

100:   Use `KSPViewFromOptions()` to allow the user to select many different `PetscViewerType` and formats from the options database.

102:   In the debugger you can do call `KSPView(ksp,0)` to display the `KSP`. (The same holds for any PETSc object viewer).

104: .seealso: [](ch_ksp), `KSP`, `PetscViewer`, `PCView()`, `PetscViewerASCIIOpen()`, `KSPViewFromOptions()`
105: @*/
106: PetscErrorCode KSPView(KSP ksp, PetscViewer viewer)
107: {
108:   PetscBool isascii, isbinary, isdraw, isstring;
109: #if defined(PETSC_HAVE_SAWS)
110:   PetscBool issaws;
111: #endif

113:   PetscFunctionBegin;
115:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp), &viewer));
117:   PetscCheckSameComm(ksp, 1, viewer, 2);

119:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
120:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
121:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
122:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
123: #if defined(PETSC_HAVE_SAWS)
124:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
125: #endif
126:   if (isascii) {
127:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ksp, viewer));
128:     PetscCall(PetscViewerASCIIPushTab(viewer));
129:     PetscTryTypeMethod(ksp, view, viewer);
130:     PetscCall(PetscViewerASCIIPopTab(viewer));
131:     if (ksp->guess_zero) {
132:       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT ", initial guess is zero\n", ksp->max_it));
133:     } else {
134:       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT ", nonzero initial guess\n", ksp->max_it));
135:     }
136:     if (ksp->min_it) PetscCall(PetscViewerASCIIPrintf(viewer, "  minimum iterations=%" PetscInt_FMT "\n", ksp->min_it));
137:     if (ksp->guess_knoll) PetscCall(PetscViewerASCIIPrintf(viewer, "  using preconditioner applied to right-hand side for initial guess\n"));
138:     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%g, absolute=%g, divergence=%g\n", (double)ksp->rtol, (double)ksp->abstol, (double)ksp->divtol));
139:     if (ksp->pc_side == PC_RIGHT) {
140:       PetscCall(PetscViewerASCIIPrintf(viewer, "  right preconditioning\n"));
141:     } else if (ksp->pc_side == PC_SYMMETRIC) {
142:       PetscCall(PetscViewerASCIIPrintf(viewer, "  symmetric preconditioning\n"));
143:     } else {
144:       PetscCall(PetscViewerASCIIPrintf(viewer, "  left preconditioning\n"));
145:     }
146:     if (ksp->guess) {
147:       PetscCall(PetscViewerASCIIPushTab(viewer));
148:       PetscCall(KSPGuessView(ksp->guess, viewer));
149:       PetscCall(PetscViewerASCIIPopTab(viewer));
150:     }
151:     if (ksp->dscale) PetscCall(PetscViewerASCIIPrintf(viewer, "  diagonally scaled system\n"));
152:     if (ksp->converged == KSPConvergedSkip || ksp->normtype == KSP_NORM_NONE) PetscCall(PetscViewerASCIIPrintf(viewer, "  not checking for convergence\n"));
153:     else PetscCall(PetscViewerASCIIPrintf(viewer, "  using %s norm type for convergence test\n", KSPNormTypes[ksp->normtype]));
154:   } else if (isbinary) {
155:     PetscInt    classid = KSP_FILE_CLASSID;
156:     MPI_Comm    comm;
157:     PetscMPIInt rank;
158:     char        type[256];

160:     PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
161:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
162:     if (rank == 0) {
163:       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
164:       PetscCall(PetscStrncpy(type, ((PetscObject)ksp)->type_name, 256));
165:       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
166:     }
167:     PetscTryTypeMethod(ksp, view, viewer);
168:   } else if (isstring) {
169:     const char *type;
170:     PetscCall(KSPGetType(ksp, &type));
171:     PetscCall(PetscViewerStringSPrintf(viewer, " KSPType: %-7.7s", type));
172:     PetscTryTypeMethod(ksp, view, viewer);
173:   } else if (isdraw) {
174:     PetscDraw draw;
175:     char      str[36];
176:     PetscReal x, y, bottom, h;
177:     PetscBool flg;

179:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
180:     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
181:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg));
182:     if (!flg) {
183:       PetscCall(PetscStrncpy(str, "KSP: ", sizeof(str)));
184:       PetscCall(PetscStrlcat(str, ((PetscObject)ksp)->type_name, sizeof(str)));
185:       PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
186:       bottom = y - h;
187:     } else {
188:       bottom = y;
189:     }
190:     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
191: #if defined(PETSC_HAVE_SAWS)
192:   } else if (issaws) {
193:     PetscMPIInt rank;
194:     const char *name;

196:     PetscCall(PetscObjectGetName((PetscObject)ksp, &name));
197:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
198:     if (!((PetscObject)ksp)->amsmem && rank == 0) {
199:       char dir[1024];

201:       PetscCall(PetscObjectViewSAWs((PetscObject)ksp, viewer));
202:       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/its", name));
203:       PetscCallSAWs(SAWs_Register, (dir, &ksp->its, 1, SAWs_READ, SAWs_INT));
204:       if (!ksp->res_hist) PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE));
205:       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/res_hist", name));
206:       PetscCallSAWs(SAWs_Register, (dir, ksp->res_hist, 10, SAWs_READ, SAWs_DOUBLE));
207:     }
208: #endif
209:   } else PetscTryTypeMethod(ksp, view, viewer);
210:   if (ksp->pc) PetscCall(PCView(ksp->pc, viewer));
211:   if (isdraw) {
212:     PetscDraw draw;
213:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
214:     PetscCall(PetscDrawPopCurrentPoint(draw));
215:   }
216:   PetscFunctionReturn(PETSC_SUCCESS);
217: }

219: /*@
220:   KSPViewFromOptions - View (print) a `KSP` object based on values in the options database. Also views the `PC` and `Mat` contained by the `KSP`
221:   with `PCView()` and `MatView()`.

223:   Collective

225:   Input Parameters:
226: + A    - Krylov solver context
227: . obj  - Optional object that provides the options prefix used to query the options database
228: - name - command line option

230:   Level: intermediate

232: .seealso: [](ch_ksp), `KSP`, `KSPView()`, `PetscObjectViewFromOptions()`, `KSPCreate()`
233: @*/
234: PetscErrorCode KSPViewFromOptions(KSP A, PetscObject obj, const char name[])
235: {
236:   PetscFunctionBegin;
238:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
239:   PetscFunctionReturn(PETSC_SUCCESS);
240: }

242: /*@
243:   KSPSetNormType - Sets the type of residual norm that is used for convergence testing in `KSPSolve()` for the given `KSP` context

245:   Logically Collective

247:   Input Parameters:
248: + ksp      - Krylov solver context
249: - normtype - one of
250: .vb
251:    KSP_NORM_NONE             - skips computing the norm, this should generally only be used if you are using
252:                                the Krylov method as a smoother with a fixed small number of iterations.
253:                                Implicitly sets `KSPConvergedSkip()` as the `KSP` convergence test.
254:                                Note that certain algorithms such as `KSPGMRES` ALWAYS require the norm calculation,
255:                                for these methods the norms are still computed, they are just not used in
256:                                the convergence test.
257:    KSP_NORM_PRECONDITIONED   - the default for left-preconditioned solves, uses the 2-norm
258:                                of the preconditioned residual  $B^{-1}(b - A x)$.
259:    KSP_NORM_UNPRECONDITIONED - uses the 2-norm of the true $b - Ax$ residual.
260:    KSP_NORM_NATURAL          - uses the $A$ norm of the true $b - Ax$ residual; supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`
261: .ve

263:   Options Database Key:
264: . -ksp_norm_type (none|preconditioned|unpreconditioned|natural) - set `KSP` norm type

266:   Level: advanced

268:   Notes:
269:   The norm is always of the equations residual $\| b - A x^n \|$  (or an approximation to that norm), they are never a norm of the error in the equation.

271:   Not all combinations of preconditioner side (see `KSPSetPCSide()`) and norm types are supported by all Krylov methods.
272:   If only one is set, PETSc tries to automatically change the other to find a compatible pair.  If no such combination
273:   is supported, PETSc will generate an error.

275:   Developer Note:
276:   Supported combinations of norm and preconditioner side are set using `KSPSetSupportedNorm()` for each `KSPType`.

278: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetCheckNormIteration()`, `KSPSetPCSide()`, `KSPGetPCSide()`, `KSPNormType`
279: @*/
280: PetscErrorCode KSPSetNormType(KSP ksp, KSPNormType normtype)
281: {
282:   PetscFunctionBegin;
285:   ksp->normtype = ksp->normtype_set = normtype;
286:   PetscFunctionReturn(PETSC_SUCCESS);
287: }

289: /*@
290:   KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
291:   computed and used in the convergence test of `KSPSolve()` for the given `KSP` context

293:   Logically Collective

295:   Input Parameters:
296: + ksp - Krylov solver context
297: - it  - use -1 to check at all iterations

299:   Level: advanced

301:   Notes:
302:   Currently only works with `KSPCG`, `KSPBCGS` and `KSPIBCGS`

304:   Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm

306:   On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
307:   `-ksp_monitor` the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).

309:   Certain methods such as `KSPGMRES` always compute the residual norm, this routine will not change that computation, but it will
310:   prevent the computed norm from being checked.

312: .seealso: [](ch_ksp), `KSP`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetLagNorm()`
313: @*/
314: PetscErrorCode KSPSetCheckNormIteration(KSP ksp, PetscInt it)
315: {
316:   PetscFunctionBegin;
319:   ksp->chknorm = it;
320:   PetscFunctionReturn(PETSC_SUCCESS);
321: }

323: /*@
324:   KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the `MPI_Allreduce()` used for
325:   computing the inner products needed for the next iteration.

327:   Logically Collective

329:   Input Parameters:
330: + ksp - Krylov solver context
331: - flg - `PETSC_TRUE` or `PETSC_FALSE`

333:   Options Database Key:
334: . -ksp_lag_norm - lag the calculated residual norm

336:   Level: advanced

338:   Notes:
339:   Currently only works with `KSPIBCGS`.

341:   This can reduce communication costs at the expense of doing
342:   one additional iteration because the norm used in the convergence test of `KSPSolve()` is one iteration behind the actual
343:   current residual norm (which has not yet been computed due to the lag).

345:   Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm

347:   If you lag the norm and run with, for example, `-ksp_monitor`, the residual norm reported will be the lagged one.

349:   `KSPSetCheckNormIteration()` is an alternative way of avoiding the expense of computing the residual norm at each iteration.

351: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetCheckNormIteration()`
352: @*/
353: PetscErrorCode KSPSetLagNorm(KSP ksp, PetscBool flg)
354: {
355:   PetscFunctionBegin;
358:   ksp->lagnorm = flg;
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: /*@
363:   KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a `KSPType`

365:   Logically Collective

367:   Input Parameters:
368: + ksp      - Krylov method
369: . normtype - supported norm type of the type `KSPNormType`
370: . pcside   - preconditioner side, of the type `PCSide` that can be used with this `KSPNormType`
371: - priority - positive integer preference for this combination; larger values have higher priority

373:   Level: developer

375:   Notes:
376:   This function should be called from the implementation files `KSPCreate_XXX()` to declare
377:   which norms and preconditioner sides are supported. Users should not call this
378:   function.

380:   This function can be called multiple times for each combination of `KSPNormType` and `PCSide`
381:   the `KSPType` supports

383: .seealso: [](ch_ksp), `KSP`, `KSPNormType`, `PCSide`, `KSPSetNormType()`, `KSPSetPCSide()`
384: @*/
385: PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType normtype, PCSide pcside, PetscInt priority)
386: {
387:   PetscFunctionBegin;
389:   ksp->normsupporttable[normtype][pcside] = priority;
390:   PetscFunctionReturn(PETSC_SUCCESS);
391: }

393: static PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
394: {
395:   PetscFunctionBegin;
396:   PetscCall(PetscMemzero(ksp->normsupporttable, sizeof(ksp->normsupporttable)));
397:   ksp->pc_side  = ksp->pc_side_set;
398:   ksp->normtype = ksp->normtype_set;
399:   PetscFunctionReturn(PETSC_SUCCESS);
400: }

402: PetscErrorCode KSPSetUpNorms_Private(KSP ksp, PetscBool errorifnotsupported, KSPNormType *normtype, PCSide *pcside)
403: {
404:   PetscInt i, j, best, ibest = 0, jbest = 0;

406:   PetscFunctionBegin;
407:   best = 0;
408:   for (i = 0; i < KSP_NORM_MAX; i++) {
409:     for (j = 0; j < PC_SIDE_MAX; j++) {
410:       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && ksp->normsupporttable[i][j] > best) {
411:         best  = ksp->normsupporttable[i][j];
412:         ibest = i;
413:         jbest = j;
414:       }
415:     }
416:   }
417:   if (best < 1 && errorifnotsupported) {
418:     PetscCheck(ksp->normtype != KSP_NORM_DEFAULT || ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "The %s KSP implementation did not call KSPSetSupportedNorm()", ((PetscObject)ksp)->type_name);
419:     PetscCheck(ksp->normtype != KSP_NORM_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support preconditioner side %s", ((PetscObject)ksp)->type_name, PCSides[ksp->pc_side]);
420:     PetscCheck(ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype]);
421:     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s with preconditioner side %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype], PCSides[ksp->pc_side]);
422:   }
423:   if (normtype) *normtype = (KSPNormType)ibest;
424:   if (pcside) *pcside = (PCSide)jbest;
425:   PetscFunctionReturn(PETSC_SUCCESS);
426: }

428: /*@
429:   KSPGetNormType - Gets the `KSPNormType` that is used for convergence testing during `KSPSolve()` for this `KSP` context

431:   Not Collective

433:   Input Parameter:
434: . ksp - Krylov solver context

436:   Output Parameter:
437: . normtype - the `KSPNormType` that is used for convergence testing

439:   Level: advanced

441: .seealso: [](ch_ksp), `KSPNormType`, `KSPSetNormType()`, `KSPConvergedSkip()`
442: @*/
443: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
444: {
445:   PetscFunctionBegin;
447:   PetscAssertPointer(normtype, 2);
448:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
449:   *normtype = ksp->normtype;
450:   PetscFunctionReturn(PETSC_SUCCESS);
451: }

453: #if defined(PETSC_HAVE_SAWS)
454: #include <petscviewersaws.h>
455: #endif

457: /*@
458:   KSPSetOperators - Sets the matrix associated with the linear system
459:   and a (possibly) different one from which the preconditioner will be built into the `KSP` context. The matrix will then be used during `KSPSolve()`

461:   Collective

463:   Input Parameters:
464: + ksp  - the `KSP` context
465: . Amat - the matrix that defines the linear system
466: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`.

468:   Level: beginner

470:   Notes:
471: .vb
472:   KSPSetOperators(ksp, Amat, Pmat);
473: .ve
474:   is the same as
475: .vb
476:   KSPGetPC(ksp, &pc);
477:   PCSetOperators(pc, Amat, Pmat);
478: .ve
479:   and is equivalent to
480: .vb
481:   PCCreate(PetscObjectComm((PetscObject)ksp), &pc);
482:   PCSetOperators(pc, Amat, Pmat);
483:   KSPSetPC(ksp, pc);
484: .ve

486:   If you know the operator `Amat` has a null space you can use `MatSetNullSpace()` and `MatSetTransposeNullSpace()` to supply the null
487:   space to `Amat` and the `KSP` solvers will automatically use that null space as needed during the solution process.

489:   All future calls to `KSPSetOperators()` must use the same size matrices, unless `KSPReset()` is called!

491:   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently being used from the `KSP` context.

493:   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
494:   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
495:   on it and then pass it back in your call to `KSPSetOperators()`.

497:   Developer Notes:
498:   If the operators have NOT been set with `KSPSetOperators()` then the operators
499:   are created in the `PC` and returned to the user. In this case, if both operators
500:   mat and pmat are requested, two DIFFERENT operators will be returned. If
501:   only one is requested both operators in the `PC` will be the same (i.e. as
502:   if one had called `KSPSetOperators()` with the same argument for both `Mat`s).
503:   The user must set the sizes of the returned matrices and their type etc just
504:   as if the user created them with `MatCreate()`. For example,

506: .vb
507:          KSPGetOperators(ksp/pc,&mat,NULL); is equivalent to
508:            set size, type, etc of mat

510:          MatCreate(comm,&mat);
511:          KSP/PCSetOperators(ksp/pc,mat,mat);
512:          PetscObjectDereference((PetscObject)mat);
513:            set size, type, etc of mat

515:      and

517:          KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
518:            set size, type, etc of mat and pmat

520:          MatCreate(comm,&mat);
521:          MatCreate(comm,&pmat);
522:          KSP/PCSetOperators(ksp/pc,mat,pmat);
523:          PetscObjectDereference((PetscObject)mat);
524:          PetscObjectDereference((PetscObject)pmat);
525:            set size, type, etc of mat and pmat
526: .ve

528:   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
529:   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
530:   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
531:   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
532:   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a `KSP`
533:   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
534:   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
535:   it can be created for you?

537: .seealso: [](ch_ksp), `KSP`, `Mat`, `KSPSolve()`, `KSPGetPC()`, `PCGetOperators()`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPSetComputeRHS()`
538: @*/
539: PetscErrorCode KSPSetOperators(KSP ksp, Mat Amat, Mat Pmat)
540: {
541:   PetscFunctionBegin;
545:   if (Amat) PetscCheckSameComm(ksp, 1, Amat, 2);
546:   if (Pmat) PetscCheckSameComm(ksp, 1, Pmat, 3);
547:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
548:   PetscCall(PCSetOperators(ksp->pc, Amat, Pmat));
549:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
550:   PetscFunctionReturn(PETSC_SUCCESS);
551: }

553: /*@
554:   KSPGetOperators - Gets the matrix associated with the linear system
555:   and a (possibly) different one used to construct the preconditioner from the `KSP` context

557:   Collective

559:   Input Parameter:
560: . ksp - the `KSP` context

562:   Output Parameters:
563: + Amat - the matrix that defines the linear system
564: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`.

566:   Level: intermediate

568:   Notes:
569:   If `KSPSetOperators()` has not been called then the `KSP` object will attempt to automatically create the matrix `Amat` and return it

571:   Use `KSPGetOperatorsSet()` to determine if matrices have been provided.

573:   DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.

575: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetPC()`, `PCSetOperators()`, `KSPSetOperators()`, `KSPGetOperatorsSet()`
576: @*/
577: PetscErrorCode KSPGetOperators(KSP ksp, Mat *Amat, Mat *Pmat)
578: {
579:   PetscFunctionBegin;
581:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
582:   PetscCall(PCGetOperators(ksp->pc, Amat, Pmat));
583:   PetscFunctionReturn(PETSC_SUCCESS);
584: }

586: /*@
587:   KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
588:   possibly a different one from which the preconditioner will be built have been set in the `KSP` with `KSPSetOperators()`

590:   Not Collective, though the results on all processes will be the same

592:   Input Parameter:
593: . ksp - the `KSP` context

595:   Output Parameters:
596: + mat  - the matrix associated with the linear system was set
597: - pmat - matrix from which the preconditioner will be built, usually the same as `mat` was set

599:   Level: intermediate

601:   Note:
602:   This routine exists because if you call `KSPGetOperators()` on a `KSP` that does not yet have operators they are
603:   automatically created in the call.

605: .seealso: [](ch_ksp), `KSP`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`, `PCGetOperatorsSet()`
606: @*/
607: PetscErrorCode KSPGetOperatorsSet(KSP ksp, PetscBool *mat, PetscBool *pmat)
608: {
609:   PetscFunctionBegin;
611:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
612:   PetscCall(PCGetOperatorsSet(ksp->pc, mat, pmat));
613:   PetscFunctionReturn(PETSC_SUCCESS);
614: }

616: /*@C
617:   KSPSetPreSolve - Sets a function that is called at the beginning of each `KSPSolve()`. Used in conjunction with `KSPSetPostSolve()`.

619:   Logically Collective

621:   Input Parameters:
622: + ksp      - the solver object
623: . presolve - the function to call before the solve, see` KSPPSolveFn`
624: - ctx      - an optional context needed by the function

626:   Level: developer

628:   Notes:
629:   The function provided here `presolve` is used to modify the right hand side, and possibly the matrix, of the linear system to be solved.
630:   The function provided with `KSPSetPostSolve()` then modifies the resulting solution of that linear system to obtain the correct solution
631:   to the initial linear system.

633:   The functions `PCPreSolve()` and `PCPostSolve()` provide a similar functionality and are used, for example with `PCEISENSTAT`.

635: .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPostSolve()`, `PCEISENSTAT`, `PCPreSolve()`, `PCPostSolve()`
636: @*/
637: PetscErrorCode KSPSetPreSolve(KSP ksp, KSPPSolveFn *presolve, PetscCtx ctx)
638: {
639:   PetscFunctionBegin;
641:   ksp->presolve = presolve;
642:   ksp->prectx   = ctx;
643:   PetscFunctionReturn(PETSC_SUCCESS);
644: }

646: /*@C
647:   KSPSetPostSolve - Sets a function that is called at the end of each `KSPSolve()` (whether it converges or not). Used in conjunction with `KSPSetPreSolve()`.

649:   Logically Collective

651:   Input Parameters:
652: + ksp       - the solver object
653: . postsolve - the function to call after the solve, see` KSPPSolveFn`
654: - ctx       - an optional context needed by the function

656:   Level: developer

658: .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPreSolve()`, `PCEISENSTAT`
659: @*/
660: PetscErrorCode KSPSetPostSolve(KSP ksp, KSPPSolveFn *postsolve, PetscCtx ctx)
661: {
662:   PetscFunctionBegin;
664:   ksp->postsolve = postsolve;
665:   ksp->postctx   = ctx;
666:   PetscFunctionReturn(PETSC_SUCCESS);
667: }

669: /*@
670:   KSPSetNestLevel - sets the amount of nesting the `KSP` has. That is the number of levels of `KSP` above this `KSP` in a linear solve.

672:   Collective

674:   Input Parameters:
675: + ksp   - the `KSP`
676: - level - the nest level

678:   Level: developer

680:   Note:
681:   For example, the `KSP` in each block of a `KSPBJACOBI` has a level of 1, while the outer `KSP` has a level of 0.

683: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
684: @*/
685: PetscErrorCode KSPSetNestLevel(KSP ksp, PetscInt level)
686: {
687:   PetscFunctionBegin;
690:   ksp->nestlevel = level;
691:   PetscFunctionReturn(PETSC_SUCCESS);
692: }

694: /*@
695:   KSPGetNestLevel - gets the amount of nesting the `KSP` has

697:   Not Collective

699:   Input Parameter:
700: . ksp - the `KSP`

702:   Output Parameter:
703: . level - the nest level

705:   Level: developer

707: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
708: @*/
709: PetscErrorCode KSPGetNestLevel(KSP ksp, PetscInt *level)
710: {
711:   PetscFunctionBegin;
713:   PetscAssertPointer(level, 2);
714:   *level = ksp->nestlevel;
715:   PetscFunctionReturn(PETSC_SUCCESS);
716: }

718: /*@
719:   KSPCreate - Creates the `KSP` context. This `KSP` context is used in PETSc to solve linear systems with `KSPSolve()`

721:   Collective

723:   Input Parameter:
724: . comm - MPI communicator

726:   Output Parameter:
727: . inksp - location to put the `KSP` context

729:   Level: beginner

731:   Note:
732:   The default `KSPType` is `KSPGMRES` with a restart of 30, using modified Gram-Schmidt orthogonalization. The `KSPType` may be
733:   changed with `KSPSetType()`

735: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetType()`
736: @*/
737: PetscErrorCode KSPCreate(MPI_Comm comm, KSP *inksp)
738: {
739:   KSP      ksp;
740:   PetscCtx ctx;

742:   PetscFunctionBegin;
743:   PetscAssertPointer(inksp, 2);
744:   PetscCall(KSPInitializePackage());

746:   PetscCall(PetscHeaderCreate(ksp, KSP_CLASSID, "KSP", "Krylov Method", "KSP", comm, KSPDestroy, KSPView));
747:   ksp->default_max_it = ksp->max_it = 10000;
748:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;

750:   ksp->default_rtol = ksp->rtol = 1.e-5;
751:   ksp->default_abstol = ksp->abstol = PetscDefined(USE_REAL_SINGLE) ? 1.e-25 : 1.e-50;
752:   ksp->default_divtol = ksp->divtol = 1.e4;

754:   ksp->chknorm  = -1;
755:   ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
756:   ksp->rnorm                        = 0.0;
757:   ksp->its                          = 0;
758:   ksp->guess_zero                   = PETSC_TRUE;
759:   ksp->calc_sings                   = PETSC_FALSE;
760:   ksp->res_hist                     = NULL;
761:   ksp->res_hist_alloc               = NULL;
762:   ksp->res_hist_len                 = 0;
763:   ksp->res_hist_max                 = 0;
764:   ksp->res_hist_reset               = PETSC_TRUE;
765:   ksp->err_hist                     = NULL;
766:   ksp->err_hist_alloc               = NULL;
767:   ksp->err_hist_len                 = 0;
768:   ksp->err_hist_max                 = 0;
769:   ksp->err_hist_reset               = PETSC_TRUE;
770:   ksp->numbermonitors               = 0;
771:   ksp->numberreasonviews            = 0;
772:   ksp->setfromoptionscalled         = 0;
773:   ksp->nmax                         = PETSC_DECIDE;

775:   PetscCall(KSPConvergedDefaultCreate(&ctx));
776:   PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy));
777:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
778:   ksp->ops->buildresidual = KSPBuildResidualDefault;

780:   ksp->vec_sol    = NULL;
781:   ksp->vec_rhs    = NULL;
782:   ksp->pc         = NULL;
783:   ksp->data       = NULL;
784:   ksp->nwork      = 0;
785:   ksp->work       = NULL;
786:   ksp->reason     = KSP_CONVERGED_ITERATING;
787:   ksp->setupstage = KSP_SETUP_NEW;

789:   PetscCall(KSPNormSupportTableReset_Private(ksp));

791:   *inksp = ksp;
792:   PetscFunctionReturn(PETSC_SUCCESS);
793: }

795: /*@
796:   KSPSetType - Sets the algorithm/method to be used to solve the linear system with the given `KSP`

798:   Logically Collective

800:   Input Parameters:
801: + ksp  - the Krylov space context
802: - type - a known method

804:   Options Database Key:
805: . -ksp_type type - Sets the method; see `KSPType`

807:   Level: intermediate

809:   Notes:
810:   See `KSPType` for available methods (for instance, `KSPCG` or `KSPGMRES`).

812:   Normally, it is best to use the `KSPSetFromOptions()` command and
813:   then set the `KSP` type from the options database rather than by using
814:   this routine.  Using the options database provides the user with
815:   maximum flexibility in evaluating the many different Krylov methods.
816:   The `KSPSetType()` routine is provided for those situations where it
817:   is necessary to set the iterative solver independently of the command
818:   line or options database.  This might be the case, for example, when
819:   the choice of iterative solver changes during the execution of the
820:   program, and the user's application is taking responsibility for
821:   choosing the appropriate method.  In other words, this routine is
822:   not for beginners.

824:   Developer Note:
825:   `KSPRegister()` is used to add Krylov types to `KSPList` from which they are accessed by `KSPSetType()`.

827: .seealso: [](ch_ksp), `PCSetType()`, `KSPType`, `KSPRegister()`, `KSPCreate()`, `KSP`
828: @*/
829: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
830: {
831:   PetscBool match;
832:   PetscErrorCode (*r)(KSP);

834:   PetscFunctionBegin;
836:   PetscAssertPointer(type, 2);

838:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, type, &match));
839:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

841:   PetscCall(PetscFunctionListFind(KSPList, type, &r));
842:   PetscCheck(r, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested KSP type %s", type);
843:   /* Destroy the previous private KSP context */
844:   PetscTryTypeMethod(ksp, destroy);

846:   /* Reinitialize function pointers in KSPOps structure */
847:   PetscCall(PetscMemzero(ksp->ops, sizeof(struct _KSPOps)));
848:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
849:   ksp->ops->buildresidual = KSPBuildResidualDefault;
850:   PetscCall(KSPNormSupportTableReset_Private(ksp));
851:   ksp->converged_neg_curve = PETSC_FALSE; // restore default
852:   ksp->setupnewmatrix      = PETSC_FALSE; // restore default (setup not called in case of new matrix)
853:   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
854:   ksp->setupstage     = KSP_SETUP_NEW;
855:   ksp->guess_not_read = PETSC_FALSE; // restore default
856:   PetscCall((*r)(ksp));
857:   PetscCall(PetscObjectChangeTypeName((PetscObject)ksp, type));
858:   PetscFunctionReturn(PETSC_SUCCESS);
859: }

861: /*@
862:   KSPGetType - Gets the `KSP` type as a string from the `KSP` object.

864:   Not Collective

866:   Input Parameter:
867: . ksp - Krylov context

869:   Output Parameter:
870: . type - name of the `KSP` method

872:   Level: intermediate

874: .seealso: [](ch_ksp), `KSPType`, `KSP`, `KSPSetType()`
875: @*/
876: PetscErrorCode KSPGetType(KSP ksp, KSPType *type)
877: {
878:   PetscFunctionBegin;
880:   PetscAssertPointer(type, 2);
881:   *type = ((PetscObject)ksp)->type_name;
882:   PetscFunctionReturn(PETSC_SUCCESS);
883: }

885: /*@C
886:   KSPRegister -  Adds a method, `KSPType`, to the Krylov subspace solver package.

888:   Not Collective, No Fortran Support

890:   Input Parameters:
891: + sname    - name of a new user-defined solver
892: - function - routine to create method

894:   Level: advanced

896:   Note:
897:   `KSPRegister()` may be called multiple times to add several user-defined solvers.

899:   Example Usage:
900: .vb
901:    KSPRegister("my_solver", MySolverCreate);
902: .ve

904:   Then, your solver can be chosen with the procedural interface via
905: .vb
906:   KSPSetType(ksp, "my_solver")
907: .ve
908:   or at runtime via the option `-ksp_type my_solver`

910: .seealso: [](ch_ksp), `KSP`, `KSPType`, `KSPSetType`, `KSPRegisterAll()`
911: @*/
912: PetscErrorCode KSPRegister(const char sname[], PetscErrorCode (*function)(KSP))
913: {
914:   PetscFunctionBegin;
915:   PetscCall(KSPInitializePackage());
916:   PetscCall(PetscFunctionListAdd(&KSPList, sname, function));
917:   PetscFunctionReturn(PETSC_SUCCESS);
918: }

920: PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
921: {
922:   PetscFunctionBegin;
923:   PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN));
924:   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
925:   PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN));
926:   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
927:   PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN));
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: /*@C
932:   KSPMonitorRegister -  Registers a Krylov subspace solver monitor routine that may be accessed with `KSPMonitorSetFromOptions()`

934:   Not Collective

936:   Input Parameters:
937: + name    - name of a new monitor type
938: . vtype   - A `PetscViewerType` for the output
939: . format  - A `PetscViewerFormat` for the output
940: . monitor - Monitor routine, see `KSPMonitorRegisterFn`
941: . create  - Creation routine, or `NULL`
942: - destroy - Destruction routine, or `NULL`

944:   Level: advanced

946:   Notes:
947:   `KSPMonitorRegister()` may be called multiple times to add several user-defined monitors.

949:   The calling sequence for the given function matches the calling sequence used by `KSPMonitorFn` functions passed to `KSPMonitorSet()` with the additional
950:   requirement that its final argument be a `PetscViewerAndFormat`.

952:   Example Usage:
953: .vb
954:   KSPMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
955: .ve

957:   Then, your monitor can be chosen with the procedural interface via
958: .vb
959:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_my_monitor", "my_monitor", NULL)
960: .ve
961:   or at runtime via the option `-ksp_monitor_my_monitor`

963: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegisterAll()`, `KSPMonitorSetFromOptions()`
964: @*/
965: PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format, KSPMonitorRegisterFn *monitor, KSPMonitorRegisterCreateFn *create, KSPMonitorRegisterDestroyFn *destroy)
966: {
967:   char key[PETSC_MAX_PATH_LEN];

969:   PetscFunctionBegin;
970:   PetscCall(KSPInitializePackage());
971:   PetscCall(KSPMonitorMakeKey_Internal(name, vtype, format, key));
972:   PetscCall(PetscFunctionListAdd(&KSPMonitorList, key, monitor));
973:   if (create) PetscCall(PetscFunctionListAdd(&KSPMonitorCreateList, key, create));
974:   if (destroy) PetscCall(PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }