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: /*@C
 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: /*@C
 71:   KSPView - Prints the `KSP` data structure.

 73:   Collective

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

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

 82:   Level: beginner

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

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

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

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

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

110:   PetscFunctionBegin;
112:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp), &viewer));
114:   PetscCheckSameComm(ksp, 1, viewer, 2);

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

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

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

192:     PetscCall(PetscObjectGetName((PetscObject)ksp, &name));
193:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
194:     if (!((PetscObject)ksp)->amsmem && rank == 0) {
195:       char dir[1024];

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

215: /*@C
216:   KSPViewFromOptions - View a `KSP` object based on values in the options database

218:   Collective

220:   Input Parameters:
221: + A    - Krylov solver context
222: . obj  - Optional object
223: - name - command line option

225:   Level: intermediate

227: .seealso: [](ch_ksp), `KSP`, `KSPView`, `PetscObjectViewFromOptions()`, `KSPCreate()`
228: @*/
229: PetscErrorCode KSPViewFromOptions(KSP A, PetscObject obj, const char name[])
230: {
231:   PetscFunctionBegin;
233:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
234:   PetscFunctionReturn(PETSC_SUCCESS);
235: }

237: /*@
238:   KSPSetNormType - Sets the norm that is used for convergence testing.

240:   Logically Collective

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

258:   Options Database Key:
259: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural> - set `KSP` norm type

261:   Level: advanced

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

268:   Developer Note:
269:   Supported combinations of norm and preconditioner side are set using `KSPSetSupportedNorm()`.

271: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetCheckNormIteration()`, `KSPSetPCSide()`, `KSPGetPCSide()`, `KSPNormType`
272: @*/
273: PetscErrorCode KSPSetNormType(KSP ksp, KSPNormType normtype)
274: {
275:   PetscFunctionBegin;
278:   ksp->normtype = ksp->normtype_set = normtype;
279:   PetscFunctionReturn(PETSC_SUCCESS);
280: }

282: /*@
283:   KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
284:   computed and used in the convergence test.

286:   Logically Collective

288:   Input Parameters:
289: + ksp - Krylov solver context
290: - it  - use -1 to check at all iterations

292:   Notes:
293:   Currently only works with `KSPCG`, `KSPBCGS` and `KSPIBCGS`

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

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

300:   Level: advanced

302: .seealso: [](ch_ksp), `KSP`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`
303: @*/
304: PetscErrorCode KSPSetCheckNormIteration(KSP ksp, PetscInt it)
305: {
306:   PetscFunctionBegin;
309:   ksp->chknorm = it;
310:   PetscFunctionReturn(PETSC_SUCCESS);
311: }

313: /*@
314:   KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the `MPI_Allreduce()` for
315:   computing the inner products for the next iteration.  This can reduce communication costs at the expense of doing
316:   one additional iteration.

318:   Logically Collective

320:   Input Parameters:
321: + ksp - Krylov solver context
322: - flg - `PETSC_TRUE` or `PETSC_FALSE`

324:   Options Database Key:
325: . -ksp_lag_norm - lag the calculated residual norm

327:   Level: advanced

329:   Notes:
330:   Currently only works with `KSPIBCGS`.

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

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

336: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetCheckNormIteration()`
337: @*/
338: PetscErrorCode KSPSetLagNorm(KSP ksp, PetscBool flg)
339: {
340:   PetscFunctionBegin;
343:   ksp->lagnorm = flg;
344:   PetscFunctionReturn(PETSC_SUCCESS);
345: }

347: /*@
348:   KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a `KSP`

350:   Logically Collective

352:   Input Parameters:
353: + ksp      - Krylov method
354: . normtype - supported norm type
355: . pcside   - preconditioner side that can be used with this norm
356: - priority - positive integer preference for this combination; larger values have higher priority

358:   Level: developer

360:   Note:
361:   This function should be called from the implementation files `KSPCreate_XXX()` to declare
362:   which norms and preconditioner sides are supported. Users should not need to call this
363:   function.

365: .seealso: [](ch_ksp), `KSP`, `KSPNormType`, `PCSide`, `KSPSetNormType()`, `KSPSetPCSide()`
366: @*/
367: PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType normtype, PCSide pcside, PetscInt priority)
368: {
369:   PetscFunctionBegin;
371:   ksp->normsupporttable[normtype][pcside] = priority;
372:   PetscFunctionReturn(PETSC_SUCCESS);
373: }

