Actual source code: itcl.c


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

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

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

 13:    Logically Collective on ksp

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

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

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

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

 37:    Level: advanced

 39: .seealso: `KSPAppendOptionsPrefix()`, `KSPGetOptionsPrefix()`
 40: @*/
 41: PetscErrorCode KSPSetOptionsPrefix(KSP ksp, const char prefix[])
 42: {
 44:   if (!ksp->pc) KSPGetPC(ksp, &ksp->pc);
 45:   PCSetOptionsPrefix(ksp->pc, prefix);
 46:   PetscObjectSetOptionsPrefix((PetscObject)ksp, prefix);
 47:   return 0;
 48: }

 50: /*@C
 51:    KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
 52:    KSP options in the database.

 54:    Logically Collective on ksp

 56:    Input Parameters:
 57: +  ksp - the Krylov context
 58: -  prefix - the prefix string to prepend to all KSP option requests

 60:    Notes:
 61:    A hyphen (-) must NOT be given at the beginning of the prefix name.
 62:    The first character of all runtime options is AUTOMATICALLY the hyphen.

 64:    Level: advanced

 66: .seealso: `KSPSetOptionsPrefix()`, `KSPGetOptionsPrefix()`
 67: @*/
 68: PetscErrorCode KSPAppendOptionsPrefix(KSP ksp, const char prefix[])
 69: {
 71:   if (!ksp->pc) KSPGetPC(ksp, &ksp->pc);
 72:   PCAppendOptionsPrefix(ksp->pc, prefix);
 73:   PetscObjectAppendOptionsPrefix((PetscObject)ksp, prefix);
 74:   return 0;
 75: }

 77: /*@
 78:    KSPSetUseFischerGuess - Use the Paul Fischer algorithm or its variants

 80:    Logically Collective on ksp

 82:    Input Parameters:
 83: +  ksp - the Krylov context
 84: .  model - use model 1, model 2, model 3, or any other number to turn it off
 85: -  size - size of subspace used to generate initial guess

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

 90:    Level: advanced

 92: .seealso: `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetUseFischerGuess()`, `KSPSetGuess()`, `KSPGetGuess()`
 93: @*/
 94: PetscErrorCode KSPSetUseFischerGuess(KSP ksp, PetscInt model, PetscInt size)
 95: {
 96:   KSPGuess guess;

101:   KSPGetGuess(ksp, &guess);
102:   KSPGuessSetType(guess, KSPGUESSFISCHER);
103:   KSPGuessFischerSetModel(guess, model, size);
104:   return 0;
105: }

107: /*@
108:    KSPSetGuess - Set the initial guess object

110:    Logically Collective on ksp

112:    Input Parameters:
113: +  ksp - the Krylov context
114: -  guess - the object created with KSPGuessCreate()

116:    Level: advanced

118:    Notes:
119:     this allows a single KSP to be used with several different initial guess generators (likely for different linear
120:           solvers, see KSPSetPC()).

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

125: .seealso: `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetUseFischerGuess()`, `KSPGetGuess()`
126: @*/
127: PetscErrorCode KSPSetGuess(KSP ksp, KSPGuess guess)
128: {
131:   PetscObjectReference((PetscObject)guess);
132:   KSPGuessDestroy(&ksp->guess);
133:   ksp->guess      = guess;
134:   ksp->guess->ksp = ksp;
135:   return 0;
136: }

138: /*@
139:    KSPGetGuess - Gets the initial guess generator for the KSP.

141:    Not Collective

143:    Input Parameters:
144: .  ksp - the Krylov context

146:    Output Parameters:
147: .   guess - the object

149:    Level: developer

151: .seealso: `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `KSPSetUseFischerGuess()`, `KSPSetGuess()`
152: @*/
153: PetscErrorCode KSPGetGuess(KSP ksp, KSPGuess *guess)
154: {
157:   if (!ksp->guess) {
158:     const char *prefix;

160:     KSPGuessCreate(PetscObjectComm((PetscObject)ksp), &ksp->guess);
161:     PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix);
162:     if (prefix) PetscObjectSetOptionsPrefix((PetscObject)ksp->guess, prefix);
163:     ksp->guess->ksp = ksp;
164:   }
165:   *guess = ksp->guess;
166:   return 0;
167: }

169: /*@C
170:    KSPGetOptionsPrefix - Gets the prefix used for searching for all
171:    KSP options in the database.

173:    Not Collective

175:    Input Parameters:
176: .  ksp - the Krylov context

178:    Output Parameters:
179: .  prefix - pointer to the prefix string used is returned

181:    Notes:
182:     On the fortran side, the user should pass in a string 'prefix' of
183:    sufficient length to hold the prefix.

185:    Level: advanced

187: .seealso: `KSPSetOptionsPrefix()`, `KSPAppendOptionsPrefix()`
188: @*/
189: PetscErrorCode KSPGetOptionsPrefix(KSP ksp, const char *prefix[])
190: {
192:   PetscObjectGetOptionsPrefix((PetscObject)ksp, prefix);
193:   return 0;
194: }

196: static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
197: {
198:   PetscViewerAndFormatCreate(viewer, format, vf);
199:   (*vf)->data = ctx;
200:   return 0;
201: }

203: /*@C
204:    KSPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

206:    Collective on ksp

208:    Input Parameters:
209: +  ksp  - KSP object you wish to monitor
210: .  opt  - the command line option for this monitor
211: .  name - the monitor type one is seeking
212: -  ctx  - An optional user context for the monitor, or NULL

214:    Level: developer

216: .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
217:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
218:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
219:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
220:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
221:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
222:           `PetscOptionsFList()`, `PetscOptionsEList()`
223: @*/
224: PetscErrorCode KSPMonitorSetFromOptions(KSP ksp, const char opt[], const char name[], void *ctx)
225: {
226:   PetscErrorCode (*mfunc)(KSP, PetscInt, PetscReal, void *);
227:   PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
228:   PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
229:   PetscViewerAndFormat *vf;
230:   PetscViewer           viewer;
231:   PetscViewerFormat     format;
232:   PetscViewerType       vtype;
233:   char                  key[PETSC_MAX_PATH_LEN];
234:   PetscBool             all, flg;
235:   const char           *prefix = NULL;

237:   PetscStrcmp(opt, "-all_ksp_monitor", &all);
238:   if (!all) PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix);
239:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->options, prefix, opt, &viewer, &format, &flg);
240:   if (!flg) return 0;

242:   PetscViewerGetType(viewer, &vtype);
243:   KSPMonitorMakeKey_Internal(name, vtype, format, key);
244:   PetscFunctionListFind(KSPMonitorList, key, &mfunc);
245:   PetscFunctionListFind(KSPMonitorCreateList, key, &cfunc);
246:   PetscFunctionListFind(KSPMonitorDestroyList, key, &dfunc);
247:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
248:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

250:   (*cfunc)(viewer, format, ctx, &vf);
251:   PetscObjectDereference((PetscObject)viewer);
252:   KSPMonitorSet(ksp, mfunc, vf, (PetscErrorCode(*)(void **))dfunc);
253:   return 0;
254: }

