Actual source code: itcl.c

  1: /*
  2:     Code for setting KSP options from the options database.
  3: */

  5: #include <petsc/private/kspimpl.h>
  6: #include <petscdraw.h>

  8: /*@
  9:   KSPSetOptionsPrefix - Sets the prefix used for searching for all
 10:   `KSP` options in the database.

 12:   Logically Collective

 14:   Input Parameters:
 15: + ksp    - the Krylov context
 16: - prefix - the prefix string to prepend to all `KSP` option requests

 18:   Level: intermediate

 20:   Notes:
 21:   A hyphen (-) must NOT be given at the beginning of the prefix name.
 22:   The first character of all runtime options is AUTOMATICALLY the
 23:   hyphen.

 25:   For example, to distinguish between the runtime options for two
 26:   different `KSP` contexts, one could call
 27: .vb
 28:       KSPSetOptionsPrefix(ksp1,"sys1_")
 29:       KSPSetOptionsPrefix(ksp2,"sys2_")
 30: .ve

 32:   This would enable use of different options for each system, such as
 33: .vb
 34:       -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3
 35:       -sys2_ksp_type bcgs  -sys2_ksp_rtol 1.e-4
 36: .ve

 38: .seealso: [](ch_ksp), `KSP`, `KSPAppendOptionsPrefix()`, `KSPGetOptionsPrefix()`, `KSPSetFromOptions()`
 39: @*/
 40: PetscErrorCode KSPSetOptionsPrefix(KSP ksp, const char prefix[])
 41: {
 42:   PetscBool ispcmpi;

 44:   PetscFunctionBegin;
 46:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
 47:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &ispcmpi));
 48:   if (ispcmpi) {
 49:     size_t     len;
 50:     const char suffix[] = "mpi_linear_solver_server_";
 51:     char      *newprefix;

 53:     PetscCall(PetscStrlen(prefix, &len));
 54:     PetscCall(PetscMalloc1(len + sizeof(suffix) + 1, &newprefix));
 55:     PetscCall(PetscStrncpy(newprefix, prefix, len + sizeof(suffix)));
 56:     PetscCall(PetscStrlcat(newprefix, suffix, len + sizeof(suffix)));
 57:     PetscCall(PCSetOptionsPrefix(ksp->pc, newprefix));
 58:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ksp, newprefix));
 59:     PetscCall(PetscFree(newprefix));
 60:   } else {
 61:     PetscCall(PCSetOptionsPrefix(ksp->pc, prefix));
 62:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ksp, prefix));
 63:   }
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*@
 68:   KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
 69:   `KSP` options in the database.

 71:   Logically Collective

 73:   Input Parameters:
 74: + ksp    - the Krylov context
 75: - prefix - the prefix string to prepend to all `KSP` option requests

 77:   Level: intermediate

 79:   Note:
 80:   A hyphen (-) must NOT be given at the beginning of the prefix name.
 81:   The first character of all runtime options is AUTOMATICALLY the hyphen.

 83: .seealso: [](ch_ksp), `KSP`, `KSPSetOptionsPrefix()`, `KSPGetOptionsPrefix()`, `KSPSetFromOptions()`
 84: @*/
 85: PetscErrorCode KSPAppendOptionsPrefix(KSP ksp, const char prefix[])
 86: {
 87:   PetscFunctionBegin;
 89:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
 90:   PetscCall(PCAppendOptionsPrefix(ksp->pc, prefix));
 91:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)ksp, prefix));
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: /*@
 96:   KSPSetUseFischerGuess - Use the Paul Fischer algorithm or its variants to compute initial guesses for a set of solves with related right-hand sides

 98:   Logically Collective

100:   Input Parameters:
101: + ksp   - the Krylov context
102: . model - use model 1, model 2, model 3, or any other number to turn it off
103: - size  - size of subspace used to generate initial guess

105:   Options Database Key:
106: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves

108:   Level: advanced

110: .seealso: [](ch_ksp), `KSP`, `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetGuess()`, `KSPGetGuess()`, `KSPGuess`
111: @*/
112: PetscErrorCode KSPSetUseFischerGuess(KSP ksp, PetscInt model, PetscInt size)
113: {
114:   KSPGuess guess;

116:   PetscFunctionBegin;
120:   PetscCall(KSPGetGuess(ksp, &guess));
121:   PetscCall(KSPGuessSetType(guess, KSPGUESSFISCHER));
122:   PetscCall(KSPGuessFischerSetModel(guess, model, size));
123:   PetscFunctionReturn(PETSC_SUCCESS);
124: }

