petsc-3.4.5 2014-06-29

KSPSolve

Solves linear system.

Synopsis

#include "petscksp.h" 
PetscErrorCode  KSPSolve(KSP ksp,Vec b,Vec x)
Collective on KSP

Parameter

ksp - iterative context obtained from KSPCreate()
b - the right hand side vector
x - the solution (this may be the same vector as b, then b will be overwritten with answer)

Options Database Keys

-ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
-ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
-ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and useing LAPACK
-ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
-ksp_view_mat binary - save matrix to the default binary viewer
-ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
-ksp_view_rhs binary - save right hand side vector to the default binary viewer
-ksp_view_solution binary - save computed solution vector to the default binary viewer (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
-ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
-ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
-ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
-ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
-ksp_view - print the ksp data structure at the end of the system solution

Notes

If one uses KSPSetDM() then x or b need not be passed. Use KSPGetSolution() to access the solution in this case.

The operator is specified with KSPSetOperators().

Call KSPGetConvergedReason() to determine if the solver converged or failed and why. The number of iterations can be obtained from KSPGetIterationNumber().

If using a direct method (e.g., via the KSP solver KSPPREONLY and a preconditioner such as PCLU/PCILU), then its=1. See KSPSetTolerances() and KSPDefaultConverged() for more details.

Understanding Convergence

The routines KSPMonitorSet(), KSPComputeEigenvalues(), and KSPComputeEigenvaluesExplicitly() provide information on additional options to monitor convergence and print eigenvalue information.

Keywords

KSP, solve, linear system

See Also

KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
KSPSolveTranspose(), KSPGetIterationNumber()

Level:beginner
Location:
src/ksp/ksp/interface/itfunc.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages src/ksp/pc/examples/tutorials/ex1.c.html
src/ksp/pc/examples/tutorials/ex2.c.html
src/ksp/ksp/examples/tutorials/ex1.c.html
src/ksp/ksp/examples/tutorials/ex2.c.html
src/ksp/ksp/examples/tutorials/ex3.c.html
src/ksp/ksp/examples/tutorials/ex4.c.html
src/ksp/ksp/examples/tutorials/ex5.c.html
src/ksp/ksp/examples/tutorials/ex7.c.html
src/ksp/ksp/examples/tutorials/ex8.c.html
src/ksp/ksp/examples/tutorials/ex8g.c.html