256: /*@
257:    KSPSetFromOptions - Sets KSP options from the options database.
258:    This routine must be called before KSPSetUp() if the user is to be
259:    allowed to set the Krylov type.

261:    Collective on ksp

263:    Input Parameter:
264: .  ksp - the Krylov space context

266:    Options Database Keys:
267: +   -ksp_max_it - maximum number of linear iterations
268: .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
269:                 if residual norm decreases by this factor than convergence is declared
270: .   -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
271:                 norm is less than this then convergence is declared
272: .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
273: .   -ksp_converged_use_initial_residual_norm - see KSPConvergedDefaultSetUIRNorm()
274: .   -ksp_converged_use_min_initial_residual_norm - see KSPConvergedDefaultSetUMIRNorm()
275: .   -ksp_converged_maxits - see KSPConvergedDefaultSetConvergedMaxits()
276: .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
277:                        convergence test (say you always want to run with 5 iterations) to
278:                        save on communication overhead
279:                     preconditioned - default for left preconditioning
280:                     unpreconditioned - see KSPSetNormType()
281:                     natural - see KSPSetNormType()
282: .   -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
283:        works only for PCBCGS, PCIBCGS and and PCCG
284: .   -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use
285:        the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations.
286:        This will require 1 more iteration of the solver than usual.
287: .   -ksp_guess_type - Type of initial guess generator for repeated linear solves
288: .   -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
289: .   -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
290: .   -ksp_test_null_space - tests the null space set with MatSetNullSpace() to see if it truly is a null space
291: .   -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
292: .   -ksp_monitor_cancel - cancel all previous convergene monitor routines set
293: .   -ksp_monitor - print residual norm at each iteration
294: .   -ksp_monitor draw::draw_lg - plot residual norm at each iteration
295: .   -ksp_monitor_true_residual - print true residual norm at each iteration
296: .   -all_ksp_monitor <optional filename> - print residual norm at each iteration for ALL KSP solves, regardless of their prefix. This is
297:                                            useful for PCFIELDSPLIT, PCMG, etc that have inner solvers and you wish to track the convergence of all the solvers
298: .   -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration
299: .   -ksp_monitor_singular_value - monitor extreme singular values at each iteration
300: .   -ksp_converged_reason - view the convergence state at the end of the solve
301: .   -ksp_use_explicittranspose - transpose the system explicitly in KSPSolveTranspose
302: .   -ksp_error_if_not_converged - stop the program as soon as an error is detected in a KSPSolve(), KSP_DIVERGED_ITS is not treated as an error on inner KSPSolves
303: -   -ksp_converged_rate - view the convergence rate at the end of the solve

305:    Notes:
306:    To see all options, run your program with the -help option
307:    or consult Users-Manual: ch_ksp

309:    Level: beginner

311: .seealso: `KSPSetOptionsPrefix()`, `KSPResetFromOptions()`, `KSPSetUseFischerGuess()`

313: @*/
314: PetscErrorCode KSPSetFromOptions(KSP ksp)
315: {
316:   const char *convtests[] = {"default", "skip", "lsqr"}, *prefix;
317:   char        type[256], guesstype[256], monfilename[PETSC_MAX_PATH_LEN];
318:   PetscBool   flg, flag, reuse, set;
319:   PetscInt    indx, model[2] = {0, 0}, nmax;
320:   KSPNormType normtype;
321:   PCSide      pcside;
322:   void       *ctx;
323:   MPI_Comm    comm;

326:   PetscObjectGetComm((PetscObject)ksp, &comm);
327:   PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix);
328:   if (!ksp->skippcsetfromoptions) {
329:     if (!ksp->pc) KSPGetPC(ksp, &ksp->pc);
330:     PCSetFromOptions(ksp->pc);
331:   }

333:   KSPRegisterAll();
334:   PetscObjectOptionsBegin((PetscObject)ksp);
335:   PetscOptionsFList("-ksp_type", "Krylov method", "KSPSetType", KSPList, (char *)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES), type, 256, &flg);
336:   if (flg) KSPSetType(ksp, type);
337:   /*
338:     Set the type if it was never set.
339:   */
340:   if (!((PetscObject)ksp)->type_name) KSPSetType(ksp, KSPGMRES);

342:   KSPResetViewers(ksp);

344:   /* Cancels all monitors hardwired into code before call to KSPSetFromOptions() */
345:   PetscOptionsBool("-ksp_monitor_cancel", "Remove any hardwired monitor routines", "KSPMonitorCancel", PETSC_FALSE, &flg, &set);
346:   if (set && flg) KSPMonitorCancel(ksp);
347:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor", "preconditioned_residual", NULL);
348:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_short", "preconditioned_residual_short", NULL);
349:   KSPMonitorSetFromOptions(ksp, "-all_ksp_monitor", "preconditioned_residual", NULL);
350:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_range", "preconditioned_residual_range", NULL);
351:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_true_residual", "true_residual", NULL);
352:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_max", "true_residual_max", NULL);
353:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_solution", "solution", NULL);
354:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_singular_value", "singular_value", ksp);
355:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_error", "error", ksp);
356:   PetscOptionsBool("-ksp_monitor_pause_final", "Pauses all draw monitors at the final iterate", "KSPMonitorPauseFinal_Internal", PETSC_FALSE, &ksp->pauseFinal, NULL);
357:   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);
358:   if (flg) KSPSetInitialGuessNonzero(ksp, flag);

