Skip to content

Conversation

pv
Copy link
Member

@pv pv commented Feb 10, 2018

The sparse.linalg iterative solvers treat the tol parameter in an inconsistent way. This PR rectifies the situation preserving backward compatibility, and adds deprecation warnings preparing the way for more sane defaults in future Scipy releases.

The fortran-based solvers have the 'legacy' behavior: tol is considered as the tolerance relative to norm(b), and also as the absolute tolerance for the initial residual (but not subsequent ones!). For lgmres, gcrotmk, the legacy behavior is to consider tol both as an absolute and relative tolerance. These legacy behaviors are wonky, and for most cases not what is wanted. The usually needed interpretation for tol would be the relative tolerance only.

This PR makes it possible to specify both relative and absolute tolerances, by adding an atol= keyword argument. If only tol= is specified, we emit a DeprecationWarning and follow the legacy behavior. In the future, we may want to change the default behavior to atol=0.

Moreover, this PR fixes the behavior of the tolerance parameters when there is left-preconditioning. Previously, gmres applied these inconsistently in the stop condition and inner loop stop condition. The new behavior is to consider the tolerances to apply to the non-preconditioned residual.

The inner loop still only knows about the preconditioned residual, so an adaptive tolerance control mechanism is added (also to lgmres). This kicks in only when there is a mismatch between the two conditions (i.e. possible with preconditioning, or when running into rounding error).

I also had to adjust the default tolerances ARPACK uses in iterative inversion --- it wants to run GMRES at tol ~ eps, but it's not necessarily possible to achieve such relative tolerance due to rounding error. So bump the tolerance upwards.

There may be changes in behavior when running at very small tolerances where the solver runs into rounding error, or when dealing with preconditioned problems. Previously, these solvers could exit when the required tolerance was not satisfied, and with the stricter checking they may now exit signaling non-convergence in such cases. However, this is my opinion bug.

Closes gh-7624

@rc: ping

pv added 7 commits February 11, 2018 17:40
Use an explicit docstring rather than a templated one for minres,
because it's not that similar to the other solvers.

Also make the description of the `tol` parameter accurate. It's only the
relative tolerance for this solver.
Also emit deprecationwarnings for calls to gcrotmk, lgmres without a
specified value for atol parameter, to prepare for changing it in the
future.
Add a way to specify absolute tolerance in iterative fortran solvers.

Previously, these solvers had only the 'legacy' behavior: `tol` is the
tolerance relative to norm(b) and also considered as the absolute
tolerance on the first iteration (!). The behavior originates
from the checks of the form

    IF ( <rc>NRM2( N, WORK(1,R), 1 ).LT.TOL) GO TO 200

This commit changes the internal workings of the fortran solvers
so that:

- The solvers internally deal only with the absolute 2-norm convergence
  criterion. The threshold for this criterion is computed before entering
  the solver.

- All BNRM2 computations are moved python-side.

- STOPTEST2 is removed, as it is easily written in python.

- The strange 'legacy' behavior of the tol= argument is preserved, but
  with a deprecation warning.

- Change the first-iteration check in the fortran codes to check for
  ``RESID .LE. TOL`` instead of ``RESID .LT. TOL`` so that zero rhs
  vectors are handled correctly also for tol=0.
…nverses

The 'legacy' convergence criterion used by gmres is more of a bug, and
ARPACK can simply use only the relative-tolerance convergence criterion.
@pv
Copy link
Member Author

pv commented Feb 12, 2018

Not ready yet --- forgot to check with the left-preconditioning...

@pv pv changed the title Fix tolerance specification in sparse.linalg iterative solvers WIP: Fix tolerance specification in sparse.linalg iterative solvers Feb 12, 2018
@pv pv changed the title WIP: Fix tolerance specification in sparse.linalg iterative solvers Fix tolerance specification in sparse.linalg iterative solvers Feb 17, 2018
@pv
Copy link
Member Author

pv commented Feb 17, 2018

ok, should be ready to go now

@pv pv force-pushed the isolve-tol branch 4 times, most recently from 233ad18 to 8defb02 Compare February 17, 2018 16:10
@pv
Copy link
Member Author

pv commented Feb 17, 2018

An arguable point here is whether gmres should simply check only the preconditioned residual for termination (the other solvers don't) rather than trying to satisfy the requirement for the unpreconditioned result. In this case, only relative tolerance is useful, however.

pv added 2 commits February 17, 2018 20:58
…ditioner

Make changes to ensure tol= and atol= kwargs in gmres() apply to the
residual rather than the preconditioned residual.

Be explicit in the GMRESREVCOM that tolerances apply to the
left-preconditioned residual.

When inner loop residual decreases below tolerance in the GMRESREVCOM,
break out from the inner loop (similarly as during restart) rather than
exiting with success. This ensures _stopcheck is run, which checks the
termination condition for the non-preconditioned residual before
exiting.

As the termination conditions for the outer and inner loops may
mismatch, control the inner tolerance adaptively.
The termination condition may not be satisfied for the actual residual
even though it is satisfied for the preconditioned one, so control the
inner loop tolerance separately.
pv added 2 commits February 17, 2018 20:58
The iterative inversion via GMRES cannot necessarily achieve machine
precision, so run it at tolerance similar to that expected from
computation of vector norms.

Also ensure the eps for the correct dtype is used.
@rgommers
Copy link
Member

This should go into 1.1.x I think. @rc if you have time to review, that would be great. If not, I'll give it a go

@rc
Copy link
Contributor

rc commented Mar 26, 2018

I have checked the new code with my preconditioning problem (see #7624) and it works like a charm, thanks! I have also looked at the changes, and they seems fine, although I did not analyze every line.

+1 for merging.

@person142
Copy link
Member

it works like a charm

Great! I also looked through the changes a while back and was happy, so it looks like we're at merge.

@person142 person142 merged commit 9952987 into scipy:master Mar 26, 2018
@rgommers
Copy link
Member

Awesome!

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
maintenance Items related to regular maintenance tasks scipy.sparse.linalg
Projects
None yet
Development

Successfully merging this pull request may close these issues.

allow setting both absolute and relative tolerance of sparse iterative solvers
4 participants