126: /*@
127:   KSPSetGuess - Set the initial guess object `KSPGuess` to be used by the `KSP` object to generate initial guesses

129:   Logically Collective

131:   Input Parameters:
132: + ksp   - the Krylov context
133: - guess - the object created with `KSPGuessCreate()`

135:   Level: advanced

137:   Notes:
138:   this allows a single `KSP` to be used with several different initial guess generators (likely for different linear
139:   solvers, see `KSPSetPC()`).

141:   This increases the reference count of the guess object, you must destroy the object with `KSPGuessDestroy()`
142:   before the end of the program.

144: .seealso: [](ch_ksp), `KSP`, `KSPGuess`, `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetUseFischerGuess()`, `KSPGetGuess()`
145: @*/
146: PetscErrorCode KSPSetGuess(KSP ksp, KSPGuess guess)
147: {
148:   PetscFunctionBegin;
151:   PetscCall(PetscObjectReference((PetscObject)guess));
152:   PetscCall(KSPGuessDestroy(&ksp->guess));
153:   ksp->guess      = guess;
154:   ksp->guess->ksp = ksp;
155:   PetscFunctionReturn(PETSC_SUCCESS);
156: }

158: /*@
159:   KSPGetGuess - Gets the initial guess generator for the `KSP`.

161:   Not Collective

163:   Input Parameter:
164: . ksp - the Krylov context

166:   Output Parameter:
167: . guess - the object

169:   Level: developer

171: .seealso: [](ch_ksp), `KSPGuess`, `KSP`, `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetUseFischerGuess()`, `KSPSetGuess()`
172: @*/
173: PetscErrorCode KSPGetGuess(KSP ksp, KSPGuess *guess)
174: {
175:   PetscFunctionBegin;
177:   PetscAssertPointer(guess, 2);
178:   if (!ksp->guess) {
179:     const char *prefix;

181:     PetscCall(KSPGuessCreate(PetscObjectComm((PetscObject)ksp), &ksp->guess));
182:     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
183:     if (prefix) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ksp->guess, prefix));
184:     ksp->guess->ksp = ksp;
185:   }
186:   *guess = ksp->guess;
187:   PetscFunctionReturn(PETSC_SUCCESS);
188: }

190: /*@
191:   KSPGetOptionsPrefix - Gets the prefix used for searching for all
192:   `KSP` options in the database.

194:   Not Collective

196:   Input Parameter:
197: . ksp - the Krylov context

199:   Output Parameter:
200: . prefix - pointer to the prefix string used is returned

202:   Level: advanced

204:   Fortran Note:
205:   Pass in a string 'prefix' of
206:   sufficient length to hold the prefix.

208: .seealso: [](ch_ksp), `KSP`, `KSPSetFromOptions()`, `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`
209: @*/
210: PetscErrorCode KSPGetOptionsPrefix(KSP ksp, const char *prefix[])
211: {
212:   PetscFunctionBegin;
214:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, prefix));
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
219: {
220:   PetscFunctionBegin;
221:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
222:   (*vf)->data = ctx;
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: /*@C
227:   KSPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user in the options database

229:   Collective

231:   Input Parameters:
232: + ksp  - `KSP` object you wish to monitor
233: . opt  - the command line option for this monitor
234: . name - the monitor type one is seeking
235: - ctx  - An optional user context for the monitor, or `NULL`

237:   Level: developer

239: .seealso: [](ch_ksp), `KSPMonitorRegister()`, `KSPMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
240:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
241:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
242:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
243:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
244:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
245:           `PetscOptionsFList()`, `PetscOptionsEList()`
246: @*/
247: PetscErrorCode KSPMonitorSetFromOptions(KSP ksp, const char opt[], const char name[], void *ctx)
248: {
249:   PetscErrorCode (*mfunc)(KSP, PetscInt, PetscReal, void *);
250:   PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
251:   PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
252:   PetscViewerAndFormat *vf;
253:   PetscViewer           viewer;
254:   PetscViewerFormat     format;
255:   PetscViewerType       vtype;
256:   char                  key[PETSC_MAX_PATH_LEN];
257:   PetscBool             all, flg;
258:   const char           *prefix = NULL;

260:   PetscFunctionBegin;
261:   PetscCall(PetscStrcmp(opt, "-all_ksp_monitor", &all));
262:   if (!all) PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
263:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->options, prefix, opt, &viewer, &format, &flg));
264:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