360:   PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg);
361:   if (flg) {
362:     KSPGetReusePreconditioner(ksp, &reuse);
363:     PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL);
364:     KSPSetReusePreconditioner(ksp, reuse);
365:     PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL);
366:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view);
367:     flg = PETSC_FALSE;
368:     PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set);
369:     if (set && flg) KSPConvergedReasonViewCancel(ksp);
370:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat);
371:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat);
372:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs);
373:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol);
374:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp);
375:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes);
376:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp);
377:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale);

379:     KSPGetDiagonalScale(ksp, &flag);
380:     PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg);
381:     if (flg) KSPSetDiagonalScale(ksp, flag);
382:     KSPGetDiagonalScaleFix(ksp, &flag);
383:     PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg);
384:     if (flg) KSPSetDiagonalScaleFix(ksp, flag);
385:     nmax = ksp->nmax;
386:     PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL);
387:     PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg);
388:     if (flg) KSPSetMatSolveBatchSize(ksp, nmax);
389:     goto skipoptions;
390:   }

392:   PetscOptionsInt("-ksp_max_it", "Maximum number of iterations", "KSPSetTolerances", ksp->max_it, &ksp->max_it, NULL);
393:   PetscOptionsReal("-ksp_rtol", "Relative decrease in residual norm", "KSPSetTolerances", ksp->rtol, &ksp->rtol, NULL);
394:   PetscOptionsReal("-ksp_atol", "Absolute value of residual norm", "KSPSetTolerances", ksp->abstol, &ksp->abstol, NULL);
395:   PetscOptionsReal("-ksp_divtol", "Residual norm increase cause divergence", "KSPSetTolerances", ksp->divtol, &ksp->divtol, NULL);

397:   PetscOptionsBool("-ksp_converged_use_initial_residual_norm", "Use initial residual norm for computing relative convergence", "KSPConvergedDefaultSetUIRNorm", PETSC_FALSE, &flag, &set);
398:   if (set && flag) KSPConvergedDefaultSetUIRNorm(ksp);
399:   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);
400:   if (set && flag) KSPConvergedDefaultSetUMIRNorm(ksp);
401:   PetscOptionsBool("-ksp_converged_maxits", "Declare convergence if the maximum number of iterations is reached", "KSPConvergedDefaultSetConvergedMaxits", PETSC_FALSE, &flag, &set);
402:   if (set) KSPConvergedDefaultSetConvergedMaxits(ksp, flag);
403:   KSPGetReusePreconditioner(ksp, &reuse);
404:   PetscOptionsBool("-ksp_reuse_preconditioner", "Use initial preconditioner and don't ever compute a new one", "KSPReusePreconditioner", reuse, &reuse, NULL);
405:   KSPSetReusePreconditioner(ksp, reuse);

407:   PetscOptionsBool("-ksp_knoll", "Use preconditioner applied to b for initial guess", "KSPSetInitialGuessKnoll", ksp->guess_knoll, &ksp->guess_knoll, NULL);
408:   PetscOptionsBool("-ksp_error_if_not_converged", "Generate error if solver does not converge", "KSPSetErrorIfNotConverged", ksp->errorifnotconverged, &ksp->errorifnotconverged, NULL);
409:   PetscOptionsFList("-ksp_guess_type", "Initial guess in Krylov method", NULL, KSPGuessList, NULL, guesstype, 256, &flg);
410:   if (flg) {
411:     KSPGetGuess(ksp, &ksp->guess);
412:     KSPGuessSetType(ksp->guess, guesstype);
413:     KSPGuessSetFromOptions(ksp->guess);
414:   } else { /* old option for KSP */
415:     nmax = 2;
416:     PetscOptionsIntArray("-ksp_fischer_guess", "Use Paul Fischer's algorithm or its variants for initial guess", "KSPSetUseFischerGuess", model, &nmax, &flag);
417:     if (flag) {
419:       KSPSetUseFischerGuess(ksp, model[0], model[1]);
420:     }
421:   }

423:   PetscOptionsEList("-ksp_convergence_test", "Convergence test", "KSPSetConvergenceTest", convtests, 3, "default", &indx, &flg);
424:   if (flg) {
425:     switch (indx) {
426:     case 0:
427:       KSPConvergedDefaultCreate(&ctx);
428:       KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy);
429:       break;
430:     case 1:
431:       KSPSetConvergenceTest(ksp, KSPConvergedSkip, NULL, NULL);
432:       break;
433:     case 2:
434:       KSPConvergedDefaultCreate(&ctx);
435:       KSPSetConvergenceTest(ksp, KSPLSQRConvergedDefault, ctx, KSPConvergedDefaultDestroy);
436:       break;
437:     }
438:   }