375: static PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
376: {
377:   PetscFunctionBegin;
378:   PetscCall(PetscMemzero(ksp->normsupporttable, sizeof(ksp->normsupporttable)));
379:   ksp->pc_side  = ksp->pc_side_set;
380:   ksp->normtype = ksp->normtype_set;
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

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

388:   PetscFunctionBegin;
389:   best = 0;
390:   for (i = 0; i < KSP_NORM_MAX; i++) {
391:     for (j = 0; j < PC_SIDE_MAX; j++) {
392:       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
393:         best  = ksp->normsupporttable[i][j];
394:         ibest = i;
395:         jbest = j;
396:       }
397:     }
398:   }
399:   if (best < 1 && errorifnotsupported) {
400:     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);
401:     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]);
402:     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]);
403:     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]);
404:   }
405:   if (normtype) *normtype = (KSPNormType)ibest;
406:   if (pcside) *pcside = (PCSide)jbest;
407:   PetscFunctionReturn(PETSC_SUCCESS);
408: }

410: /*@
411:   KSPGetNormType - Gets the norm that is used for convergence testing.

413:   Not Collective

415:   Input Parameter:
416: . ksp - Krylov solver context

418:   Output Parameter:
419: . normtype - norm that is used for convergence testing

421:   Level: advanced

423: .seealso: [](ch_ksp), `KSPNormType`, `KSPSetNormType()`, `KSPConvergedSkip()`
424: @*/
425: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
426: {
427:   PetscFunctionBegin;
429:   PetscAssertPointer(normtype, 2);
430:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
431:   *normtype = ksp->normtype;
432:   PetscFunctionReturn(PETSC_SUCCESS);
433: }

435: #if defined(PETSC_HAVE_SAWS)
436: #include <petscviewersaws.h>
437: #endif