266:   PetscCall(PetscViewerGetType(viewer, &vtype));
267:   PetscCall(KSPMonitorMakeKey_Internal(name, vtype, format, key));
268:   PetscCall(PetscFunctionListFind(KSPMonitorList, key, &mfunc));
269:   PetscCall(PetscFunctionListFind(KSPMonitorCreateList, key, &cfunc));
270:   PetscCall(PetscFunctionListFind(KSPMonitorDestroyList, key, &dfunc));
271:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
272:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

274:   PetscCall((*cfunc)(viewer, format, ctx, &vf));
275:   PetscCall(PetscViewerDestroy(&viewer));
276:   PetscCall(KSPMonitorSet(ksp, mfunc, vf, (PetscErrorCode (*)(void **))dfunc));
277:   PetscFunctionReturn(PETSC_SUCCESS);
278: }

280: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP);

282: /*@
283:   KSPSetFromOptions - Sets `KSP` options from the options database.
284:   This routine must be called before `KSPSetUp()` if the user is to be
285:   allowed to set the Krylov type.

287:   Collective

289:   Input Parameter:
290: . ksp - the Krylov space context

292:   Options Database Keys:
293: + -ksp_rtol rtol                                                          - relative tolerance used in default determination of convergence, i.e.
294:                                                                             if residual norm decreases by this factor than convergence is declared
295: . -ksp_atol abstol                                                        - absolute tolerance used in default convergence test, i.e. if residual
296:                                                                             norm is less than this then convergence is declared
297: . -ksp_divtol tol                                                         - if residual norm increases by this factor than divergence is declared
298: . -ksp_max_it                                                             - maximum number of linear iterations
299: . -ksp_min_it                                                             - minimum number of linear iterations to use, defaults to zero

301: . -ksp_reuse_preconditioner <true,false>                                  - reuse the previously computed preconditioner

303: . -ksp_converged_use_initial_residual_norm                                - see `KSPConvergedDefaultSetUIRNorm()`
304: . -ksp_converged_use_min_initial_residual_norm                            - see `KSPConvergedDefaultSetUMIRNorm()`
305: . -ksp_converged_maxits                                                   - see `KSPConvergedDefaultSetConvergedMaxits()`
306: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>           - see `KSPSetNormType()`
307: . -ksp_check_norm_iteration it                                            - do not compute residual norm until iteration number it (does compute at 0th iteration)
308:                                                                             works only for `KSPBCGS`, `KSPIBCGS`, and `KSPCG`
309: . -ksp_lag_norm                                                           - compute the norm of the residual for the ith iteration on the i+1 iteration;
310:                                                                             this means that one can use the norm of the residual for convergence test WITHOUT
311:                                                                             an extra `MPI_Allreduce()` limiting global synchronizations.
312:                                                                             This will require 1 more iteration of the solver than usual.
313: . -ksp_guess_type                                                         - Type of initial guess generator for repeated linear solves
314: . -ksp_fischer_guess <model,size>                                         - uses the Fischer initial guess generator for repeated linear solves
315: . -ksp_constant_null_space                                                - assume the operator (matrix) has the constant vector in its null space
316: . -ksp_test_null_space                                                    - tests the null space set with `MatSetNullSpace()` to see if it truly is a null space
317: . -ksp_knoll                                                              - compute initial guess by applying the preconditioner to the right-hand side
318: . -ksp_monitor_cancel                                                     - cancel all previous convergene monitor routines set
319: . -ksp_monitor                                                            - print residual norm at each iteration
320: . -ksp_monitor draw::draw_lg                                              - plot residual norm at each iteration, see `KSPMonitorResidual()`
321: . -ksp_monitor_true_residual                                              - print the true l2 residual norm at each iteration, see `KSPMonitorTrueResidual()`
322: . -all_ksp_monitor <optional filename>                                    - print residual norm at each iteration for ALL KSP solves, regardless of their prefix. This is
323:                                                                             useful for `PCFIELDSPLIT`, `PCMG`, etc that have inner solvers and
324:                                                                             you wish to track the convergence of all the solvers
325: . -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration
326: . -ksp_monitor_singular_value                                             - monitor extreme singular values at each iteration
327: . -ksp_converged_reason                                                   - view the convergence state at the end of the solve
328: . -ksp_use_explicittranspose                                              - transpose the system explicitly in `KSPSolveTranspose()`
329: . -ksp_error_if_not_converged                                             - stop the program as soon as an error is detected in a `KSPSolve()`, `KSP_DIVERGED_ITS`
330:                                                                             is not treated as an error on inner solves
331: - -ksp_converged_rate                                                     - view the convergence rate at the end of the solve

333:   Level: beginner

335:   Note:
336:   To see all options, run your program with the `-help` option or consult [](ch_ksp)

338: .seealso: [](ch_ksp), `KSP`, `KSPSetOptionsPrefix()`, `KSPResetFromOptions()`, `KSPSetUseFischerGuess()`
339: @*/
340: PetscErrorCode KSPSetFromOptions(KSP ksp)
341: {
342:   const char *convtests[] = {"default", "skip", "lsqr"}, *prefix;
343:   char        type[256], guesstype[256], monfilename[PETSC_MAX_PATH_LEN];
344:   PetscBool   flg, flag, reuse, set;
345:   PetscInt    indx, model[2] = {0, 0}, nmax, max_it;
346:   KSPNormType normtype;
347:   PCSide      pcside;
348:   void       *ctx;
349:   MPI_Comm    comm;
350:   PetscReal   rtol, abstol, divtol;

352:   PetscFunctionBegin;

355:   PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
356:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
357:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
358:   if (!ksp->skippcsetfromoptions) PetscCall(PCSetFromOptions(ksp->pc));

360:   PetscCall(KSPRegisterAll());
361:   PetscObjectOptionsBegin((PetscObject)ksp);
362:   PetscCall(PetscOptionsFList("-ksp_type", "Krylov method", "KSPSetType", KSPList, (char *)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES), type, 256, &flg));
363:   if (flg) PetscCall(KSPSetType(ksp, type));
364:   /*
365:     Set the type if it was never set.
366:   */
367:   if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp, KSPGMRES));

369:   PetscCall(KSPResetViewers(ksp));

371:   /* Cancels all monitors hardwired into code before call to KSPSetFromOptions() */
372:   PetscCall(PetscOptionsBool("-ksp_monitor_cancel", "Remove any hardwired monitor routines", "KSPMonitorCancel", PETSC_FALSE, &flg, &set));
373:   if (set && flg) PetscCall(KSPMonitorCancel(ksp));
374:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor", "preconditioned_residual", NULL));
375:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_short", "preconditioned_residual_short", NULL));
376:   PetscCall(KSPMonitorSetFromOptions(ksp, "-all_ksp_monitor", "preconditioned_residual", NULL));
377:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_range", "preconditioned_residual_range", NULL));
378:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_true_residual", "true_residual", NULL));
379:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_max", "true_residual_max", NULL));
380:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_solution", "solution", NULL));
381:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_singular_value", "singular_value", ksp));
382:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_error", "error", ksp));
383:   PetscCall(PetscOptionsBool("-ksp_monitor_pause_final", "Pauses all draw monitors at the final iterate", "KSPMonitorPauseFinal_Internal", PETSC_FALSE, &ksp->pauseFinal, NULL));
384:   PetscCall(PetscOptionsBool("-ksp_initial_guess_nonzero", "Use the contents of the solution vector for initial guess", "KSPSetInitialNonzero", ksp->guess_zero ? PETSC_FALSE : PETSC_TRUE, &flag, &flg));
385:   if (flg) PetscCall(KSPSetInitialGuessNonzero(ksp, flag));

387:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg));
388:   if (flg) {
389:     PetscCall(KSPGetReusePreconditioner(ksp, &reuse));
390:     PetscCall(PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL));
391:     PetscCall(KSPSetReusePreconditioner(ksp, reuse));
392:     PetscCall(PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL));
393:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view));
394:     PetscCall(PetscViewerDestroy(&ksp->convergedreasonviewer));
395:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_converged_reason", &ksp->convergedreasonviewer, &ksp->convergedreasonformat, NULL));
396:     flg = PETSC_FALSE;
397:     PetscCall(PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set));
398:     if (set && flg) PetscCall(KSPConvergedReasonViewCancel(ksp));
399:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat));
400:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat));
401:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs));
402:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol));
403:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp));
404:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes));
405:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp));
406:     PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale));