440:   KSPSetUpNorms_Private(ksp, PETSC_FALSE, &normtype, NULL);
441:   PetscOptionsEnum("-ksp_norm_type", "KSP Norm type", "KSPSetNormType", KSPNormTypes, (PetscEnum)normtype, (PetscEnum *)&normtype, &flg);
442:   if (flg) KSPSetNormType(ksp, normtype);

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

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

449:   KSPGetDiagonalScale(ksp, &flag);
450:   PetscOptionsBool("-ksp_diagonal_scale", "Diagonal scale matrix before building preconditioner", "KSPSetDiagonalScale", flag, &flag, &flg);
451:   if (flg) KSPSetDiagonalScale(ksp, flag);
452:   KSPGetDiagonalScaleFix(ksp, &flag);
453:   PetscOptionsBool("-ksp_diagonal_scale_fix", "Fix diagonally scaled matrix after solve", "KSPSetDiagonalScaleFix", flag, &flag, &flg);
454:   if (flg) KSPSetDiagonalScaleFix(ksp, flag);

456:   PetscOptionsBool("-ksp_constant_null_space", "Add constant null space to Krylov solver matrix", "MatSetNullSpace", PETSC_FALSE, &flg, &set);
457:   if (set && flg) {
458:     MatNullSpace nsp;
459:     Mat          Amat = NULL;

461:     MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, &nsp);
462:     if (ksp->pc) PCGetOperators(ksp->pc, &Amat, NULL);
463:     if (Amat) {
464:       MatSetNullSpace(Amat, nsp);
465:       MatNullSpaceDestroy(&nsp);
466:     } else SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot set nullspace, matrix has not yet been provided");
467:   }

469:   flg = PETSC_FALSE;
470:   if (ksp->pc) {
471:     PetscObjectTypeCompare((PetscObject)ksp->pc, PCKSP, &flg);
472:     if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc, PCBJACOBI, &flg);
473:     if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc, PCDEFLATION, &flg);
474:   }

476:   if (flg) {
477:     /* A hack for using dynamic tolerance in preconditioner */
478:     PetscOptionsString("-sub_ksp_dynamic_tolerance", "Use dynamic tolerance for PC if PC is a KSP", "KSPMonitorDynamicTolerance", "stdout", monfilename, sizeof(monfilename), &flg);
479:     if (flg) {
480:       KSPDynTolCtx *scale;
481:       PetscMalloc1(1, &scale);
482:       scale->bnrm = -1.0;
483:       scale->coef = 1.0;
484:       PetscOptionsReal("-sub_ksp_dynamic_tolerance_param", "Parameter of dynamic tolerance for inner PCKSP", "KSPMonitorDynamicToleranceParam", scale->coef, &scale->coef, &flg);
485:       KSPMonitorSet(ksp, KSPMonitorDynamicTolerance, scale, KSPMonitorDynamicToleranceDestroy);
486:     }
487:   }

489:   /*
490:    Calls Python function
491:   */
492:   PetscOptionsString("-ksp_monitor_python", "Use Python function", "KSPMonitorSet", NULL, monfilename, sizeof(monfilename), &flg);
493:   if (flg) PetscPythonMonitorSet((PetscObject)ksp, monfilename);
494:   /*
495:     Graphically plots preconditioned residual norm and range of residual element values
496:   */
497:   PetscOptionsBool("-ksp_monitor_lg_range", "Monitor graphically range of preconditioned residual norm", "KSPMonitorSet", PETSC_FALSE, &flg, &set);
498:   if (set && flg) {
499:     PetscViewer ctx;

501:     PetscViewerDrawOpen(comm, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &ctx);
502:     KSPMonitorSet(ksp, KSPMonitorLGRange, ctx, (PetscErrorCode(*)(void **))PetscViewerDestroy);
503:   }
504:   /* TODO Do these show up in help? */
505:   PetscOptionsHasName(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &flg);
506:   if (flg) {
507:     const char *RateTypes[] = {"default", "residual", "error", "PetscRateType", "RATE_", NULL};
508:     PetscEnum   rtype       = (PetscEnum)1;

510:     PetscOptionsGetEnum(((PetscObject)ksp)->options, prefix, "-ksp_converged_rate_type", RateTypes, &rtype, &flg);
511:     if (rtype == (PetscEnum)0 || rtype == (PetscEnum)1) KSPSetResidualHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE);
512:     if (rtype == (PetscEnum)0 || rtype == (PetscEnum)2) KSPSetErrorHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE);
513:   }
514:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view", &ksp->viewer, &ksp->format, &ksp->view);
515:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pre", &ksp->viewerPre, &ksp->formatPre, &ksp->viewPre);

