1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPLSQR.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>KSPLSQR</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.4.2 2013-07-02</b></div>
<A NAME="KSPLSQR"><H1>KSPLSQR</H1></A>
This implements LSQR
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_lsqr_set_standard_error </B></TD><TD>- Set Standard Error Estimates of Solution see KSPLSQRSetStandardErrorVec()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_lsqr_monitor </B></TD><TD>- Monitor residual norm and norm of residual of normal equations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>see <A HREF="../KSP/KSPSolve.html#KSPSolve">KSPSolve</A>()</B></TD><TD>-
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
This varient, when applied with no preconditioning is identical to the original algorithm in exact arithematic; however, in practice, with no preconditioning
due to inexact arithematic, it can converge differently. Hence when no preconditioner is used (<A HREF="../PC/PCType.html#PCType">PCType</A> <A HREF="../PC/PCNONE.html#PCNONE">PCNONE</A>) it automatically reverts to the original algorithm.
<P>
With the PETSc built-in preconditioners, such as ICC, one should call <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>(ksp,A,A'*A,...) since the preconditioner needs to work
for the normal equations A'*A.
<P>
Supports only left preconditioning.
<P>
References:The original unpreconditioned algorithm can be found in Paige and Saunders, ACM Transactions on Mathematical Software, Vol 8, pp 43-71, 1982.
In exact arithmetic the LSQR method (with no preconditioning) is identical to the <A HREF="../KSP/KSPCG.html#KSPCG">KSPCG</A> algorithm applied to the normal equations.
The preconditioned varient was implemented by Bas van't Hof and is essentially a left preconditioning for the Normal Equations. It appears the implementation with preconditioner
track the true norm of the residual and uses that in the convergence test.
<P>
Developer Notes: How is this related to the <A HREF="../KSP/KSPCGNE.html#KSPCGNE">KSPCGNE</A> implementation? One difference is that <A HREF="../KSP/KSPCGNE.html#KSPCGNE">KSPCGNE</A> applies
the preconditioner transpose times the preconditioner, so one does not need to pass A'*A as the third argument to <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>().
<P>
<P>
For least squares problems without a zero to A*x = b, there are additional convergence tests for the residual of the normal equations, A'*(b - Ax), see <A HREF="../KSP/KSPLSQRDefaultConverged.html#KSPLSQRDefaultConverged">KSPLSQRDefaultConverged</A>()
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../KSP/KSPCreate.html#KSPCreate">KSPCreate</A>(), <A HREF="../KSP/KSPSetType.html#KSPSetType">KSPSetType</A>(), <A HREF="../KSP/KSPType.html#KSPType">KSPType</A> (for list of available types), <A HREF="../KSP/KSP.html#KSP">KSP</A>, <A HREF="../KSP/KSPLSQRDefaultConverged.html#KSPLSQRDefaultConverged">KSPLSQRDefaultConverged</A>()
<BR>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>beginner
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/ksp/impls/lsqr/lsqr.c.html#KSPLSQR">src/ksp/ksp/impls/lsqr/lsqr.c</A>
<BR><A HREF="./index.html">Index of all KSP routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|