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: /*@C
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: /*@C
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
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: /*@C
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()`, `PetscOptionsGetViewer()`, `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(PetscOptionsGetViewer(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(PetscObjectDereference((PetscObject)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_max_it - maximum number of linear iterations
294: . -ksp_min_it - minimum number of linear iterations to use, defaults to zero
295: . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
296: if residual norm decreases by this factor than convergence is declared
297: . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
298: norm is less than this then convergence is declared
299: . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
300: . -ksp_converged_use_initial_residual_norm - see `KSPConvergedDefaultSetUIRNorm()`
301: . -ksp_converged_use_min_initial_residual_norm - see `KSPConvergedDefaultSetUMIRNorm()`
302: . -ksp_converged_maxits - see `KSPConvergedDefaultSetConvergedMaxits()`
303: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural> - see `KSPSetNormType()`
304: . -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
305: works only for `KSPBCGS`, `KSPIBCGS`, and `KSPCG`
306: . -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration;
307: this means that one can use the norm of the residual for convergence test WITHOUT
308: an extra `MPI_Allreduce()` limiting global synchronizations.
309: This will require 1 more iteration of the solver than usual.
310: . -ksp_guess_type - Type of initial guess generator for repeated linear solves
311: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
312: . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
313: . -ksp_test_null_space - tests the null space set with `MatSetNullSpace()` to see if it truly is a null space
314: . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
315: . -ksp_monitor_cancel - cancel all previous convergene monitor routines set
316: . -ksp_monitor - print residual norm at each iteration
317: . -ksp_monitor draw::draw_lg - plot residual norm at each iteration
318: . -ksp_monitor_true_residual - print true residual norm at each iteration
319: . -all_ksp_monitor <optional filename> - print residual norm at each iteration for ALL KSP solves, regardless of their prefix. This is
320: useful for `PCFIELDSPLIT`, `PCMG`, etc that have inner solvers and
321: you wish to track the convergence of all the solvers
322: . -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration
323: . -ksp_monitor_singular_value - monitor extreme singular values at each iteration
324: . -ksp_converged_reason - view the convergence state at the end of the solve
325: . -ksp_use_explicittranspose - transpose the system explicitly in KSPSolveTranspose
326: . -ksp_error_if_not_converged - stop the program as soon as an error is detected in a `KSPSolve()`, `KSP_DIVERGED_ITS`
327: is not treated as an error on inner solves
328: - -ksp_converged_rate - view the convergence rate at the end of the solve
330: Level: beginner
332: Note:
333: To see all options, run your program with the `-help` option or consult [](ch_ksp)
335: .seealso: [](ch_ksp), `KSP`, `KSPSetOptionsPrefix()`, `KSPResetFromOptions()`, `KSPSetUseFischerGuess()`
336: @*/
337: PetscErrorCode KSPSetFromOptions(KSP ksp)
338: {
339: const char *convtests[] = {"default", "skip", "lsqr"}, *prefix;
340: char type[256], guesstype[256], monfilename[PETSC_MAX_PATH_LEN];
341: PetscBool flg, flag, reuse, set;
342: PetscInt indx, model[2] = {0, 0}, nmax;
343: KSPNormType normtype;
344: PCSide pcside;
345: void *ctx;
346: MPI_Comm comm;
348: PetscFunctionBegin;
350: PetscCall(KSPCheckPCMPI(ksp));
352: PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
353: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
354: if (!ksp->skippcsetfromoptions) {
355: if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
356: PetscCall(PCSetFromOptions(ksp->pc));
357: }
359: PetscCall(KSPRegisterAll());
360: PetscObjectOptionsBegin((PetscObject)ksp);
361: PetscCall(PetscOptionsFList("-ksp_type", "Krylov method", "KSPSetType", KSPList, (char *)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES), type, 256, &flg));
362: if (flg) PetscCall(KSPSetType(ksp, type));
363: /*
364: Set the type if it was never set.
365: */
366: if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp, KSPGMRES));
368: PetscCall(KSPResetViewers(ksp));
370: /* Cancels all monitors hardwired into code before call to KSPSetFromOptions() */
371: PetscCall(PetscOptionsBool("-ksp_monitor_cancel", "Remove any hardwired monitor routines", "KSPMonitorCancel", PETSC_FALSE, &flg, &set));
372: if (set && flg) PetscCall(KSPMonitorCancel(ksp));
373: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor", "preconditioned_residual", NULL));
374: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_short", "preconditioned_residual_short", NULL));
375: PetscCall(KSPMonitorSetFromOptions(ksp, "-all_ksp_monitor", "preconditioned_residual", NULL));
376: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_range", "preconditioned_residual_range", NULL));
377: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_true_residual", "true_residual", NULL));
378: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_max", "true_residual_max", NULL));
379: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_solution", "solution", NULL));
380: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_singular_value", "singular_value", ksp));
381: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_monitor_error", "error", ksp));
382: PetscCall(PetscOptionsBool("-ksp_monitor_pause_final", "Pauses all draw monitors at the final iterate", "KSPMonitorPauseFinal_Internal", PETSC_FALSE, &ksp->pauseFinal, NULL));
383: 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));
384: if (flg) PetscCall(KSPSetInitialGuessNonzero(ksp, flag));
386: PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg));
387: if (flg) {
388: PetscCall(KSPGetReusePreconditioner(ksp, &reuse));
389: PetscCall(PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL));
390: PetscCall(KSPSetReusePreconditioner(ksp, reuse));
391: PetscCall(PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL));
392: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view));
393: flg = PETSC_FALSE;
394: PetscCall(PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set));
395: if (set && flg) PetscCall(KSPConvergedReasonViewCancel(ksp));
396: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat));
397: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat));
398: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs));
399: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol));
400: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp));
401: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes));
402: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp));
403: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale));
405: PetscCall(KSPGetDiagonalScale(ksp, &flag));
406: PetscCall(PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg));
407: if (flg) PetscCall(KSPSetDiagonalScale(ksp, flag));
408: PetscCall(KSPGetDiagonalScaleFix(ksp, &flag));
409: PetscCall(PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg));
410: if (flg) PetscCall(KSPSetDiagonalScaleFix(ksp, flag));
411: nmax = ksp->nmax;
412: PetscCall(PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL));
413: PetscCall(PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg));
414: if (flg) PetscCall(KSPSetMatSolveBatchSize(ksp, nmax));
415: goto skipoptions;
416: }
418: PetscCall(PetscOptionsBoundedInt("-ksp_max_it", "Maximum number of iterations", "KSPSetTolerances", ksp->max_it, &ksp->max_it, NULL, 0));
419: PetscCall(PetscOptionsRangeInt("-ksp_min_it", "Minimum number of iterations", "KSPSetMinimumIterations", ksp->min_it, &ksp->min_it, NULL, 0, ksp->max_it));
420: PetscCall(PetscOptionsReal("-ksp_rtol", "Relative decrease in residual norm", "KSPSetTolerances", ksp->rtol, &ksp->rtol, NULL));
421: PetscCall(PetscOptionsReal("-ksp_atol", "Absolute value of residual norm", "KSPSetTolerances", ksp->abstol, &ksp->abstol, NULL));
422: PetscCall(PetscOptionsReal("-ksp_divtol", "Residual norm increase cause divergence", "KSPSetTolerances", ksp->divtol, &ksp->divtol, NULL));
424: PetscCall(PetscOptionsBool("-ksp_converged_use_initial_residual_norm", "Use initial residual norm for computing relative convergence", "KSPConvergedDefaultSetUIRNorm", PETSC_FALSE, &flag, &set));
425: if (set && flag) PetscCall(KSPConvergedDefaultSetUIRNorm(ksp));
426: 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));
427: if (set && flag) PetscCall(KSPConvergedDefaultSetUMIRNorm(ksp));
428: PetscCall(PetscOptionsBool("-ksp_converged_maxits", "Declare convergence if the maximum number of iterations is reached", "KSPConvergedDefaultSetConvergedMaxits", PETSC_FALSE, &flag, &set));
429: if (set) PetscCall(KSPConvergedDefaultSetConvergedMaxits(ksp, flag));
430: PetscCall(KSPGetConvergedNegativeCurvature(ksp, &flag));
431: PetscCall(PetscOptionsBool("-ksp_converged_neg_curve", "Declare convergence if negative curvature is detected", "KSPConvergedNegativeCurvature", flag, &flag, &set));
432: if (set) PetscCall(KSPSetConvergedNegativeCurvature(ksp, flag));
433: PetscCall(KSPGetReusePreconditioner(ksp, &reuse));
434: PetscCall(PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL));
435: PetscCall(KSPSetReusePreconditioner(ksp, reuse));
437: PetscCall(PetscOptionsBool("-ksp_knoll", "Use preconditioner applied to b for initial guess", "KSPSetInitialGuessKnoll", ksp->guess_knoll, &ksp->guess_knoll, NULL));
438: PetscCall(PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL));
439: PetscCall(PetscOptionsFList("-ksp_guess_type", "Initial guess in Krylov method", NULL, KSPGuessList, NULL, guesstype, 256, &flg));
440: if (flg) {
441: PetscCall(KSPGetGuess(ksp, &ksp->guess));
442: PetscCall(KSPGuessSetType(ksp->guess, guesstype));
443: PetscCall(KSPGuessSetFromOptions(ksp->guess));
444: } else { /* old option for KSP */
445: nmax = 2;
446: PetscCall(PetscOptionsIntArray("-ksp_fischer_guess", "Use Paul Fischer's algorithm or its variants for initial guess", "KSPSetUseFischerGuess", model, &nmax, &flag));
447: if (flag) {
448: PetscCheck(nmax == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in model,size as arguments");
449: PetscCall(KSPSetUseFischerGuess(ksp, model[0], model[1]));
450: }
451: }
453: PetscCall(PetscOptionsEList("-ksp_convergence_test", "Convergence test", "KSPSetConvergenceTest", convtests, 3, "default", &indx, &flg));
454: if (flg) {
455: switch (indx) {
456: case 0:
457: PetscCall(KSPConvergedDefaultCreate(&ctx));
458: PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy));
459: break;
460: case 1:
461: PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedSkip, NULL, NULL));
462: break;
463: case 2:
464: PetscCall(KSPConvergedDefaultCreate(&ctx));
465: PetscCall(KSPSetConvergenceTest(ksp, KSPLSQRConvergedDefault, ctx, KSPConvergedDefaultDestroy));
466: break;
467: }
468: }
470: PetscCall(KSPSetUpNorms_Private(ksp, PETSC_FALSE, &normtype, NULL));
471: PetscCall(PetscOptionsEnum("-ksp_norm_type", "KSP Norm type", "KSPSetNormType", KSPNormTypes, (PetscEnum)normtype, (PetscEnum *)&normtype, &flg));
472: if (flg) PetscCall(KSPSetNormType(ksp, normtype));
474: PetscCall(PetscOptionsInt("-ksp_check_norm_iteration", "First iteration to compute residual norm", "KSPSetCheckNormIteration", ksp->chknorm, &ksp->chknorm, NULL));
476: PetscCall(PetscOptionsBool("-ksp_lag_norm", "Lag the calculation of the residual norm", "KSPSetLagNorm", ksp->lagnorm, &flag, &flg));
477: if (flg) PetscCall(KSPSetLagNorm(ksp, flag));
479: PetscCall(KSPGetDiagonalScale(ksp, &flag));
480: PetscCall(PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg));
481: if (flg) PetscCall(KSPSetDiagonalScale(ksp, flag));
482: PetscCall(KSPGetDiagonalScaleFix(ksp, &flag));
483: PetscCall(PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg));
484: if (flg) PetscCall(KSPSetDiagonalScaleFix(ksp, flag));
486: PetscCall(PetscOptionsBool("-ksp_constant_null_space", "Add constant null space to Krylov solver matrix", "MatSetNullSpace", PETSC_FALSE, &flg, &set));
487: if (set && flg) {
488: MatNullSpace nsp;
489: Mat Amat = NULL;
491: PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, &nsp));
492: if (ksp->pc) PetscCall(PCGetOperators(ksp->pc, &Amat, NULL));
493: PetscCheck(Amat, comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot set nullspace, matrix has not yet been provided");
494: PetscCall(MatSetNullSpace(Amat, nsp));
495: PetscCall(MatNullSpaceDestroy(&nsp));
496: }
498: flg = PETSC_FALSE;
499: if (ksp->pc) {
500: PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCKSP, &flg));
501: if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCBJACOBI, &flg));
502: if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCDEFLATION, &flg));
503: }
505: if (flg) {
506: /* Using dynamic tolerance in preconditioner */
507: PetscCall(PetscOptionsString("-sub_ksp_dynamic_tolerance", "Use dynamic tolerance for inner PC", "KSPMonitorDynamicTolerance", "stdout", monfilename, sizeof(monfilename), &flg));
508: if (flg) {
509: void *scale;
510: PetscReal coeff = 1.0;
512: PetscCall(KSPMonitorDynamicToleranceCreate(&scale));
513: PetscCall(PetscOptionsReal("-sub_ksp_dynamic_tolerance", "Coefficient of dynamic tolerance for inner PC", "KSPMonitorDynamicTolerance", coeff, &coeff, &flg));
514: if (flg) PetscCall(KSPMonitorDynamicToleranceSetCoefficient(scale, coeff));
515: PetscCall(KSPMonitorSet(ksp, KSPMonitorDynamicTolerance, scale, KSPMonitorDynamicToleranceDestroy));
516: }
517: }
519: /*
520: Calls Python function
521: */
522: PetscCall(PetscOptionsString("-ksp_monitor_python", "Use Python function", "KSPMonitorSet", NULL, monfilename, sizeof(monfilename), &flg));
523: if (flg) PetscCall(PetscPythonMonitorSet((PetscObject)ksp, monfilename));
524: /*
525: Graphically plots preconditioned residual norm and range of residual element values
526: */
527: PetscCall(PetscOptionsBool("-ksp_monitor_lg_range", "Monitor graphically range of preconditioned residual norm", "KSPMonitorSet", PETSC_FALSE, &flg, &set));
528: if (set && flg) {
529: PetscViewer ctx;
531: PetscCall(PetscViewerDrawOpen(comm, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &ctx));
532: PetscCall(KSPMonitorSet(ksp, KSPMonitorLGRange, ctx, (PetscErrorCode(*)(void **))PetscViewerDestroy));
533: }
534: /* TODO Do these show up in help? */
535: PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &flg));
536: if (flg) {
537: const char *RateTypes[] = {"default", "residual", "error", "PetscRateType", "RATE_", NULL};
538: PetscEnum rtype = (PetscEnum)1;
540: PetscCall(PetscOptionsGetEnum(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate_type", RateTypes, &rtype, &flg));
541: if (rtype == (PetscEnum)0 || rtype == (PetscEnum)1) PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE));
542: if (rtype == (PetscEnum)0 || rtype == (PetscEnum)2) PetscCall(KSPSetErrorHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE));
543: }
544: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view));
545: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pre", &ksp->viewerPre, &ksp->formatPre, &ksp->viewPre));
547: flg = PETSC_FALSE;
548: PetscCall(PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set));
549: if (set && flg) PetscCall(KSPConvergedReasonViewCancel(ksp));
550: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &ksp->viewerRate, &ksp->formatRate, &ksp->viewRate));
551: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat));
552: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat));
553: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs));
554: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol));
555: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp));
556: PetscCall(PetscOptionsDeprecated("-ksp_compute_eigenvalues", "-ksp_view_eigenvalues", "3.9", NULL));
557: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues", &ksp->viewerEV, &ksp->formatEV, &ksp->viewEV));
558: PetscCall(PetscOptionsDeprecated("-ksp_compute_singularvalues", "-ksp_view_singularvalues", "3.9", NULL));
559: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_singularvalues", &ksp->viewerSV, &ksp->formatSV, &ksp->viewSV));
560: PetscCall(PetscOptionsDeprecated("-ksp_compute_eigenvalues_explicitly", "-ksp_view_eigenvalues_explicit", "3.9", NULL));
561: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues_explicit", &ksp->viewerEVExp, &ksp->formatEVExp, &ksp->viewEVExp));
562: PetscCall(PetscOptionsDeprecated("-ksp_final_residual", "-ksp_view_final_residual", "3.9", NULL));
563: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes));
564: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp));
565: PetscCall(PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale));
567: /* Deprecated options */
568: if (!ksp->viewEV) {
569: /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigenvalues",...) below because the argument handling is different */
570: PetscCall(PetscOptionsDeprecated("-ksp_plot_eigenvalues", NULL, "3.9", "Use -ksp_view_eigenvalues draw"));
571: PetscCall(PetscOptionsName("-ksp_plot_eigenvalues", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw]", "KSPView", &ksp->viewEV));
572: if (ksp->viewEV) {
573: ksp->formatEV = PETSC_VIEWER_DEFAULT;
574: ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
575: PetscCall(PetscObjectReference((PetscObject)ksp->viewerEV));
576: }
577: }
578: if (!ksp->viewEV) {
579: PetscCall(PetscOptionsDeprecated("-ksp_plot_eigenvalues_explicitly", NULL, "3.9", "Use -ksp_view_eigenvalues_explicit draw"));
580: /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigencontours",...) below because the argument handling is different */
581: PetscCall(PetscOptionsName("-ksp_plot_eigencontours", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw::draw_contour]", "KSPView", &ksp->viewEV));
582: if (ksp->viewEV) {
583: ksp->formatEV = PETSC_VIEWER_DRAW_CONTOUR;
584: ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
585: PetscCall(PetscObjectReference((PetscObject)ksp->viewerEV));
586: }
587: }
588: if (!ksp->viewEVExp) {
589: /* Cannot remove the what otherwise would be redundant call to PetscOptionsName("-ksp_plot_eigencontours_explicitly",...) below because the argument handling is different */
590: PetscCall(PetscOptionsName("-ksp_plot_eigenvalues_explicitly", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues_explicit draw]", "KSPView", &ksp->viewEVExp));
591: if (ksp->viewEVExp) {
592: ksp->formatEVExp = PETSC_VIEWER_DEFAULT;
593: ksp->viewerEVExp = PETSC_VIEWER_DRAW_(comm);
594: PetscCall(PetscObjectReference((PetscObject)ksp->viewerEVExp));
595: }
596: }
598: #if defined(PETSC_HAVE_SAWS)
599: /*
600: Publish convergence information using AMS
601: */
602: PetscCall(PetscOptionsBool("-ksp_monitor_saws", "Publish KSP progress using SAWs", "KSPMonitorSet", PETSC_FALSE, &flg, &set));
603: if (set && flg) {
604: void *ctx;
605: PetscCall(KSPMonitorSAWsCreate(ksp, &ctx));
606: PetscCall(KSPMonitorSet(ksp, KSPMonitorSAWs, ctx, KSPMonitorSAWsDestroy));
607: PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
608: }
609: #endif
611: PetscCall(KSPSetUpNorms_Private(ksp, PETSC_FALSE, NULL, &pcside));
612: PetscCall(PetscOptionsEnum("-ksp_pc_side", "KSP preconditioner side", "KSPSetPCSide", PCSides, (PetscEnum)pcside, (PetscEnum *)&pcside, &flg));
613: if (flg) PetscCall(KSPSetPCSide(ksp, pcside));
615: if (ksp->viewSV || ksp->viewEV) PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
617: #if defined(PETSC_HAVE_SAWS)
618: {
619: PetscBool set;
620: flg = PETSC_FALSE;
621: PetscCall(PetscOptionsBool("-ksp_saws_block", "Block for SAWs at end of KSPSolve", "PetscObjectSAWsBlock", ((PetscObject)ksp)->amspublishblock, &flg, &set));
622: if (set) PetscCall(PetscObjectSAWsSetBlock((PetscObject)ksp, flg));
623: }
624: #endif
626: nmax = ksp->nmax;
627: PetscCall(PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL));
628: PetscCall(PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg));
629: if (flg) PetscCall(KSPSetMatSolveBatchSize(ksp, nmax));
631: flg = PETSC_FALSE;
632: PetscCall(PetscOptionsBool("-ksp_use_explicittranspose", "Explicitly transpose the system in KSPSolveTranspose", "KSPSetUseExplicitTranspose", ksp->transpose.use_explicittranspose, &flg, &set));
633: if (set) PetscCall(KSPSetUseExplicitTranspose(ksp, flg));
635: PetscTryTypeMethod(ksp, setfromoptions, PetscOptionsObject);
636: skipoptions:
637: /* process any options handlers added with PetscObjectAddOptionsHandler() */
638: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)ksp, PetscOptionsObject));
639: PetscOptionsEnd();
640: ksp->setfromoptionscalled++;
641: PetscFunctionReturn(PETSC_SUCCESS);
642: }
644: /*@
645: KSPResetFromOptions - Sets `KSP` parameters from user options ONLY if the `KSP` was previously set from options
647: Collective
649: Input Parameter:
650: . ksp - the `KSP` context
652: Level: advanced
654: .seealso: [](ch_ksp), `KSPSetFromOptions()`, `KSPSetOptionsPrefix()`
655: @*/
656: PetscErrorCode KSPResetFromOptions(KSP ksp)
657: {
658: PetscFunctionBegin;
659: if (ksp->setfromoptionscalled) PetscCall(KSPSetFromOptions(ksp));
660: PetscFunctionReturn(PETSC_SUCCESS);
661: }