517:   flg = PETSC_FALSE;
518:   PetscOptionsBool("-ksp_converged_reason_view_cancel", "Cancel all the converged reason view functions set using KSPConvergedReasonViewSet", "KSPConvergedReasonViewCancel", PETSC_FALSE, &flg, &set);
519:   if (set && flg) KSPConvergedReasonViewCancel(ksp);
520:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_converged_rate", &ksp->viewerRate, &ksp->formatRate, &ksp->viewRate);
521:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat", &ksp->viewerMat, &ksp->formatMat, &ksp->viewMat);
522:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_pmat", &ksp->viewerPMat, &ksp->formatPMat, &ksp->viewPMat);
523:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_rhs", &ksp->viewerRhs, &ksp->formatRhs, &ksp->viewRhs);
524:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_solution", &ksp->viewerSol, &ksp->formatSol, &ksp->viewSol);
525:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_mat_explicit", &ksp->viewerMatExp, &ksp->formatMatExp, &ksp->viewMatExp);
526:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues", &ksp->viewerEV, &ksp->formatEV, &ksp->viewEV);
527:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_singularvalues", &ksp->viewerSV, &ksp->formatSV, &ksp->viewSV);
528:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_eigenvalues_explicit", &ksp->viewerEVExp, &ksp->formatEVExp, &ksp->viewEVExp);
529:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes);
530:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_preconditioned_operator_explicit", &ksp->viewerPOpExp, &ksp->formatPOpExp, &ksp->viewPOpExp);
531:   PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_view_diagonal_scale", &ksp->viewerDScale, &ksp->formatDScale, &ksp->viewDScale);

533:   /* Deprecated options */
534:   if (!ksp->viewEV) {
535:     PetscOptionsDeprecated("-ksp_compute_eigenvalues", NULL, "3.9", "Use -ksp_view_eigenvalues");
536:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_compute_eigenvalues", &ksp->viewerEV, &ksp->formatEV, &ksp->viewEV);
537:   }
538:   if (!ksp->viewEV) {
539:     PetscOptionsDeprecated("-ksp_plot_eigenvalues", NULL, "3.9", "Use -ksp_view_eigenvalues draw");
540:     PetscOptionsName("-ksp_plot_eigenvalues", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw]", "KSPView", &ksp->viewEV);
541:     if (ksp->viewEV) {
542:       ksp->formatEV = PETSC_VIEWER_DEFAULT;
543:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
544:       PetscObjectReference((PetscObject)ksp->viewerEV);
545:     }
546:   }
547:   if (!ksp->viewEV) {
548:     PetscOptionsDeprecated("-ksp_plot_eigencontours", NULL, "3.9", "Use -ksp_view_eigenvalues draw::draw_contour");
549:     PetscOptionsName("-ksp_plot_eigencontours", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw::draw_contour]", "KSPView", &ksp->viewEV);
550:     if (ksp->viewEV) {
551:       ksp->formatEV = PETSC_VIEWER_DRAW_CONTOUR;
552:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
553:       PetscObjectReference((PetscObject)ksp->viewerEV);
554:     }
555:   }
556:   if (!ksp->viewEVExp) {
557:     PetscOptionsDeprecated("-ksp_compute_eigenvalues_explicitly", NULL, "3.9", "Use -ksp_view_eigenvalues_explicit");
558:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_compute_eigenvalues_explicitly", &ksp->viewerEVExp, &ksp->formatEVExp, &ksp->viewEVExp);
559:   }
560:   if (!ksp->viewEVExp) {
561:     PetscOptionsDeprecated("-ksp_plot_eigenvalues_explicitly", NULL, "3.9", "Use -ksp_view_eigenvalues_explicit draw");
562:     PetscOptionsName("-ksp_plot_eigenvalues_explicitly", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues_explicit draw]", "KSPView", &ksp->viewEVExp);
563:     if (ksp->viewEVExp) {
564:       ksp->formatEVExp = PETSC_VIEWER_DEFAULT;
565:       ksp->viewerEVExp = PETSC_VIEWER_DRAW_(comm);
566:       PetscObjectReference((PetscObject)ksp->viewerEVExp);
567:     }
568:   }
569:   if (!ksp->viewSV) {
570:     PetscOptionsDeprecated("-ksp_compute_singularvalues", NULL, "3.9", "Use -ksp_view_singularvalues");
571:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_compute_singularvalues", &ksp->viewerSV, &ksp->formatSV, &ksp->viewSV);
572:   }
573:   if (!ksp->viewFinalRes) {
574:     PetscOptionsDeprecated("-ksp_final_residual", NULL, "3.9", "Use -ksp_view_final_residual");
575:     PetscOptionsGetViewer(comm, ((PetscObject)ksp)->options, prefix, "-ksp_final_residual", &ksp->viewerFinalRes, &ksp->formatFinalRes, &ksp->viewFinalRes);
576:   }

