2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "gmx_fatal.h"
46 #include "md_support.h"
47 #include "md_logging.h"
48 #include "types/iteratedconstraints.h"
51 #define CONVERGEITER 0.000000001
52 #define CLOSE_ENOUGH 0.000001000
54 #define CONVERGEITER 0.0001
55 #define CLOSE_ENOUGH 0.0050
58 /* we want to keep track of the close calls. If there are too many, there might be some other issues.
59 so we make sure that it's either less than some predetermined number, or if more than that number,
60 only some small fraction of the total. */
61 #define MAX_NUMBER_CLOSE 50
62 #define FRACTION_CLOSE 0.001
64 /* maximum length of cyclic traps to check, emerging from limited numerical precision */
67 void gmx_iterate_init(gmx_iterate_t *iterate,gmx_bool bIterate)
72 iterate->bIterate = bIterate;
73 iterate->num_close = 0;
74 for (i=0;i<MAXITERCONST+2;i++)
76 iterate->allrelerr[i] = 0;
80 gmx_bool done_iterating(const t_commrec *cr,FILE *fplog, int nsteps, gmx_iterate_t *iterate, gmx_bool bFirstIterate, real fom, real *newf)
82 /* monitor convergence, and use a secant search to propose new
85 The secant method computes x_{i+1} = x_{i} - f(x_{i}) * ---------------------
88 The function we are trying to zero is fom-x, where fom is the
89 "figure of merit" which is the pressure (or the veta value) we
90 would get by putting in an old value of the pressure or veta into
91 the incrementor function for the step or half step. I have
92 verified that this gives the same answer as self consistent
93 iteration, usually in many fewer steps, especially for small tau_p.
95 We could possibly eliminate an iteration with proper use
96 of the value from the previous step, but that would take a bit
97 more bookkeeping, especially for veta, since tests indicate the
98 function of veta on the last step is not sufficiently close to
99 guarantee convergence this step. This is
100 good enough for now. On my tests, I could use tau_p down to
101 0.02, which is smaller that would ever be necessary in
102 practice. Generally, 3-5 iterations will be sufficient */
104 real relerr,err,xmin;
111 iterate->f = fom-iterate->x;
118 iterate->f = fom-iterate->x; /* we want to zero this difference */
119 if ((iterate->iter_i > 1) && (iterate->iter_i < MAXITERCONST))
121 if (iterate->f==iterate->fprev)
127 *newf = iterate->x - (iterate->x-iterate->xprev)*(iterate->f)/(iterate->f-iterate->fprev);
132 /* just use self-consistent iteration the first step to initialize, or
133 if it's not converging (which happens occasionally -- need to investigate why) */
137 /* Consider a slight shortcut allowing us to exit one sooner -- we check the
138 difference between the closest of x and xprev to the new
139 value. To be 100% certain, we should check the difference between
140 the last result, and the previous result, or
142 relerr = (fabs((x-xprev)/fom));
144 but this is pretty much never necessary under typical conditions.
145 Checking numerically, it seems to lead to almost exactly the same
146 trajectories, but there are small differences out a few decimal
147 places in the pressure, and eventually in the v_eta, but it could
150 if (fabs(*newf-x) < fabs(*newf - xprev)) { xmin = x;} else { xmin = xprev;}
151 relerr = (fabs((*newf-xmin) / *newf));
154 err = fabs((iterate->f-iterate->fprev));
155 relerr = fabs(err/fom);
157 iterate->allrelerr[iterate->iter_i] = relerr;
159 if (iterate->iter_i > 0)
163 fprintf(debug,"Iterating NPT constraints: %6i %20.12f%14.6g%20.12f\n",
164 iterate->iter_i,fom,relerr,*newf);
167 if ((relerr < CONVERGEITER) || (err < CONVERGEITER) || (fom==0) || ((iterate->x == iterate->xprev) && iterate->iter_i > 1))
169 iterate->bIterate = FALSE;
172 fprintf(debug,"Iterating NPT constraints: CONVERGED\n");
176 if (iterate->iter_i > MAXITERCONST)
178 if (relerr < CLOSE_ENOUGH)
181 for (i=1;i<CYCLEMAX;i++) {
182 if ((iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-1]) &&
183 (iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-(1+2*i)])) {
187 fprintf(debug,"Exiting from an NPT iterating cycle of length %d\n",i);
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.
198 iterate->num_close++;
203 /* Step #2: test if we are reasonably close for other reasons, then monitor the number. If not, die */
205 /* how many close calls have we had? If less than a few, we're OK */
206 if (iterate->num_close < MAX_NUMBER_CLOSE)
208 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);
209 iterate->num_close++;
211 /* if more than a few, check the total fraction. If too high, die. */
212 } else if (iterate->num_close/(double)nsteps > FRACTION_CLOSE) {
213 gmx_fatal(FARGS,"Could not converge NPT constraints, too many exceptions (%d%%\n",iterate->num_close/(double)nsteps);
219 gmx_fatal(FARGS,"Could not converge NPT constraints\n");
224 iterate->xprev = iterate->x;
226 iterate->fprev = iterate->f;