KSPSetFromOptions#

Sets KSP options from the options database. This routine must be called before KSPSetUp() if the user is to be allowed to set the Krylov type.

Synopsis#

#include "petscksp.h" 
PetscErrorCode KSPSetFromOptions(KSP ksp)

Collective

Input Parameter#

  • ksp - the Krylov space context

Options Database Keys#

  • -ksp_max_it - maximum number of linear iterations

  • -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e. if residual norm decreases by this factor than convergence is declared

  • -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual norm is less than this then convergence is declared

  • -ksp_divtol tol - if residual norm increases by this factor than divergence is declared

  • -ksp_converged_use_initial_residual_norm - see KSPConvergedDefaultSetUIRNorm()

  • -ksp_converged_use_min_initial_residual_norm - see KSPConvergedDefaultSetUMIRNorm()

  • -ksp_converged_maxits - see KSPConvergedDefaultSetConvergedMaxits()

  • -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using convergence test (say you always want to run with 5 iterations) to save on communication overhead preconditioned - default for left preconditioning unpreconditioned - see KSPSetNormType() natural - see KSPSetNormType()

  • -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration) works only for KSPBCGS, KSPIBCGS and and KSPCG

  • -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations. This will require 1 more iteration of the solver than usual.

  • -ksp_guess_type - Type of initial guess generator for repeated linear solves

  • -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves

  • -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space

  • -ksp_test_null_space - tests the null space set with MatSetNullSpace() to see if it truly is a null space

  • -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side

  • -ksp_monitor_cancel - cancel all previous convergene monitor routines set

  • -ksp_monitor - print residual norm at each iteration

  • -ksp_monitor draw::draw_lg - plot residual norm at each iteration

  • -ksp_monitor_true_residual - print true residual norm at each iteration

  • -all_ksp_monitor - print residual norm at each iteration for ALL KSP solves, regardless of their prefix. This is useful for PCFIELDSPLIT, PCMG, etc that have inner solvers and you wish to track the convergence of all the solvers

  • -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration

  • -ksp_monitor_singular_value - monitor extreme singular values at each iteration

  • -ksp_converged_reason - view the convergence state at the end of the solve

  • -ksp_use_explicittranspose - transpose the system explicitly in KSPSolveTranspose

  • -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 solves

  • -ksp_converged_rate - view the convergence rate at the end of the solve

Notes#

To see all options, run your program with the -help option or consult KSP: Linear System Solvers

See Also#

KSP: Linear System Solvers, KSP, KSPSetOptionsPrefix(), KSPResetFromOptions(), KSPSetUseFischerGuess()

Level#

beginner

Location#

src/ksp/ksp/interface/itcl.c

Examples#

src/dm/impls/stag/tutorials/ex1.c
src/dm/impls/stag/tutorials/ex2.c
src/dm/impls/stag/tutorials/ex3.c
src/dm/impls/stag/tutorials/ex4.c
src/dm/impls/stag/tutorials/ex8.c
src/dm/impls/swarm/tutorials/ex1.c
src/dm/impls/swarm/tutorials/ex1f90.F90
src/ksp/ksp/tutorials/bench_kspsolve.c
src/ksp/ksp/tutorials/ex1.c
src/ksp/ksp/tutorials/ex10.c
src/ksp/ksp/tutorials/ex100.c

Implementations#

KSPSetFromOptions_BCGS in src/ksp/ksp/impls/bcgs/bcgs.c
KSPSetFromOptions_BCGSL in src/ksp/ksp/impls/bcgsl/bcgsl.c
KSPSetFromOptions_CG in src/ksp/ksp/impls/cg/cg.c
KSPSetFromOptions_PIPELCG in src/ksp/ksp/impls/cg/pipelcg/pipelcg.c
KSPSetFromOptions_PIPEPRCG in src/ksp/ksp/impls/cg/pipeprcg/pipeprcg.c
KSPSetFromOptions_Chebyshev in src/ksp/ksp/impls/cheby/cheby.c
KSPSetFromOptions_FCG in src/ksp/ksp/impls/fcg/fcg.c
KSPSetFromOptions_PIPEFCG in src/ksp/ksp/impls/fcg/pipefcg/pipefcg.c
KSPSetFromOptions_FETIDP in src/ksp/ksp/impls/fetidp/fetidp.c
KSPSetFromOptions_GCR in src/ksp/ksp/impls/gcr/gcr.c
KSPSetFromOptions_PIPEGCR in src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c
KSPSetFromOptions_AGMRES in src/ksp/ksp/impls/gmres/agmres/agmres.c
KSPSetFromOptions_DGMRES in src/ksp/ksp/impls/gmres/dgmres/dgmres.c
KSPSetFromOptions_FGMRES in src/ksp/ksp/impls/gmres/fgmres/fgmres.c
KSPSetFromOptions_GMRES in src/ksp/ksp/impls/gmres/gmres.c
KSPSetFromOptions_LGMRES in src/ksp/ksp/impls/gmres/lgmres/lgmres.c
KSPSetFromOptions_PGMRES in src/ksp/ksp/impls/gmres/pgmres/pgmres.c
KSPSetFromOptions_PIPEFGMRES in src/ksp/ksp/impls/gmres/pipefgmres/pipefgmres.c
KSPSetFromOptions_HPDDM in src/ksp/ksp/impls/hpddm/hpddm.cxx
KSPSetFromOptions_LCD in src/ksp/ksp/impls/lcd/lcd.c
KSPSetFromOptions_LSQR in src/ksp/ksp/impls/lsqr/lsqr.c
KSPSetFromOptions_MINRES in src/ksp/ksp/impls/minres/minres.c
KSPSetFromOptions_QCG in src/ksp/ksp/impls/qcg/qcg.c
KSPSetFromOptions_Richardson in src/ksp/ksp/impls/rich/rich.c
KSPSetFromOptions_TSIRM in src/ksp/ksp/impls/tsirm/tsirm.c


Edit on GitLab

Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages