[SciPy-User] fmin_bfgs stuck in infinite loop

Johann Cohen-Tanugi johann.cohen-tanugi@lupm.univ-montp2...
Mon Oct 24 15:46:44 CDT 2011


I posted a comment in 
https://github.com/scipy/scipy/commit/a31acbf966f50665997d1fd1ced41d3a8b3e4187
but I agree that a comprehensive set of checks at optimize.py level is 
warranted, rather than a kludge in a sub-method.
J

On 10/24/2011 10:34 PM, josef.pktd@gmail.com wrote:
> On Mon, Oct 24, 2011 at 4:14 PM,<josef.pktd@gmail.com>  wrote:
>> On Mon, Oct 24, 2011 at 3:59 PM,<josef.pktd@gmail.com>  wrote:
>>> tricky things these reply to all, forwarding to list
>>>
>>> ---------- Forwarded message ----------
>>> From:<josef.pktd@gmail.com>
>>> Date: Mon, Oct 24, 2011 at 3:52 PM
>>> Subject: Re: [SciPy-User] fmin_bfgs stuck in infinite loop
>>> To: johann.cohen-tanugi@lupm.in2p3.fr
>>>
>>>
>>> On Mon, Oct 24, 2011 at 3:43 PM, Johann Cohen-Tanugi
>>> <johann.cohen-tanugi@lupm.univ-montp2.fr>  wrote:
>>>> indeed, see the email I just sent : for nan linesearch_wolfe1 does not exit
>>>> gracefully, so that Ralf's patch is never encountered.
>>> I'm not sure what's going on,
>>>
>>> I just copied the few lines from
>>> https://github.com/scipy/scipy/commit/a31acbf into my scipy 0.9 and
>>> the original example stops and I'm not able to produce an endless loop
>>> anymore when I try to change around with any of the examples, even
>>> when I start with a nan, it stops immediately.  I only tried the one
>>> parameter example.
>> Nope, still endless in linesearch with -np.exp(-x/2.) -np.exp(-x**2)
>> the third iteration starts with nan as parameter and never returns
>> from the line search.
>>
>> So something like your original fix to exit a nan linesearch looks necessary.
> In my example the linesearch is getting pk=np.nan
>
> So I tried the fix still in the optimization loop, at the end of the loop after
>
>          Hk = numpy.dot(A1,numpy.dot(Hk,A2)) + rhok * sk[:,numpy.newaxis] \
>                   * sk[numpy.newaxis,:]
>
> I added the pk check
> #add this
>          pk = -numpy.dot(Hk,gfk)
>          if not numpy.isfinite(pk):
>              #make sure we will have a valid search direction
>              warnflag = 2
>              break
>
> Now my example stops with the previous iteration, when there is no new
> search direction in the next iteration loop.
> I'm duplicating now the pk fix, so it's not really a fix yet, but it
> seems to work.
>
> (I would be glad if bfgs is more robust so I dare to use it also.)
> bfgs is missing a max function evaluation limit, I build it into my
> objective function to avoid the endless loop.
>
> Josef
>
>> Josef
>>
>>
>>> Can you check that you actually run the code that has this fix?
>>>
>>> Josef
>>>
>>>
>>>> Johann
>>>>
>>>> On 10/24/2011 09:26 PM, josef.pktd@gmail.com wrote:
>>>>> On Mon, Oct 24, 2011 at 2:39 PM, Johann Cohen-Tanugi
>>>>> <johann.cohentanugi@gmail.com>    wrote:
>>>>>> Dear Josef
>>>>>> On 10/24/2011 07:58 PM, josef.pktd@gmail.com wrote:
>>>>>>> On Mon, Oct 24, 2011 at 1:56 PM,<josef.pktd@gmail.com>      wrote:
>>>>>>>> On Mon, Oct 24, 2011 at 1:50 PM, Johann Cohen-Tanugi
>>>>>>>> <johann.cohentanugi@gmail.com>      wrote:
>>>>>>>>> Hello,
>>>>>>>>> the OP is a colleague of mine and I looked quickly at the code. The
>>>>>>>>> infinite
>>>>>>>>> loop in the OP's illustrating script comes from the "while 1" loop in
>>>>>>>>> l.144
>>>>>>>>> of linesearch.py : becuse phi0 is np.nan, phi1 is returned as np.nan
>>>>>>>>> as
>>>>>>>>> well, and the break condition is never met. There is an easy fix :
>>>>>>>>>      while 1:
>>>>>>>>>          stp, phi1, derphi1, task = minpack2.dcsrch(alpha1, phi1,
>>>>>>>>> derphi1,
>>>>>>>>>                                                     c1, c2, xtol, task,
>>>>>>>>>                                                     amin, amax, isave,
>>>>>>>>> dsave)
>>>>>>>>>          if task[:2] == asbytes('FG') and not np.isnan(phi1):
>>>>>>>>>              alpha1 = stp
>>>>>>>>>              phi1 = phi(stp)
>>>>>>>>>              derphi1 = derphi(stp)
>>>>>>>>>          else:
>>>>>>>>>              break
>>>>>>>>>
>>>>>>>>> but it is not a nice kludge.... Is there a better way to secure this
>>>>>>>>> while 1
>>>>>>>>> loop? I am sure I am not covering all possible pathological cases with
>>>>>>>>> adding "not np.isnan(phi1)" in the code above.
>>>>>>>> Is this still a problem with 0.10 ?
>>>>>>>> I thought this fixed it, https://github.com/scipy/scipy/commit/a31acbf
>>>>>> Well I am a complete newby with github, but I think I went to the head of
>>>>>> master before testing, and the problem is still there. I can see the code
>>>>>> snippet from https://github.com/scipy/scipy/commit/a31acbf in my local
>>>>>> copy,
>>>>>> and this is testing against +/-inf, not nan. Changing the OP's code to
>>>>>> test
>>>>>> against inf yields :
>>>>>> In [1]: run test_bgfs.py
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/optimize.py:303:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   if (max(numpy.ravel(abs(sim[1:] - sim[0])))<= xtol \
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/optimize.py:308:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   xr = (1 + rho)*xbar - rho*sim[-1]
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/optimize.py:350:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   sim[j] = sim[0] + sigma*(sim[j] - sim[0])
>>>>>> fmin works [ inf]
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/linesearch.py:132:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   alpha1 = min(1.0, 1.01*2*(phi0 - old_phi0)/derphi0)
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/linesearch.py:308:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   alpha1 = min(1.0, 1.01*2*(phi0 - old_phi0)/derphi0)
>>>>>>
>>>>>> /home/cohen/sources/python/scipydev/lib/python2.6/site-packages/scipy/optimize/linesearch.py:417:
>>>>>> RuntimeWarning: invalid value encountered in subtract
>>>>>>   B = (fb-D-C*db)/(db*db)
>>>>>> fmin_bfgs gets stuck in a loop [ nan]
>>>>>>
>>>>>> so it looks like your code solves the inf situation, but not the nan.
>>>>> It's not my fix (I'm still on scipy 0.9 and avoid bfgs because I don't
>>>>> want to have to kill my interpreter session)
>>>>>
>>>>> isfinite also checks for nans
>>>>>
>>>>>>>> np.isfinite(np.nan)
>>>>> False
>>>>>
>>>>> so there should be another reason that the linesearch doesn't return.
>>>>>
>>>>> Josef
>>>>>
>>>>>
>>>>>
>>>>>
>>>>>
>>>>>>> Is http://projects.scipy.org/scipy/ticket/1542 the same?
>>>>>> yes it looks like a duplicate
>>>>>>> josef
>>>>>>>
>>>>>>>> Josef
>>>>>>>>
>>>>>>>>
>>>>>>>>> thoughts?
>>>>>>>>> Johann
>>>>>>>>>
>>>>>>>>> On 08/14/2011 01:38 AM, b9o2jnbm tsd71eam wrote:
>>>>>>>>>
>>>>>>>>> I have run into a frustrating problem where scipy.optimize.fmin_bfgs
>>>>>>>>> will
>>>>>>>>> get stuck in an infinite loop.
>>>>>>>>>
>>>>>>>>> I have submitted a bug report:
>>>>>>>>>
>>>>>>>>> http://projects.scipy.org/scipy/ticket/1494
>>>>>>>>>
>>>>>>>>> but would also like to see if anybody on this list has any suggestions
>>>>>>>>> or
>>>>>>>>> feedback.
>>>>>>>>>
>>>>>>>>> Thanks,
>>>>>>>>>
>>>>>>>>> --
>>>>>>>>> This message has been scanned for viruses and
>>>>>>>>> dangerous content by MailScanner, and is
>>>>>>>>> believed to be clean.
>>>>>>>>>
>>>>>>>>> _______________________________________________
>>>>>>>>> SciPy-User mailing list
>>>>>>>>> SciPy-User@scipy.org
>>>>>>>>> http://mail.scipy.org/mailman/listinfo/scipy-user
>>>>>>>>>
>>>>>>>>> _______________________________________________
>>>>>>>>> SciPy-User mailing list
>>>>>>>>> SciPy-User@scipy.org
>>>>>>>>> http://mail.scipy.org/mailman/listinfo/scipy-user
>>>>>>>>>
>>>>>>>>>
> _______________________________________________
> SciPy-User mailing list
> SciPy-User@scipy.org
> http://mail.scipy.org/mailman/listinfo/scipy-user
>


More information about the SciPy-User mailing list