Remove unnecessary config.h includes
[alexxy/gromacs.git] / src / gromacs / mdlib / iteratedconstraints.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include <math.h>
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/utility/fatalerror.h"
42 #include "gromacs/legacyheaders/mdrun.h"
43 #include "gromacs/legacyheaders/md_support.h"
44 #include "gromacs/legacyheaders/md_logging.h"
45 #include "gromacs/legacyheaders/types/iteratedconstraints.h"
46
47 #ifdef GMX_DOUBLE
48 #define CONVERGEITER  0.000000001
49 #define CLOSE_ENOUGH  0.000001000
50 #else
51 #define CONVERGEITER  0.0001
52 #define CLOSE_ENOUGH  0.0050
53 #endif
54
55 /* we want to keep track of the close calls.  If there are too many, there might be some other issues.
56    so we make sure that it's either less than some predetermined number, or if more than that number,
57    only some small fraction of the total. */
58 #define MAX_NUMBER_CLOSE        50
59 #define FRACTION_CLOSE       0.001
60
61 /* maximum length of cyclic traps to check, emerging from limited numerical precision  */
62 #define CYCLEMAX            20
63
64 void gmx_iterate_init(gmx_iterate_t *iterate, gmx_bool bSetIterationActive)
65 {
66     int i;
67
68     iterate->iter_i           = 0;
69     iterate->bIterationActive = bSetIterationActive;
70     iterate->num_close        = 0;
71     for (i = 0; i < MAXITERCONST+2; i++)
72     {
73         iterate->allrelerr[i] = 0;
74     }
75 }
76
77 gmx_bool done_iterating(const t_commrec *cr, FILE *fplog, int nsteps, gmx_iterate_t *iterate, gmx_bool bFirstIterate, real fom, real *newf)
78 {
79     /* monitor convergence, and use a secant search to propose new
80        values.
81                                                                   x_{i} - x_{i-1}
82        The secant method computes x_{i+1} = x_{i} - f(x_{i}) * ---------------------
83                                                                 f(x_{i}) - f(x_{i-1})
84
85        The function we are trying to zero is fom-x, where fom is the
86        "figure of merit" which is the pressure (or the veta value) we
87        would get by putting in an old value of the pressure or veta into
88        the incrementor function for the step or half step.  I have
89        verified that this gives the same answer as self consistent
90        iteration, usually in many fewer steps, especially for small tau_p.
91
92        We could possibly eliminate an iteration with proper use
93        of the value from the previous step, but that would take a bit
94        more bookkeeping, especially for veta, since tests indicate the
95        function of veta on the last step is not sufficiently close to
96        guarantee convergence this step. This is
97        good enough for now.  On my tests, I could use tau_p down to
98        0.02, which is smaller that would ever be necessary in
99        practice. Generally, 3-5 iterations will be sufficient */
100
101     real     relerr, err, xmin;
102     int      i;
103     gmx_bool incycle;
104
105     if (bFirstIterate)
106     {
107         iterate->x     = fom;
108         iterate->f     = fom-iterate->x;
109         iterate->xprev = 0;
110         iterate->fprev = 0;
111         *newf          = fom;
112     }
113     else
114     {
115         iterate->f = fom-iterate->x; /* we want to zero this difference */
116         if ((iterate->iter_i > 1) && (iterate->iter_i < MAXITERCONST))
117         {
118             if (iterate->f == iterate->fprev)
119             {
120                 *newf = iterate->f;
121             }
122             else
123             {
124                 *newf = iterate->x - (iterate->x-iterate->xprev)*(iterate->f)/(iterate->f-iterate->fprev);
125             }
126         }
127         else
128         {
129             /* just use self-consistent iteration the first step to initialize, or
130                if it's not converging (which happens occasionally -- need to investigate why) */
131             *newf = fom;
132         }
133     }
134     /* Consider a slight shortcut allowing us to exit one sooner -- we check the
135        difference between the closest of x and xprev to the new
136        value. To be 100% certain, we should check the difference between
137        the last result, and the previous result, or
138
139        relerr = (fabs((x-xprev)/fom));
140
141        but this is pretty much never necessary under typical conditions.
142        Checking numerically, it seems to lead to almost exactly the same
143        trajectories, but there are small differences out a few decimal
144        places in the pressure, and eventually in the v_eta, but it could
145        save an interation.
146
147        if (fabs(*newf-x) < fabs(*newf - xprev)) { xmin = x;} else { xmin = xprev;}
148        relerr = (fabs((*newf-xmin) / *newf));
149      */
150
151     err    = fabs((iterate->f-iterate->fprev));
152     relerr = fabs(err/fom);
153
154     iterate->allrelerr[iterate->iter_i] = relerr;
155
156     if (iterate->iter_i > 0)
157     {
158         if (debug)
159         {
160             fprintf(debug, "Iterating NPT constraints: %6i %20.12f%14.6g%20.12f\n",
161                     iterate->iter_i, fom, relerr, *newf);
162         }
163
164         if ((relerr < CONVERGEITER) || (err < CONVERGEITER) || (fom == 0) || ((iterate->x == iterate->xprev) && iterate->iter_i > 1))
165         {
166             iterate->bIterationActive = FALSE;
167             if (debug)
168             {
169                 fprintf(debug, "Iterating NPT constraints: CONVERGED\n");
170             }
171             return TRUE;
172         }
173         if (iterate->iter_i > MAXITERCONST)
174         {
175             if (relerr < CLOSE_ENOUGH)
176             {
177                 incycle = FALSE;
178                 for (i = 1; i < CYCLEMAX; i++)
179                 {
180                     if ((iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-1]) &&
181                         (iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-(1+2*i)]))
182                     {
183                         incycle = TRUE;
184                         if (debug)
185                         {
186                             fprintf(debug, "Exiting from an NPT iterating cycle of length %d\n", i);
187                         }
188                         break;
189                     }
190                 }
191
192                 if (incycle)
193                 {
194                     /* step 1: trapped in a numerical attractor */
195                     /* we are trapped in a numerical attractor, and can't converge any more, and are close to the final result.
196                        Better to give up convergence here than have the simulation die.
197                      */
198                     iterate->num_close++;
199                     iterate->bIterationActive = FALSE;
200                     return TRUE;
201                 }
202                 else
203                 {
204                     /* Step #2: test if we are reasonably close for other reasons, then monitor the number.  If not, die */
205
206                     /* how many close calls have we had?  If less than a few, we're OK */
207                     if (iterate->num_close < MAX_NUMBER_CLOSE)
208                     {
209                         md_print_warn(cr, fplog, "Slight numerical convergence deviation with NPT at step %d, relative error only %10.5g, likely not a problem, continuing\n", nsteps, relerr);
210                         iterate->num_close++;
211                         iterate->bIterationActive = FALSE;
212                         return TRUE;
213                         /* if more than a few, check the total fraction.  If too high, die. */
214                     }
215                     else if (iterate->num_close/(double)nsteps > FRACTION_CLOSE)
216                     {
217                         gmx_fatal(FARGS, "Could not converge NPT constraints, too many exceptions (%d%%\n", iterate->num_close/(double)nsteps);
218                     }
219                 }
220             }
221             else
222             {
223                 gmx_fatal(FARGS, "Could not converge NPT constraints\n");
224             }
225         }
226     }
227
228     iterate->xprev = iterate->x;
229     iterate->x     = *newf;
230     iterate->fprev = iterate->f;
231     iterate->iter_i++;
232
233     return FALSE;
234 }