408:     PetscCall(KSPGetDiagonalScale(ksp, &flag));
409:     PetscCall(PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg));
410:     if (flg) PetscCall(KSPSetDiagonalScale(ksp, flag));
411:     PetscCall(KSPGetDiagonalScaleFix(ksp, &flag));
412:     PetscCall(PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg));
413:     if (flg) PetscCall(KSPSetDiagonalScaleFix(ksp, flag));
414:     nmax = ksp->nmax;
415:     PetscCall(PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL));
416:     PetscCall(PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg));
417:     if (flg) PetscCall(KSPSetMatSolveBatchSize(ksp, nmax));
418:     goto skipoptions;
419:   }

421:   rtol   = ksp->rtol;
422:   abstol = ksp->abstol;
423:   divtol = ksp->divtol;
424:   max_it = ksp->max_it;
425:   PetscCall(PetscOptionsReal("-ksp_rtol", "Relative decrease in residual norm", "KSPSetTolerances", ksp->rtol, &rtol, NULL));
426:   PetscCall(PetscOptionsReal("-ksp_atol", "Absolute value of residual norm", "KSPSetTolerances", ksp->abstol, &abstol, NULL));
427:   PetscCall(PetscOptionsReal("-ksp_divtol", "Residual norm increase cause divergence", "KSPSetTolerances", ksp->divtol, &divtol, NULL));
428:   PetscCall(PetscOptionsInt("-ksp_max_it", "Maximum number of iterations", "KSPSetTolerances", ksp->max_it, &max_it, &flg));
429:   PetscCall(KSPSetTolerances(ksp, rtol, abstol, divtol, max_it));
430:   PetscCall(PetscOptionsRangeInt("-ksp_min_it", "Minimum number of iterations", "KSPSetMinimumIterations", ksp->min_it, &ksp->min_it, NULL, 0, ksp->max_it));

432:   PetscCall(PetscOptionsBool("-ksp_converged_use_initial_residual_norm", "Use initial residual norm for computing relative convergence", "KSPConvergedDefaultSetUIRNorm", PETSC_FALSE, &flag, &set));
433:   if (set && flag) PetscCall(KSPConvergedDefaultSetUIRNorm(ksp));
434:   PetscCall(PetscOptionsBool("-ksp_converged_use_min_initial_residual_norm", "Use minimum of initial residual norm and b for computing relative convergence", "KSPConvergedDefaultSetUMIRNorm", PETSC_FALSE, &flag, &set));
435:   if (set && flag) PetscCall(KSPConvergedDefaultSetUMIRNorm(ksp));
436:   PetscCall(PetscOptionsBool("-ksp_converged_maxits", "Declare convergence if the maximum number of iterations is reached", "KSPConvergedDefaultSetConvergedMaxits", PETSC_FALSE, &flag, &set));
437:   if (set) PetscCall(KSPConvergedDefaultSetConvergedMaxits(ksp, flag));
438:   PetscCall(KSPGetConvergedNegativeCurvature(ksp, &flag));
439:   PetscCall(PetscOptionsBool("-ksp_converged_neg_curve", "Declare convergence if negative curvature is detected", "KSPConvergedNegativeCurvature", flag, &flag, &set));
440:   if (set) PetscCall(KSPSetConvergedNegativeCurvature(ksp, flag));
441:   PetscCall(KSPGetReusePreconditioner(ksp, &reuse));
442:   PetscCall(PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL));
443:   PetscCall(KSPSetReusePreconditioner(ksp, reuse));