578: #if defined(PETSC_HAVE_SAWS)
579:   /*
580:     Publish convergence information using AMS
581:   */
582:   PetscOptionsBool("-ksp_monitor_saws", "Publish KSP progress using SAWs", "KSPMonitorSet", PETSC_FALSE, &flg, &set);
583:   if (set && flg) {
584:     void *ctx;
585:     KSPMonitorSAWsCreate(ksp, &ctx);
586:     KSPMonitorSet(ksp, KSPMonitorSAWs, ctx, KSPMonitorSAWsDestroy);
587:     KSPSetComputeSingularValues(ksp, PETSC_TRUE);
588:   }
589: #endif

591:   /* -----------------------------------------------------------------------*/
592:   KSPSetUpNorms_Private(ksp, PETSC_FALSE, NULL, &pcside);
593:   PetscOptionsEnum("-ksp_pc_side", "KSP preconditioner side", "KSPSetPCSide", PCSides, (PetscEnum)pcside, (PetscEnum *)&pcside, &flg);
594:   if (flg) KSPSetPCSide(ksp, pcside);

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

598: #if defined(PETSC_HAVE_SAWS)
599:   {
600:     PetscBool set;
601:     flg = PETSC_FALSE;
602:     PetscOptionsBool("-ksp_saws_block", "Block for SAWs at end of KSPSolve", "PetscObjectSAWsBlock", ((PetscObject)ksp)->amspublishblock, &flg, &set);
603:     if (set) PetscObjectSAWsSetBlock((PetscObject)ksp, flg);
604:   }
605: #endif

607:   nmax = ksp->nmax;
608:   PetscOptionsDeprecated("-ksp_matsolve_block_size", "-ksp_matsolve_batch_size", "3.15", NULL);
609:   PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg);
610:   if (flg) KSPSetMatSolveBatchSize(ksp, nmax);

612:   flg = PETSC_FALSE;
613:   PetscOptionsBool("-ksp_use_explicittranspose", "Explicitly transpose the system in KSPSolveTranspose", "KSPSetUseExplicitTranspose", ksp->transpose.use_explicittranspose, &flg, &set);
614:   if (set) KSPSetUseExplicitTranspose(ksp, flg);

616:   PetscTryTypeMethod(ksp, setfromoptions, PetscOptionsObject);
617: skipoptions:
618:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
619:   PetscObjectProcessOptionsHandlers((PetscObject)ksp, PetscOptionsObject);
620:   PetscOptionsEnd();
621:   ksp->setfromoptionscalled++;
622:   return 0;
623: }

625: /*@
626:    KSPResetFromOptions - Sets various KSP parameters from user options ONLY if the KSP was previously set from options

628:    Collective on ksp

630:    Input Parameter:
631: .  ksp - the KSP context

633:    Level: beginner

635: .seealso: `KSPSetFromOptions()`, `KSPSetOptionsPrefix()`
636: @*/
637: PetscErrorCode KSPResetFromOptions(KSP ksp)
638: {
639:   if (ksp->setfromoptionscalled) KSPSetFromOptions(ksp);
640:   return 0;
641: }