439: /*@
440:   KSPSetOperators - Sets the matrix associated with the linear system
441:   and a (possibly) different one from which the preconditioner will be built

443:   Collective

445:   Input Parameters:
446: + ksp  - the `KSP` context
447: . Amat - the matrix that defines the linear system
448: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

450:   Level: beginner

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

456:   All future calls to `KSPSetOperators()` must use the same size matrices!

458:   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.

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

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

473: .vb
474:          KSPGetOperators(ksp/pc,&mat,NULL); is equivalent to
475:            set size, type, etc of mat

477:          MatCreate(comm,&mat);
478:          KSP/PCSetOperators(ksp/pc,mat,mat);
479:          PetscObjectDereference((PetscObject)mat);
480:            set size, type, etc of mat

482:      and

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

487:          MatCreate(comm,&mat);
488:          MatCreate(comm,&pmat);
489:          KSP/PCSetOperators(ksp/pc,mat,pmat);
490:          PetscObjectDereference((PetscObject)mat);
491:          PetscObjectDereference((PetscObject)pmat);
492:            set size, type, etc of mat and pmat
493: .ve

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

504: .seealso: [](ch_ksp), `KSP`, `Mat`, `KSPSolve()`, `KSPGetPC()`, `PCGetOperators()`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPSetComputeRHS()`
505: @*/
506: PetscErrorCode KSPSetOperators(KSP ksp, Mat Amat, Mat Pmat)
507: {
508:   PetscFunctionBegin;
512:   if (Amat) PetscCheckSameComm(ksp, 1, Amat, 2);
513:   if (Pmat) PetscCheckSameComm(ksp, 1, Pmat, 3);
514:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
515:   PetscCall(PCSetOperators(ksp->pc, Amat, Pmat));
516:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: /*@
521:   KSPGetOperators - Gets the matrix associated with the linear system
522:   and a (possibly) different one used to construct the preconditioner.

524:   Collective

526:   Input Parameter:
527: . ksp - the `KSP` context

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

533:   Level: intermediate

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

538: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetPC()`, `PCGetOperators()`, `PCSetOperators()`, `KSPSetOperators()`, `KSPGetOperatorsSet()`
539: @*/
540: PetscErrorCode KSPGetOperators(KSP ksp, Mat *Amat, Mat *Pmat)
541: {
542:   PetscFunctionBegin;
544:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
545:   PetscCall(PCGetOperators(ksp->pc, Amat, Pmat));
546:   PetscFunctionReturn(PETSC_SUCCESS);
547: }

549: /*@C
550:   KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
551:   possibly a different one associated with the preconditioner have been set in the `KSP`.

553:   Not Collective, though the results on all processes should be the same

555:   Input Parameter:
556: . ksp - the `KSP` context

558:   Output Parameters:
559: + mat  - the matrix associated with the linear system was set
560: - pmat - matrix associated with the preconditioner was set, usually the same as `mat`

562:   Level: intermediate

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

568: .seealso: [](ch_ksp), `KSP`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`, `PCGetOperatorsSet()`
569: @*/
570: PetscErrorCode KSPGetOperatorsSet(KSP ksp, PetscBool *mat, PetscBool *pmat)
571: {
572:   PetscFunctionBegin;
574:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
575:   PetscCall(PCGetOperatorsSet(ksp->pc, mat, pmat));
576:   PetscFunctionReturn(PETSC_SUCCESS);
577: }

579: /*@C
580:   KSPSetPreSolve - Sets a function that is called at the beginning of each `KSPSolve()`

582:   Logically Collective

584:   Input Parameters:
585: + ksp      - the solver object
586: . presolve - the function to call before the solve
587: - ctx      - any context needed by the function

589:   Calling sequence of `presolve`:
590: + ksp - the `KSP` context
591: . rhs - the right-hand side vector
592: . x   - the solution vector
593: - ctx - optional user-provided context

595:   Level: developer

597: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPostSolve()`, `PCEISENSTAT`
598: @*/
599: PetscErrorCode KSPSetPreSolve(KSP ksp, PetscErrorCode (*presolve)(KSP ksp, Vec rhs, Vec x, void *ctx), void *ctx)
600: {
601:   PetscFunctionBegin;
603:   ksp->presolve = presolve;
604:   ksp->prectx   = ctx;
605:   PetscFunctionReturn(PETSC_SUCCESS);
606: }

608: /*@C
609:   KSPSetPostSolve - Sets a function that is called at the end of each `KSPSolve()` (whether it converges or not)

611:   Logically Collective

613:   Input Parameters:
614: + ksp       - the solver object
615: . postsolve - the function to call after the solve
616: - ctx       - any context needed by the function

618:   Calling sequence of `postsolve`:
619: + ksp - the `KSP` context
620: . rhs - the right-hand side vector
621: . x   - the solution vector
622: - ctx - optional user-provided context

624:   Level: developer

626: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPreSolve()`, `PCEISENSTAT`
627: @*/
628: PetscErrorCode KSPSetPostSolve(KSP ksp, PetscErrorCode (*postsolve)(KSP ksp, Vec rhs, Vec x, void *ctx), void *ctx)
629: {
630:   PetscFunctionBegin;
632:   ksp->postsolve = postsolve;
633:   ksp->postctx   = ctx;
634:   PetscFunctionReturn(PETSC_SUCCESS);
635: }

637: /*@
638:   KSPSetNestLevel - sets the amount of nesting the `KSP` has

640:   Collective

642:   Input Parameters:
643: + ksp   - the `KSP`
644: - level - the nest level

646:   Level: developer

648: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
649: @*/
650: PetscErrorCode KSPSetNestLevel(KSP ksp, PetscInt level)
651: {
652:   PetscFunctionBegin;
655:   ksp->nestlevel = level;
656:   PetscFunctionReturn(PETSC_SUCCESS);
657: }

659: /*@
660:   KSPGetNestLevel - gets the amount of nesting the `KSP` has

662:   Not Collective

664:   Input Parameter:
665: . ksp - the `KSP`

667:   Output Parameter:
668: . level - the nest level

670:   Level: developer

672: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
673: @*/
674: PetscErrorCode KSPGetNestLevel(KSP ksp, PetscInt *level)
675: {
676:   PetscFunctionBegin;
678:   PetscAssertPointer(level, 2);
679:   *level = ksp->nestlevel;
680:   PetscFunctionReturn(PETSC_SUCCESS);
681: }

683: /*@
684:   KSPCreate - Creates the `KSP` context.

686:   Collective

688:   Input Parameter:
689: . comm - MPI communicator

691:   Output Parameter:
692: . inksp - location to put the `KSP` context

694:   Level: beginner

696:   Note:
697:   The default `KSPType` is `KSPGMRES` with a restart of 30, using modified Gram-Schmidt orthogonalization.

699: .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`
700: @*/
701: PetscErrorCode KSPCreate(MPI_Comm comm, KSP *inksp)
702: {
703:   KSP   ksp;
704:   void *ctx;

706:   PetscFunctionBegin;
707:   PetscAssertPointer(inksp, 2);
708:   *inksp = NULL;
709:   PetscCall(KSPInitializePackage());

711:   PetscCall(PetscHeaderCreate(ksp, KSP_CLASSID, "KSP", "Krylov Method", "KSP", comm, KSPDestroy, KSPView));

713:   ksp->max_it  = 10000;
714:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
715:   ksp->rtol                       = 1.e-5;
716: #if defined(PETSC_USE_REAL_SINGLE)
717:   ksp->abstol = 1.e-25;
718: #else
719:   ksp->abstol = 1.e-50;
720: #endif
721:   ksp->divtol = 1.e4;

723:   ksp->chknorm  = -1;
724:   ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
725:   ksp->rnorm                        = 0.0;
726:   ksp->its                          = 0;
727:   ksp->guess_zero                   = PETSC_TRUE;
728:   ksp->calc_sings                   = PETSC_FALSE;
729:   ksp->res_hist                     = NULL;
730:   ksp->res_hist_alloc               = NULL;
731:   ksp->res_hist_len                 = 0;
732:   ksp->res_hist_max                 = 0;
733:   ksp->res_hist_reset               = PETSC_TRUE;
734:   ksp->err_hist                     = NULL;
735:   ksp->err_hist_alloc               = NULL;
736:   ksp->err_hist_len                 = 0;
737:   ksp->err_hist_max                 = 0;
738:   ksp->err_hist_reset               = PETSC_TRUE;
739:   ksp->numbermonitors               = 0;
740:   ksp->numberreasonviews            = 0;
741:   ksp->setfromoptionscalled         = 0;
742:   ksp->nmax                         = PETSC_DECIDE;

744:   PetscCall(KSPConvergedDefaultCreate(&ctx));
745:   PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy));
746:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
747:   ksp->ops->buildresidual = KSPBuildResidualDefault;

