1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
40 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
46 #include "gmx_fatal.h"
50 #define CONVERGEITER 0.000000001
51 #define CLOSE_ENOUGH 0.000001000
53 #define CONVERGEITER 0.0001
54 #define CLOSE_ENOUGH 0.0050
57 /* we want to keep track of the close calls. If there are too many, there might be some other issues.
58 so we make sure that it's either less than some predetermined number, or if more than that number,
59 only some small fraction of the total. */
60 #define MAX_NUMBER_CLOSE 50
61 #define FRACTION_CLOSE 0.001
63 /* maximum length of cyclic traps to check, emerging from limited numerical precision */
66 void gmx_iterate_init(gmx_iterate_t *iterate,gmx_bool bIterate)
71 iterate->bIterate = bIterate;
72 iterate->num_close = 0;
73 for (i=0;i<MAXITERCONST+2;i++)
75 iterate->allrelerr[i] = 0;
79 gmx_bool done_iterating(const t_commrec *cr,FILE *fplog, int nsteps, gmx_iterate_t *iterate, gmx_bool bFirstIterate, real fom, real *newf)
81 /* monitor convergence, and use a secant search to propose new
84 The secant method computes x_{i+1} = x_{i} - f(x_{i}) * ---------------------
87 The function we are trying to zero is fom-x, where fom is the
88 "figure of merit" which is the pressure (or the veta value) we
89 would get by putting in an old value of the pressure or veta into
90 the incrementor function for the step or half step. I have
91 verified that this gives the same answer as self consistent
92 iteration, usually in many fewer steps, especially for small tau_p.
94 We could possibly eliminate an iteration with proper use
95 of the value from the previous step, but that would take a bit
96 more bookkeeping, especially for veta, since tests indicate the
97 function of veta on the last step is not sufficiently close to
98 guarantee convergence this step. This is
99 good enough for now. On my tests, I could use tau_p down to
100 0.02, which is smaller that would ever be necessary in
101 practice. Generally, 3-5 iterations will be sufficient */
103 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 sprintf(buf,"Slight numerical convergence deviation with NPT at step %d, relative error only %10.5g, likely not a problem, continuing\n",nsteps,relerr);
209 md_print_warning(cr,fplog,buf);
210 iterate->num_close++;
212 /* if more than a few, check the total fraction. If too high, die. */
213 } else if (iterate->num_close/(double)nsteps > FRACTION_CLOSE) {
214 gmx_fatal(FARGS,"Could not converge NPT constraints, too many exceptions (%d%%\n",iterate->num_close/(double)nsteps);
220 gmx_fatal(FARGS,"Could not converge NPT constraints\n");
225 iterate->xprev = iterate->x;
227 iterate->fprev = iterate->f;