Sort all includes in src/gromacs
[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
41 #include "gromacs/legacyheaders/md_logging.h"
42 #include "gromacs/legacyheaders/md_support.h"
43 #include "gromacs/legacyheaders/mdrun.h"
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/legacyheaders/types/iteratedconstraints.h"
46 #include "gromacs/utility/fatalerror.h"
47
48 #ifdef GMX_DOUBLE
49 #define CONVERGEITER  0.000000001
50 #define CLOSE_ENOUGH  0.000001000
51 #else
52 #define CONVERGEITER  0.0001
53 #define CLOSE_ENOUGH  0.0050
54 #endif
55
56 /* we want to keep track of the close calls.  If there are too many, there might be some other issues.
57    so we make sure that it's either less than some predetermined number, or if more than that number,
58    only some small fraction of the total. */
59 #define MAX_NUMBER_CLOSE        50
60 #define FRACTION_CLOSE       0.001
61
62 /* maximum length of cyclic traps to check, emerging from limited numerical precision  */
63 #define CYCLEMAX            20
64
65 void gmx_iterate_init(gmx_iterate_t *iterate, gmx_bool bSetIterationActive)
66 {
67     int i;
68
69     iterate->iter_i           = 0;
70     iterate->bIterationActive = bSetIterationActive;
71     iterate->num_close        = 0;
72     for (i = 0; i < MAXITERCONST+2; i++)
73     {
74         iterate->allrelerr[i] = 0;
75     }
76 }
77
78 gmx_bool done_iterating(const t_commrec *cr, FILE *fplog, int nsteps, gmx_iterate_t *iterate, gmx_bool bFirstIterate, real fom, real *newf)
79 {
80     /* monitor convergence, and use a secant search to propose new
81        values.
82                                                                   x_{i} - x_{i-1}
83        The secant method computes x_{i+1} = x_{i} - f(x_{i}) * ---------------------
84                                                                 f(x_{i}) - f(x_{i-1})
85
86        The function we are trying to zero is fom-x, where fom is the
87        "figure of merit" which is the pressure (or the veta value) we
88        would get by putting in an old value of the pressure or veta into
89        the incrementor function for the step or half step.  I have
90        verified that this gives the same answer as self consistent
91        iteration, usually in many fewer steps, especially for small tau_p.
92
93        We could possibly eliminate an iteration with proper use
94        of the value from the previous step, but that would take a bit
95        more bookkeeping, especially for veta, since tests indicate the
96        function of veta on the last step is not sufficiently close to
97        guarantee convergence this step. This is
98        good enough for now.  On my tests, I could use tau_p down to
99        0.02, which is smaller that would ever be necessary in
100        practice. Generally, 3-5 iterations will be sufficient */
101
102     real     relerr, err, xmin;
103     int      i;
104     gmx_bool incycle;
105
106     if (bFirstIterate)
107     {
108         iterate->x     = fom;
109         iterate->f     = fom-iterate->x;
110         iterate->xprev = 0;
111         iterate->fprev = 0;
112         *newf          = fom;
113     }
114     else
115     {
116         iterate->f = fom-iterate->x; /* we want to zero this difference */
117         if ((iterate->iter_i > 1) && (iterate->iter_i < MAXITERCONST))
118         {
119             if (iterate->f == iterate->fprev)
120             {
121                 *newf = iterate->f;
122             }
123             else
124             {
125                 *newf = iterate->x - (iterate->x-iterate->xprev)*(iterate->f)/(iterate->f-iterate->fprev);
126             }
127         }
128         else
129         {
130             /* just use self-consistent iteration the first step to initialize, or
131                if it's not converging (which happens occasionally -- need to investigate why) */
132             *newf = fom;
133         }
134     }
135     /* Consider a slight shortcut allowing us to exit one sooner -- we check the
136        difference between the closest of x and xprev to the new
137        value. To be 100% certain, we should check the difference between
138        the last result, and the previous result, or
139
140        relerr = (fabs((x-xprev)/fom));
141
142        but this is pretty much never necessary under typical conditions.
143        Checking numerically, it seems to lead to almost exactly the same
144        trajectories, but there are small differences out a few decimal
145        places in the pressure, and eventually in the v_eta, but it could
146        save an interation.
147
148        if (fabs(*newf-x) < fabs(*newf - xprev)) { xmin = x;} else { xmin = xprev;}
149        relerr = (fabs((*newf-xmin) / *newf));
150      */
151
152     err    = fabs((iterate->f-iterate->fprev));
153     relerr = fabs(err/fom);
154
155     iterate->allrelerr[iterate->iter_i] = relerr;
156
157     if (iterate->iter_i > 0)
158     {
159         if (debug)
160         {
161             fprintf(debug, "Iterating NPT constraints: %6i %20.12f%14.6g%20.12f\n",
162                     iterate->iter_i, fom, relerr, *newf);
163         }
164
165         if ((relerr < CONVERGEITER) || (err < CONVERGEITER) || (fom == 0) || ((iterate->x == iterate->xprev) && iterate->iter_i > 1))
166         {
167             iterate->bIterationActive = FALSE;
168             if (debug)
169             {
170                 fprintf(debug, "Iterating NPT constraints: CONVERGED\n");
171             }
172             return TRUE;
173         }
174         if (iterate->iter_i > MAXITERCONST)
175         {
176             if (relerr < CLOSE_ENOUGH)
177             {
178                 incycle = FALSE;
179                 for (i = 1; i < CYCLEMAX; i++)
180                 {
181                     if ((iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-1]) &&
182                         (iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-(1+2*i)]))
183                     {
184                         incycle = TRUE;
185                         if (debug)
186                         {
187                             fprintf(debug, "Exiting from an NPT iterating cycle of length %d\n", i);
188                         }
189                         break;
190                     }
191                 }
192
193                 if (incycle)
194                 {
195                     /* step 1: trapped in a numerical attractor */
196                     /* we are trapped in a numerical attractor, and can't converge any more, and are close to the final result.
197                        Better to give up convergence here than have the simulation die.
198                      */
199                     iterate->num_close++;
200                     iterate->bIterationActive = FALSE;
201                     return TRUE;
202                 }
203                 else
204                 {
205                     /* Step #2: test if we are reasonably close for other reasons, then monitor the number.  If not, die */
206
207                     /* how many close calls have we had?  If less than a few, we're OK */
208                     if (iterate->num_close < MAX_NUMBER_CLOSE)
209                     {
210                         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);
211                         iterate->num_close++;
212                         iterate->bIterationActive = FALSE;
213                         return TRUE;
214                         /* if more than a few, check the total fraction.  If too high, die. */
215                     }
216                     else if (iterate->num_close/(double)nsteps > FRACTION_CLOSE)
217                     {
218                         gmx_fatal(FARGS, "Could not converge NPT constraints, too many exceptions (%d%%\n", iterate->num_close/(double)nsteps);
219                     }
220                 }
221             }
222             else
223             {
224                 gmx_fatal(FARGS, "Could not converge NPT constraints\n");
225             }
226         }
227     }
228
229     iterate->xprev = iterate->x;
230     iterate->x     = *newf;
231     iterate->fprev = iterate->f;
232     iterate->iter_i++;
233
234     return FALSE;
235 }