Merge branch 'release-4-6', adds the nbnxn functionality
[alexxy/gromacs.git] / src / gromacs / mdlib / iteratedconstraints.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
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.
15
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.
20  * 
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.
27  * 
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.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41 #include "typedefs.h"
42 #include "gmx_fatal.h"
43 #include "mdrun.h"
44 #include "md_support.h"
45 #include "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 bIterate)
65 {
66     int i;
67
68     iterate->iter_i = 0;
69     iterate->bIterate = bIterate;
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->bIterate = 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                     if ((iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-1]) &&
180                         (iterate->allrelerr[iterate->iter_i-(1+i)] == iterate->allrelerr[iterate->iter_i-(1+2*i)])) {
181                         incycle = TRUE;
182                         if (debug) 
183                         {
184                             fprintf(debug,"Exiting from an NPT iterating cycle of length %d\n",i);
185                         }
186                         break;
187                     }
188                 }
189                 
190                 if (incycle) {
191                     /* step 1: trapped in a numerical attractor */
192                     /* we are trapped in a numerical attractor, and can't converge any more, and are close to the final result.
193                        Better to give up convergence here than have the simulation die.
194                     */
195                     iterate->num_close++;
196                     return TRUE;
197                 } 
198                 else 
199                 {
200                     /* Step #2: test if we are reasonably close for other reasons, then monitor the number.  If not, die */
201                     
202                     /* how many close calls have we had?  If less than a few, we're OK */
203                     if (iterate->num_close < MAX_NUMBER_CLOSE) 
204                     {
205                         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);
206                         iterate->num_close++;
207                         return TRUE;
208                         /* if more than a few, check the total fraction.  If too high, die. */
209                     } else if (iterate->num_close/(double)nsteps > FRACTION_CLOSE) {
210                         gmx_fatal(FARGS,"Could not converge NPT constraints, too many exceptions (%d%%\n",iterate->num_close/(double)nsteps);
211                     } 
212                 }
213             }
214             else 
215             {
216                 gmx_fatal(FARGS,"Could not converge NPT constraints\n");
217             }
218         }
219     }
220     
221     iterate->xprev = iterate->x;
222     iterate->x = *newf;
223     iterate->fprev = iterate->f;
224     iterate->iter_i++;
225     
226     return FALSE;
227 }
228