445:   PetscCall(PetscOptionsBool("-ksp_knoll", "Use preconditioner applied to b for initial guess", "KSPSetInitialGuessKnoll", ksp->guess_knoll, &ksp->guess_knoll, NULL));
446:   PetscCall(PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL));
447:   PetscCall(PetscOptionsFList("-ksp_guess_type", "Initial guess in Krylov method", NULL, KSPGuessList, NULL, guesstype, 256, &flg));
448:   if (flg) {
449:     PetscCall(KSPGetGuess(ksp, &ksp->guess));
450:     PetscCall(KSPGuessSetType(ksp->guess, guesstype));
451:     PetscCall(KSPGuessSetFromOptions(ksp->guess));
452:   } else { /* old option for KSP */
453:     nmax = 2;
454:     PetscCall(PetscOptionsIntArray("-ksp_fischer_guess", "Use Paul Fischer's algorithm or its variants for initial guess", "KSPSetUseFischerGuess", model, &nmax, &flag));
455:     if (flag) {
456:       PetscCheck(nmax == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in model,size as arguments");
457:       PetscCall(KSPSetUseFischerGuess(ksp, model[0], model[1]));
458:     }
459:   }

461:   PetscCall(PetscOptionsEList("-ksp_convergence_test", "Convergence test", "KSPSetConvergenceTest", convtests, 3, "default", &indx, &flg));
462:   if (flg) {
463:     switch (indx) {
464:     case 0:
465:       PetscCall(KSPConvergedDefaultCreate(&ctx));
466:       PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy));
467:       break;
468:     case 1:
469:       PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedSkip, NULL, NULL));
470:       break;
471:     case 2:
472:       PetscCall(KSPConvergedDefaultCreate(&ctx));
473:       PetscCall(KSPSetConvergenceTest(ksp, KSPLSQRConvergedDefault, ctx, KSPConvergedDefaultDestroy));
474:       break;
475:     }
476:   }

478:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_FALSE, &normtype, NULL));
479:   PetscCall(PetscOptionsEnum("-ksp_norm_type", "KSP Norm type", "KSPSetNormType", KSPNormTypes, (PetscEnum)normtype, (PetscEnum *)&normtype, &flg));
480:   if (flg) PetscCall(KSPSetNormType(ksp, normtype));

482:   PetscCall(PetscOptionsInt("-ksp_check_norm_iteration", "First iteration to compute residual norm", "KSPSetCheckNormIteration", ksp->chknorm, &ksp->chknorm, NULL));

484:   PetscCall(PetscOptionsBool("-ksp_lag_norm", "Lag the calculation of the residual norm", "KSPSetLagNorm", ksp->lagnorm, &flag, &flg));
485:   if (flg) PetscCall(KSPSetLagNorm(ksp, flag));

487:   PetscCall(KSPGetDiagonalScale(ksp, &flag));
488:   PetscCall(PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg));
489:   if (flg) PetscCall(KSPSetDiagonalScale(ksp, flag));
490:   PetscCall(KSPGetDiagonalScaleFix(ksp, &flag));
491:   PetscCall(PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg));
492:   if (flg) PetscCall(KSPSetDiagonalScaleFix(ksp, flag));

494:   PetscCall(PetscOptionsBool("-ksp_constant_null_space", "Add constant null space to Krylov solver matrix", "MatSetNullSpace", PETSC_FALSE, &flg, &set));
495:   if (set && flg) {
496:     MatNullSpace nsp;
497:     Mat          Amat = NULL;

499:     PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, &nsp));
500:     if (ksp->pc) PetscCall(PCGetOperators(ksp->pc, &Amat, NULL));
501:     PetscCheck(Amat, comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot set nullspace, matrix has not yet been provided");
502:     PetscCall(MatSetNullSpace(Amat, nsp));
503:     PetscCall(MatNullSpaceDestroy(&nsp));
504:   }

506:   flg = PETSC_FALSE;
507:   if (ksp->pc) {
508:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCKSP, &flg));
509:     if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCBJACOBI, &flg));
510:     if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCDEFLATION, &flg));
511:   }

513:   if (flg) {
514:     /* Using dynamic tolerance in preconditioner */
515:     PetscCall(PetscOptionsString("-sub_ksp_dynamic_tolerance", "Use dynamic tolerance for inner PC", "KSPMonitorDynamicTolerance", "stdout", monfilename, sizeof(monfilename), &flg));
516:     if (flg) {
517:       void     *scale;
518:       PetscReal coeff = 1.0;

520:       PetscCall(KSPMonitorDynamicToleranceCreate(&scale));
521:       PetscCall(PetscOptionsReal("-sub_ksp_dynamic_tolerance", "Coefficient of dynamic tolerance for inner PC", "KSPMonitorDynamicTolerance", coeff, &coeff, &flg));
522:       if (flg) PetscCall(KSPMonitorDynamicToleranceSetCoefficient(scale, coeff));
523:       PetscCall(KSPMonitorSet(ksp, KSPMonitorDynamicTolerance, scale, KSPMonitorDynamicToleranceDestroy));
524:     }
525:   }

