Redefine the default boolean type to gmx_bool.
[gromacs.git] / include / types / forcerec.h
blob405b503bd0a7b39ee62af37ddc399387ef3ce914
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GRoups of Organic Molecules in ACtion for Science
36 #include "ns.h"
37 #include "genborn.h"
38 #include "qmmmrec.h"
39 #include "idef.h"
41 #ifdef __cplusplus
42 extern "C" {
43 #endif
45 /* Abstract type for PME that is defined only in the routine that use them. */
46 typedef struct gmx_pme *gmx_pme_t;
48 typedef struct {
49 real r; /* range of the table */
50 int n; /* n+1 is the number of points */
51 real scale; /* distance between two points */
52 real scale_exp; /* distance for exponential Buckingham table */
53 real *tab; /* the actual tables, per point there are 4 numbers for
54 * Coulomb, dispersion and repulsion (in total 12 numbers)
56 } t_forcetable;
58 typedef struct {
59 t_forcetable tab;
60 /* We duplicate tables for cache optimization purposes */
61 real *coultab; /* Coul only */
62 real *vdwtab; /* Vdw only */
63 /* The actual neighbor lists, short and long range, see enum above
64 * for definition of neighborlist indices.
66 t_nblist nlist_sr[eNL_NR];
67 t_nblist nlist_lr[eNL_NR];
68 } t_nblists;
70 /* macros for the cginfo data in forcerec */
71 /* The maximum cg size is 255, because we only have space for 8 bits in cginfo,
72 * this cg size entry is actually only read with domain decomposition.
74 #define MAX_CHARGEGROUP_SIZE 256
75 #define SET_CGINFO_GID(cgi,gid) (cgi) = (((cgi) & ~65535) | (gid) )
76 #define GET_CGINFO_GID(cgi) ( (cgi) & 65535)
77 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
78 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
79 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
80 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
81 #define SET_CGINFO_SOLOPT(cgi,opt) (cgi) = (((cgi) & ~(15<<18)) | ((opt)<<18))
82 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 15)
83 /* This bit is only used with bBondComm in the domain decomposition */
84 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1<<22))
85 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi) & (1<<22))
86 #define SET_CGINFO_NATOMS(cgi,opt) (cgi) = (((cgi) & ~(255<<23)) | ((opt)<<23))
87 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>23) & 255)
90 /* Value to be used in mdrun for an infinite cut-off.
91 * Since we need to compare with the cut-off squared,
92 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
94 #define GMX_CUTOFF_INF 1E+18
97 enum { egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
98 egCOUL14, egLJ14, egGB, egNR };
100 typedef struct {
101 int nener; /* The number of energy group pairs */
102 real *ener[egNR]; /* Energy terms for each pair of groups */
103 } gmx_grppairener_t;
105 typedef struct {
106 real term[F_NRE]; /* The energies for all different interaction types */
107 gmx_grppairener_t grpp;
108 double dvdl_lin; /* Contributions to dvdl with linear lam-dependence */
109 double dvdl_nonlin; /* Idem, but non-linear dependence */
110 int n_lambda;
111 double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
112 } gmx_enerdata_t;
113 /* The idea is that dvdl terms with linear lambda dependence will be added
114 * automatically to enerpart_lambda. Terms with non-linear lambda dependence
115 * should explicitly determine the energies at foreign lambda points
116 * when n_lambda > 0.
119 typedef struct {
120 int cg_start;
121 int cg_end;
122 int cg_mod;
123 int *cginfo;
124 } cginfo_mb_t;
127 /* ewald table type */
128 typedef struct ewald_tab *ewald_tab_t;
130 typedef struct {
131 /* Domain Decomposition */
132 gmx_bool bDomDec;
134 /* PBC stuff */
135 int ePBC;
136 gmx_bool bMolPBC;
137 int rc_scaling;
138 rvec posres_com;
139 rvec posres_comB;
141 gmx_bool UseOptimizedKernels;
143 /* Use special N*N kernels? */
144 gmx_bool bAllvsAll;
145 /* Private work data */
146 void *AllvsAll_work;
147 void *AllvsAll_workgb;
149 /* Cut-Off stuff.
150 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
152 real rlist,rlistlong;
154 /* Dielectric constant resp. multiplication factor for charges */
155 real zsquare,temp;
156 real epsilon_r,epsilon_rf,epsfac;
158 /* Constants for reaction fields */
159 real kappa,k_rf,c_rf;
161 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
162 double qsum[2];
163 rvec mu_tot[2];
165 /* Dispersion correction stuff */
166 int eDispCorr;
167 /* The shift of the shift or user potentials */
168 real enershiftsix;
169 real enershifttwelve;
170 /* Integrated differces for energy and virial with cut-off functions */
171 real enerdiffsix;
172 real enerdifftwelve;
173 real virdiffsix;
174 real virdifftwelve;
175 /* Constant for long range dispersion correction (average dispersion)
176 * for topology A/B ([0]/[1]) */
177 real avcsix[2];
178 /* Constant for long range repulsion term. Relative difference of about
179 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
181 real avctwelve[2];
183 /* Fudge factors */
184 real fudgeQQ;
186 /* Table stuff */
187 gmx_bool bcoultab;
188 gmx_bool bvdwtab;
189 /* The normal tables are in the nblists struct(s) below */
190 t_forcetable tab14; /* for 1-4 interactions only */
192 /* PPPM & Shifting stuff */
193 real rcoulomb_switch,rcoulomb;
194 real *phi;
196 /* VdW stuff */
197 double reppow;
198 real rvdw_switch,rvdw;
199 real bham_b_max;
201 /* Free energy ? */
202 int efep;
203 real sc_alpha;
204 int sc_power;
205 real sc_sigma6_def;
206 real sc_sigma6_min;
207 gmx_bool bSepDVDL;
209 /* NS Stuff */
210 int eeltype;
211 int vdwtype;
212 int cg0,hcg;
213 /* solvent_opt contains the enum for the most common solvent
214 * in the system, which will be optimized.
215 * It can be set to esolNO to disable all water optimization */
216 int solvent_opt;
217 int nWatMol;
218 gmx_bool bGrid;
219 gmx_bool bExcl_IntraCGAll_InterCGNone;
220 cginfo_mb_t *cginfo_mb;
221 int *cginfo;
222 rvec *cg_cm;
223 int cg_nalloc;
224 rvec *shift_vec;
226 /* The neighborlists including tables */
227 int nnblists;
228 int *gid2nblists;
229 t_nblists *nblists;
231 /* The wall tables (if used) */
232 int nwall;
233 t_forcetable **wall_tab;
235 /* This mask array of length nn determines whether or not this bit of the
236 * neighbourlists should be computed. Usually all these are true of course,
237 * but not when shells are used. During minimisation all the forces that
238 * include shells are done, then after minimsation is converged the remaining
239 * forces are computed.
241 /* gmx_bool *bMask; */
243 /* The number of charge groups participating in do_force_lowlevel */
244 int ncg_force;
245 /* The number of atoms participating in do_force_lowlevel */
246 int natoms_force;
247 /* The number of atoms participating in force and constraints */
248 int natoms_force_constr;
249 /* The allocation size of vectors of size natoms_force */
250 int nalloc_force;
252 /* Twin Range stuff, f_twin has size natoms_force */
253 gmx_bool bTwinRange;
254 int nlr;
255 rvec *f_twin;
257 /* Forces that should not enter into the virial summation:
258 * PPPM/PME/Ewald/posres
260 gmx_bool bF_NoVirSum;
261 int f_novirsum_n;
262 int f_novirsum_nalloc;
263 rvec *f_novirsum_alloc;
264 /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
265 * points to the normal force vectors wen pressure is not requested.
267 rvec *f_novirsum;
269 /* Long-range forces and virial for PPPM/PME/Ewald */
270 gmx_pme_t pmedata;
271 tensor vir_el_recip;
273 /* PME/Ewald stuff */
274 gmx_bool bEwald;
275 real ewaldcoeff;
276 ewald_tab_t ewald_table;
278 /* Virial Stuff */
279 rvec *fshift;
280 rvec vir_diag_posres;
281 dvec vir_wall_z;
283 /* Non bonded Parameter lists */
284 int ntype; /* Number of atom types */
285 gmx_bool bBHAM;
286 real *nbfp;
288 /* Energy group pair flags */
289 int *egp_flags;
291 /* xmdrun flexible constraints */
292 real fc_stepsize;
294 /* Generalized born implicit solvent */
295 gmx_bool bGB;
296 /* Generalized born stuff */
297 real gb_epsilon_solvent;
298 /* Table data for GB */
299 t_forcetable gbtab;
300 /* VdW radius for each atomtype (dim is thus ntype) */
301 real *atype_radius;
302 /* Effective radius (derived from effective volume) for each type */
303 real *atype_vol;
304 /* Implicit solvent - surface tension for each atomtype */
305 real *atype_surftens;
306 /* Implicit solvent - radius for GB calculation */
307 real *atype_gb_radius;
308 /* Implicit solvent - overlap for HCT model */
309 real *atype_S_hct;
310 /* Generalized born interaction data */
311 gmx_genborn_t *born;
313 /* Table scale for GB */
314 real gbtabscale;
315 /* Table range for GB */
316 real gbtabr;
317 /* GB neighborlists (the sr list will contain for each atom all other atoms
318 * (for use in the SA calculation) and the lr list will contain
319 * for each atom all atoms 1-4 or greater (for use in the GB calculation)
321 t_nblist gblist_sr;
322 t_nblist gblist_lr;
323 t_nblist gblist;
325 /* Inverse square root of the Born radii for implicit solvent */
326 real *invsqrta;
327 /* Derivatives of the potential with respect to the Born radii */
328 real *dvda;
329 /* Derivatives of the Born radii with respect to coordinates */
330 real *dadx;
331 real *dadx_rawptr;
332 int nalloc_dadx; /* Allocated size of dadx */
334 /* If > 0 signals Test Particle Insertion,
335 * the value is the number of atoms of the molecule to insert
336 * Only the energy difference due to the addition of the last molecule
337 * should be calculated.
339 gmx_bool n_tpi;
341 /* Neighbor searching stuff */
342 gmx_ns_t ns;
344 /* QMMM stuff */
345 gmx_bool bQMMM;
346 t_QMMMrec *qr;
348 /* QM-MM neighborlists */
349 t_nblist QMMMlist;
351 /* Limit for printing large forces, negative is don't print */
352 real print_force;
354 /* coarse load balancing time measurement */
355 double t_fnbf;
356 double t_wait;
357 int timesteps;
359 /* User determined parameters, copied from the inputrec */
360 int userint1;
361 int userint2;
362 int userint3;
363 int userint4;
364 real userreal1;
365 real userreal2;
366 real userreal3;
367 real userreal4;
368 } t_forcerec;
370 #define C6(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))]
371 #define C12(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
372 #define BHAMC(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))]
373 #define BHAMA(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
374 #define BHAMB(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]
376 #ifdef __cplusplus
378 #endif