Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxlib / inputrec.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-2010, The GROMACS development team.
6  * Copyright (c) 2012,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 "gromacs/legacyheaders/inputrec.h"
40
41 #include "gromacs/legacyheaders/macros.h"
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/utility/fatalerror.h"
44
45
46 /* The minimum number of integration steps required for reasonably accurate
47  * integration of first and second order coupling algorithms.
48  */
49 const int nstmin_berendsen_tcouple =  5;
50 const int nstmin_berendsen_pcouple = 10;
51 const int nstmin_harmonic          = 20;
52
53 static int nst_wanted(const t_inputrec *ir)
54 {
55     if (ir->nstlist > 0)
56     {
57         return ir->nstlist;
58     }
59     else
60     {
61         return 10;
62     }
63 }
64
65 int ir_optimal_nstcalcenergy(const t_inputrec *ir)
66 {
67     return nst_wanted(ir);
68 }
69
70 int tcouple_min_integration_steps(int etc)
71 {
72     int n;
73
74     switch (etc)
75     {
76         case etcNO:
77             n = 0;
78             break;
79         case etcBERENDSEN:
80         case etcYES:
81             n = nstmin_berendsen_tcouple;
82             break;
83         case etcVRESCALE:
84             /* V-rescale supports instantaneous rescaling */
85             n = 0;
86             break;
87         case etcNOSEHOOVER:
88             n = nstmin_harmonic;
89             break;
90         case etcANDERSEN:
91         case etcANDERSENMASSIVE:
92             n = 1;
93             break;
94         default:
95             gmx_incons("Unknown etc value");
96             n = 0;
97     }
98
99     return n;
100 }
101
102 int ir_optimal_nsttcouple(const t_inputrec *ir)
103 {
104     int  nmin, nwanted, n;
105     real tau_min;
106     int  g;
107
108     nmin = tcouple_min_integration_steps(ir->etc);
109
110     nwanted = nst_wanted(ir);
111
112     tau_min = 1e20;
113     if (ir->etc != etcNO)
114     {
115         for (g = 0; g < ir->opts.ngtc; g++)
116         {
117             if (ir->opts.tau_t[g] > 0)
118             {
119                 tau_min = min(tau_min, ir->opts.tau_t[g]);
120             }
121         }
122     }
123
124     if (nmin == 0 || ir->delta_t*nwanted <= tau_min)
125     {
126         n = nwanted;
127     }
128     else
129     {
130         n = (int)(tau_min/(ir->delta_t*nmin) + 0.001);
131         if (n < 1)
132         {
133             n = 1;
134         }
135         while (nwanted % n != 0)
136         {
137             n--;
138         }
139     }
140
141     return n;
142 }
143
144 int pcouple_min_integration_steps(int epc)
145 {
146     int n;
147
148     switch (epc)
149     {
150         case epcNO:
151             n = 0;
152             break;
153         case etcBERENDSEN:
154         case epcISOTROPIC:
155             n = nstmin_berendsen_pcouple;
156             break;
157         case epcPARRINELLORAHMAN:
158         case epcMTTK:
159             n = nstmin_harmonic;
160             break;
161         default:
162             gmx_incons("Unknown epc value");
163             n = 0;
164     }
165
166     return n;
167 }
168
169 int ir_optimal_nstpcouple(const t_inputrec *ir)
170 {
171     int  nmin, nwanted, n;
172
173     nmin = pcouple_min_integration_steps(ir->epc);
174
175     nwanted = nst_wanted(ir);
176
177     if (nmin == 0 || ir->delta_t*nwanted <= ir->tau_p)
178     {
179         n = nwanted;
180     }
181     else
182     {
183         n = (int)(ir->tau_p/(ir->delta_t*nmin) + 0.001);
184         if (n < 1)
185         {
186             n = 1;
187         }
188         while (nwanted % n != 0)
189         {
190             n--;
191         }
192     }
193
194     return n;
195 }
196
197 gmx_bool ir_coulomb_switched(const t_inputrec *ir)
198 {
199     return (ir->coulombtype == eelSWITCH ||
200             ir->coulombtype == eelSHIFT ||
201             ir->coulombtype == eelENCADSHIFT ||
202             ir->coulombtype == eelPMESWITCH ||
203             ir->coulombtype == eelPMEUSERSWITCH ||
204             ir->coulomb_modifier == eintmodPOTSWITCH ||
205             ir->coulomb_modifier == eintmodFORCESWITCH);
206 }
207
208 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec *ir)
209 {
210     return (ir->cutoff_scheme == ecutsVERLET ||
211             ir_coulomb_switched(ir) || ir->coulomb_modifier != eintmodNONE ||
212             ir->coulombtype == eelRF_ZERO);
213 }
214
215 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec *ir)
216 {
217     return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == eelUSER || ir->coulombtype == eelPMEUSER);
218 }
219
220 gmx_bool ir_vdw_switched(const t_inputrec *ir)
221 {
222     return (ir->vdwtype == evdwSWITCH ||
223             ir->vdwtype == evdwSHIFT ||
224             ir->vdwtype == evdwENCADSHIFT ||
225             ir->vdw_modifier == eintmodPOTSWITCH ||
226             ir->vdw_modifier == eintmodFORCESWITCH);
227 }
228
229 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec *ir)
230 {
231     return (ir->cutoff_scheme == ecutsVERLET ||
232             ir_vdw_switched(ir) || ir->vdw_modifier != eintmodNONE);
233 }
234
235 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec *ir)
236 {
237     return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == evdwUSER);
238 }