527:   /*
528:    Calls Python function
529:   */
530:   PetscCall(PetscOptionsString("-ksp_monitor_python", "Use Python function", "KSPMonitorSet", NULL, monfilename, sizeof(monfilename), &flg));
531:   if (flg) PetscCall(PetscPythonMonitorSet((PetscObject)ksp, monfilename));
532:   /*
533:     Graphically plots preconditioned residual norm and range of residual element values
534:   */
535:   PetscCall(PetscOptionsBool("-ksp_monitor_lg_range", "Monitor graphically range of preconditioned residual norm", "KSPMonitorSet", PETSC_FALSE, &flg, &set));
536:   if (set && flg) {
537:     PetscViewer ctx;

539:     PetscCall(PetscViewerDrawOpen(comm, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &ctx));
540:     PetscCall(KSPMonitorSet(ksp, KSPMonitorLGRange, ctx, (PetscErrorCode (*)(void **))PetscViewerDestroy));
541:   }
542:   /* TODO Do these show up in help? */
543:   PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &flg));
544:   if (flg) {
545:     const char *RateTypes[] = {"default", "residual", "error", "PetscRateType", "RATE_", NULL};
546:     PetscEnum   rtype       = (PetscEnum)1;

548:     PetscCall(PetscOptionsGetEnum(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate_type", RateTypes, &rtype, &flg));
549:     if (rtype == (PetscEnum)0 || rtype == (PetscEnum)1) PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE));
550:     if (rtype == (PetscEnum)0 || rtype == (PetscEnum)2) PetscCall(KSPSetErrorHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE));
551:   }
552:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view));
553:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pre", &ksp->viewerPre, &ksp->formatPre, &ksp->viewPre));

555:   PetscCall(PetscViewerDestroy(&ksp->convergedreasonviewer));
556:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_converged_reason", &ksp->convergedreasonviewer, &ksp->convergedreasonformat, NULL));
557:   flg = PETSC_FALSE;
558:   PetscCall(PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set));
559:   if (set && flg) PetscCall(KSPConvergedReasonViewCancel(ksp));
560:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &ksp->viewerRate, &ksp->formatRate, &ksp->viewRate));
561:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat));
562:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat));
563:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs));
564:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol));
565:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp));
566:   PetscCall(PetscOptionsDeprecated("-ksp_compute_eigenvalues", "-ksp_view_eigenvalues", "3.9", NULL));
567:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues", &ksp->viewerEV, &ksp->formatEV, &ksp->viewEV));
568:   PetscCall(PetscOptionsDeprecated("-ksp_compute_singularvalues", "-ksp_view_singularvalues", "3.9", NULL));
569:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_singularvalues", &ksp->viewerSV, &ksp->formatSV, &ksp->viewSV));
570:   PetscCall(PetscOptionsDeprecated("-ksp_compute_eigenvalues_explicitly", "-ksp_view_eigenvalues_explicit", "3.9", NULL));
571:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues_explicit", &ksp->viewerEVExp, &ksp->formatEVExp, &ksp->viewEVExp));
572:   PetscCall(PetscOptionsDeprecated("-ksp_final_residual", "-ksp_view_final_residual", "3.9", NULL));
573:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes));
574:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp));
575:   PetscCall(PetscOptionsCreateViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale));