749:   ksp->vec_sol    = NULL;
750:   ksp->vec_rhs    = NULL;
751:   ksp->pc         = NULL;
752:   ksp->data       = NULL;
753:   ksp->nwork      = 0;
754:   ksp->work       = NULL;
755:   ksp->reason     = KSP_CONVERGED_ITERATING;
756:   ksp->setupstage = KSP_SETUP_NEW;

758:   PetscCall(KSPNormSupportTableReset_Private(ksp));

760:   *inksp = ksp;
761:   PetscFunctionReturn(PETSC_SUCCESS);
762: }

764: /*@C
765:   KSPSetType - Builds the `KSP` data structure for a particular `KSPType`

767:   Logically Collective

769:   Input Parameters:
770: + ksp  - the Krylov space context
771: - type - a known method

773:   Options Database Key:
774: . -ksp_type  <method> - Sets the method; use `-help` for a list  of available methods (for instance, cg or gmres)

776:   Level: intermediate

778:   Notes:
779:   See "petsc/include/petscksp.h" for available methods (for instance, `KSPCG` or `KSPGMRES`).

781:   Normally, it is best to use the `KSPSetFromOptions()` command and
782:   then set the `KSP` type from the options database rather than by using
783:   this routine.  Using the options database provides the user with
784:   maximum flexibility in evaluating the many different Krylov methods.
785:   The `KSPSetType()` routine is provided for those situations where it
786:   is necessary to set the iterative solver independently of the command
787:   line or options database.  This might be the case, for example, when
788:   the choice of iterative solver changes during the execution of the
789:   program, and the user's application is taking responsibility for
790:   choosing the appropriate method.  In other words, this routine is
791:   not for beginners.

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

796: .seealso: [](ch_ksp), `PCSetType()`, `KSPType`, `KSPRegister()`, `KSPCreate()`, `KSP`
797: @*/
798: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
799: {
800:   PetscBool match;
801:   PetscErrorCode (*r)(KSP);

803:   PetscFunctionBegin;
805:   PetscAssertPointer(type, 2);

807:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, type, &match));
808:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

