Cleans up steepest descent initial potential hack
authorJohn Baxter <automata@unm.edu>
Thu, 22 Jan 2015 00:27:50 +0000 (17:27 -0700)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 23 Jan 2015 18:28:08 +0000 (19:28 +0100)
Addressed Bug #1676, in which there was a hack to
the do_steep() function in minimize.cpp where the
initial value of s_min->epot was set to
s_try->epot + 1; this patch ensures that the first
step's s_min->epot value is the true minimum
potential energy value and adds explicit checks
for the first step where necessary to ensure that
the behavior of do_steep() is unchanged.

Change-Id: I19c64a4948dfa24805345d7e593397d12ad4cc48

src/gromacs/mdlib/minimize.cpp

index a9c7500b64461f3c27b993c3fabb5f506dcddbfa..2a2bf8714ebded2a3588f1cc4598c82a9b9988bd 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -2479,7 +2479,7 @@ double do_steep(FILE *fplog, t_commrec *cr,
 
         if (count == 0)
         {
-            s_min->epot = s_try->epot + 1;
+            s_min->epot = s_try->epot;
         }
 
         /* Print it if necessary  */
@@ -2489,10 +2489,10 @@ double do_steep(FILE *fplog, t_commrec *cr,
             {
                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
-                        (s_try->epot < s_min->epot) ? '\n' : '\r');
+                        ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
             }
 
-            if (s_try->epot < s_min->epot)
+            if ( (count == 0) || (s_try->epot < s_min->epot) )
             {
                 /* Store the new (lower) energies  */
                 upd_mdebin(mdebin, FALSE, FALSE, (double)count,