577:   /* Deprecated options */
578:   if (!ksp->viewEV) {
579:     /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigenvalues",...) below because the argument handling is different */
580:     PetscCall(PetscOptionsDeprecated("-ksp_plot_eigenvalues", NULL, "3.9", "Use -ksp_view_eigenvalues draw"));
581:     PetscCall(PetscOptionsName("-ksp_plot_eigenvalues", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw]", "KSPView", &ksp->viewEV));
582:     if (ksp->viewEV) {
583:       ksp->formatEV = PETSC_VIEWER_DEFAULT;
584:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
585:       PetscCall(PetscObjectReference((PetscObject)ksp->viewerEV));
586:     }
587:   }
588:   if (!ksp->viewEV) {
589:     PetscCall(PetscOptionsDeprecated("-ksp_plot_eigenvalues_explicitly", NULL, "3.9", "Use -ksp_view_eigenvalues_explicit draw"));
590:     /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigencontours",...) below because the argument handling is different */
591:     PetscCall(PetscOptionsName("-ksp_plot_eigencontours", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw::draw_contour]", "KSPView", &ksp->viewEV));
592:     if (ksp->viewEV) {
593:       ksp->formatEV = PETSC_VIEWER_DRAW_CONTOUR;
594:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
595:       PetscCall(PetscObjectReference((PetscObject)ksp->viewerEV));
596:     }
597:   }
598:   if (!ksp->viewEVExp) {
599:     /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigencontours_explicitly",...) below because the argument handling is different */
600:     PetscCall(PetscOptionsName("-ksp_plot_eigenvalues_explicitly", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues_explicit draw]", "KSPView", &ksp->viewEVExp));
601:     if (ksp->viewEVExp) {
602:       ksp->formatEVExp = PETSC_VIEWER_DEFAULT;
603:       ksp->viewerEVExp = PETSC_VIEWER_DRAW_(comm);
604:       PetscCall(PetscObjectReference((PetscObject)ksp->viewerEVExp));
605:     }
606:   }

608: #if defined(PETSC_HAVE_SAWS)
609:   /*
610:     Publish convergence information using AMS
611:   */
612:   PetscCall(PetscOptionsBool("-ksp_monitor_saws", "Publish KSP progress using SAWs", "KSPMonitorSet", PETSC_FALSE, &flg, &set));
613:   if (set && flg) {
614:     void *ctx;
615:     PetscCall(KSPMonitorSAWsCreate(ksp, &ctx));
616:     PetscCall(KSPMonitorSet(ksp, KSPMonitorSAWs, ctx, KSPMonitorSAWsDestroy));
617:     PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
618:   }
619: #endif

621:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_FALSE, NULL, &pcside));
622:   PetscCall(PetscOptionsEnum("-ksp_pc_side", "KSP preconditioner side", "KSPSetPCSide", PCSides, (PetscEnum)pcside, (PetscEnum *)&pcside, &flg));
623:   if (flg) PetscCall(KSPSetPCSide(ksp, pcside));

625:   if (ksp->viewSV || ksp->viewEV) PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));

627: #if defined(PETSC_HAVE_SAWS)
628:   {
629:     PetscBool set;
630:     flg = PETSC_FALSE;
631:     PetscCall(PetscOptionsBool("-ksp_saws_block", "Block for SAWs at end of KSPSolve", "PetscObjectSAWsBlock", ((PetscObject)ksp)->amspublishblock, &flg, &set));
632:     if (set) PetscCall(PetscObjectSAWsSetBlock((PetscObject)ksp, flg));
633:   }
634: #endif

636:   nmax = ksp->nmax;
637:   PetscCall(PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL));
638:   PetscCall(PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg));
639:   if (flg) PetscCall(KSPSetMatSolveBatchSize(ksp, nmax));

641:   flg = PETSC_FALSE;
642:   PetscCall(PetscOptionsBool("-ksp_use_explicittranspose", "Explicitly transpose the system in KSPSolveTranspose", "KSPSetUseExplicitTranspose", ksp->transpose.use_explicittranspose, &flg, &set));
643:   if (set) PetscCall(KSPSetUseExplicitTranspose(ksp, flg));

645:   PetscTryTypeMethod(ksp, setfromoptions, PetscOptionsObject);
646: skipoptions:
647:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
648:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)ksp, PetscOptionsObject));
649:   PetscOptionsEnd();
650:   ksp->setfromoptionscalled++;
651:   PetscFunctionReturn(PETSC_SUCCESS);
652: }

654: /*@
655:   KSPResetFromOptions - Sets `KSP` parameters from user options ONLY if the `KSP` was previously set from options

657:   Collective

659:   Input Parameter:
660: . ksp - the `KSP` context

662:   Level: advanced

664: .seealso: [](ch_ksp), `KSPSetFromOptions()`, `KSPSetOptionsPrefix()`
665: @*/
666: PetscErrorCode KSPResetFromOptions(KSP ksp)
667: {
668:   PetscFunctionBegin;
669:   if (ksp->setfromoptionscalled) PetscCall(KSPSetFromOptions(ksp));
670:   PetscFunctionReturn(PETSC_SUCCESS);
671: }