810:   PetscCall(PetscFunctionListFind(KSPList, type, &r));
811:   PetscCheck(r, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested KSP type %s", type);
812:   /* Destroy the previous private KSP context */
813:   PetscTryTypeMethod(ksp, destroy);
814:   ksp->ops->destroy = NULL;

816:   /* Reinitialize function pointers in KSPOps structure */
817:   PetscCall(PetscMemzero(ksp->ops, sizeof(struct _KSPOps)));
818:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
819:   ksp->ops->buildresidual = KSPBuildResidualDefault;
820:   PetscCall(KSPNormSupportTableReset_Private(ksp));
821:   ksp->converged_neg_curve = PETSC_FALSE; // restore default
822:   ksp->setupnewmatrix      = PETSC_FALSE; // restore default (setup not called in case of new matrix)
823:   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
824:   ksp->setupstage     = KSP_SETUP_NEW;
825:   ksp->guess_not_read = PETSC_FALSE; // restore default
826:   PetscCall((*r)(ksp));
827:   PetscCall(PetscObjectChangeTypeName((PetscObject)ksp, type));
828:   PetscFunctionReturn(PETSC_SUCCESS);
829: }

831: /*@C
832:   KSPGetType - Gets the `KSP` type as a string from the `KSP` object.

834:   Not Collective

836:   Input Parameter:
837: . ksp - Krylov context

839:   Output Parameter:
840: . type - name of the `KSP` method

842:   Level: intermediate

844: .seealso: [](ch_ksp), `KSPType`, `KSP`, `KSPSetType()`
845: @*/
846: PetscErrorCode KSPGetType(KSP ksp, KSPType *type)
847: {
848:   PetscFunctionBegin;
850:   PetscAssertPointer(type, 2);
851:   *type = ((PetscObject)ksp)->type_name;
852:   PetscFunctionReturn(PETSC_SUCCESS);
853: }

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

858:   Not Collective

860:   Input Parameters:
861: + sname    - name of a new user-defined solver
862: - function - routine to create method

864:   Level: advanced

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

869:   Example Usage:
870: .vb
871:    KSPRegister("my_solver", MySolverCreate);
872: .ve

874:   Then, your solver can be chosen with the procedural interface via
875: .vb
876:   KSPSetType(ksp, "my_solver")
877: .ve
878:   or at runtime via the option `-ksp_type my_solver`

880: .seealso: [](ch_ksp), `KSP`, `KSPType`, `KSPSetType`, `KSPRegisterAll()`
881: @*/
882: PetscErrorCode KSPRegister(const char sname[], PetscErrorCode (*function)(KSP))
883: {
884:   PetscFunctionBegin;
885:   PetscCall(KSPInitializePackage());
886:   PetscCall(PetscFunctionListAdd(&KSPList, sname, function));
887:   PetscFunctionReturn(PETSC_SUCCESS);
888: }

890: PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
891: {
892:   PetscFunctionBegin;
893:   PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN));
894:   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
895:   PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN));
896:   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
897:   PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN));
898:   PetscFunctionReturn(PETSC_SUCCESS);
899: }

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

904:   Not Collective

906:   Input Parameters:
907: + name    - name of a new monitor routine
908: . vtype   - A `PetscViewerType` for the output
909: . format  - A `PetscViewerFormat` for the output
910: . monitor - Monitor routine
911: . create  - Creation routine, or `NULL`
912: - destroy - Destruction routine, or `NULL`

914:   Level: advanced

916:   Note:
917:   `KSPMonitorRegister()` may be called multiple times to add several user-defined monitors.

919:   Example Usage:
920: .vb
921:   KSPMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
922: .ve

924:   Then, your monitor can be chosen with the procedural interface via
925: .vb
926:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_my_monitor", "my_monitor", NULL)
927: .ve
928:   or at runtime via the option `-ksp_monitor_my_monitor`

930: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegisterAll()`, `KSPMonitorSetFromOptions()`
931: @*/
932: PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format, PetscErrorCode (*monitor)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*create)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*destroy)(PetscViewerAndFormat **))
933: {
934:   char key[PETSC_MAX_PATH_LEN];

936:   PetscFunctionBegin;
937:   PetscCall(KSPInitializePackage());
938:   PetscCall(KSPMonitorMakeKey_Internal(name, vtype, format, key));
939:   PetscCall(PetscFunctionListAdd(&KSPMonitorList, key, monitor));
940:   if (create) PetscCall(PetscFunctionListAdd(&KSPMonitorCreateList, key, create));
941:   if (destroy) PetscCall(PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy));
942:   PetscFunctionReturn(PETSC_SUCCESS);
943: }