3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * GRoups of Organic Molecules in ACtion for Science
40 #include "nb_verlet.h"
41 #include "interaction_const.h"
48 /* Abstract type for PME that is defined only in the routine that use them. */
49 typedef struct gmx_pme *gmx_pme_t;
52 real r; /* range of the table */
53 int n; /* n+1 is the number of points */
54 real scale; /* distance between two points */
55 real scale_exp; /* distance for exponential Buckingham table */
56 real *tab; /* the actual tables, per point there are 4 numbers for
57 * Coulomb, dispersion and repulsion (in total 12 numbers)
63 /* We duplicate tables for cache optimization purposes */
64 real *coultab; /* Coul only */
65 real *vdwtab; /* Vdw only */
66 /* The actual neighbor lists, short and long range, see enum above
67 * for definition of neighborlist indices.
69 t_nblist nlist_sr[eNL_NR];
70 t_nblist nlist_lr[eNL_NR];
73 /* macros for the cginfo data in forcerec */
74 /* The maximum cg size in cginfo is 63
75 * because we only have space for 6 bits in cginfo,
76 * this cg size entry is actually only read with domain decomposition.
77 * But there is a smaller limit due to the t_excl data structure
78 * which is defined in nblist.h.
80 #define SET_CGINFO_GID(cgi,gid) (cgi) = (((cgi) & ~65535) | (gid) )
81 #define GET_CGINFO_GID(cgi) ( (cgi) & 65535)
82 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
83 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
84 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
85 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
86 #define SET_CGINFO_SOLOPT(cgi,opt) (cgi) = (((cgi) & ~(3<<18)) | ((opt)<<18))
87 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 3)
88 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1<<20))
89 #define GET_CGINFO_CONSTR(cgi) ( (cgi) & (1<<20))
90 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1<<21))
91 #define GET_CGINFO_SETTLE(cgi) ( (cgi) & (1<<21))
92 /* This bit is only used with bBondComm in the domain decomposition */
93 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1<<22))
94 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi) & (1<<22))
95 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1<<23))
96 #define GET_CGINFO_HAS_VDW(cgi) ( (cgi) & (1<<23))
97 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1<<24))
98 #define GET_CGINFO_HAS_Q(cgi) ( (cgi) & (1<<24))
99 #define SET_CGINFO_NATOMS(cgi,opt) (cgi) = (((cgi) & ~(63<<25)) | ((opt)<<25))
100 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>25) & 63)
103 /* Value to be used in mdrun for an infinite cut-off.
104 * Since we need to compare with the cut-off squared,
105 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
107 #define GMX_CUTOFF_INF 1E+18
109 /* enums for the neighborlist type */
110 enum { enbvdwNONE,enbvdwLJ,enbvdwBHAM,enbvdwTAB,enbvdwNR};
111 /* OOR is "one over r" -- standard coul */
112 enum { enbcoulNONE,enbcoulOOR,enbcoulRF,enbcoulTAB,enbcoulGB,enbcoulFEWALD,enbcoulNR};
114 enum { egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
115 egCOUL14, egLJ14, egGB, egNR };
118 int nener; /* The number of energy group pairs */
119 real *ener[egNR]; /* Energy terms for each pair of groups */
123 real term[F_NRE]; /* The energies for all different interaction types */
124 gmx_grppairener_t grpp;
125 double dvdl_lin[efptNR]; /* Contributions to dvdl with linear lam-dependence */
126 double dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence */
128 int fep_state; /*current fep state -- just for printing */
129 double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
131 /* The idea is that dvdl terms with linear lambda dependence will be added
132 * automatically to enerpart_lambda. Terms with non-linear lambda dependence
133 * should explicitly determine the energies at foreign lambda points
145 /* ewald table type */
146 typedef struct ewald_tab *ewald_tab_t;
151 unsigned red_mask; /* Mask for marking which parts of f are filled */
154 gmx_grppairener_t grpp;
161 interaction_const_t *ic;
163 /* Domain Decomposition */
173 gmx_hw_info_t *hwinfo;
174 gmx_bool use_cpu_acceleration;
176 /* Use special N*N kernels? */
178 /* Private work data */
180 void *AllvsAll_workgb;
183 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
185 real rlist,rlistlong;
187 /* Dielectric constant resp. multiplication factor for charges */
189 real epsilon_r,epsilon_rf,epsfac;
191 /* Constants for reaction fields */
192 real kappa,k_rf,c_rf;
194 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
199 /* Dispersion correction stuff */
201 /* The shift of the shift or user potentials */
203 real enershifttwelve;
204 /* Integrated differces for energy and virial with cut-off functions */
209 /* Constant for long range dispersion correction (average dispersion)
210 * for topology A/B ([0]/[1]) */
212 /* Constant for long range repulsion term. Relative difference of about
213 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
223 /* The normal tables are in the nblists struct(s) below */
224 t_forcetable tab14; /* for 1-4 interactions only */
226 /* PPPM & Shifting stuff */
227 real rcoulomb_switch,rcoulomb;
232 real rvdw_switch,rvdw;
249 /* solvent_opt contains the enum for the most common solvent
250 * in the system, which will be optimized.
251 * It can be set to esolNO to disable all water optimization */
255 gmx_bool bExcl_IntraCGAll_InterCGNone;
256 cginfo_mb_t *cginfo_mb;
262 /* The neighborlists including tables */
267 int cutoff_scheme; /* old- or Verlet-style cutoff */
268 gmx_bool bNonbonded; /* true if nonbonded calculations are *not* turned off */
269 nonbonded_verlet_t *nbv;
271 /* The wall tables (if used) */
273 t_forcetable **wall_tab;
275 /* The number of charge groups participating in do_force_lowlevel */
277 /* The number of atoms participating in do_force_lowlevel */
279 /* The number of atoms participating in force and constraints */
280 int natoms_force_constr;
281 /* The allocation size of vectors of size natoms_force */
284 /* Twin Range stuff, f_twin has size natoms_force */
289 /* Forces that should not enter into the virial summation:
290 * PPPM/PME/Ewald/posres
292 gmx_bool bF_NoVirSum;
294 int f_novirsum_nalloc;
295 rvec *f_novirsum_alloc;
296 /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
297 * points to the normal force vectors wen pressure is not requested.
301 /* Long-range forces and virial for PPPM/PME/Ewald */
305 /* PME/Ewald stuff */
308 ewald_tab_t ewald_table;
312 rvec vir_diag_posres;
315 /* Non bonded Parameter lists */
316 int ntype; /* Number of atom types */
320 /* Energy group pair flags */
323 /* xmdrun flexible constraints */
326 /* Generalized born implicit solvent */
328 /* Generalized born stuff */
329 real gb_epsilon_solvent;
330 /* Table data for GB */
332 /* VdW radius for each atomtype (dim is thus ntype) */
334 /* Effective radius (derived from effective volume) for each type */
336 /* Implicit solvent - surface tension for each atomtype */
337 real *atype_surftens;
338 /* Implicit solvent - radius for GB calculation */
339 real *atype_gb_radius;
340 /* Implicit solvent - overlap for HCT model */
342 /* Generalized born interaction data */
345 /* Table scale for GB */
347 /* Table range for GB */
349 /* GB neighborlists (the sr list will contain for each atom all other atoms
350 * (for use in the SA calculation) and the lr list will contain
351 * for each atom all atoms 1-4 or greater (for use in the GB calculation)
357 /* Inverse square root of the Born radii for implicit solvent */
359 /* Derivatives of the potential with respect to the Born radii */
361 /* Derivatives of the Born radii with respect to coordinates */
364 int nalloc_dadx; /* Allocated size of dadx */
366 /* If > 0 signals Test Particle Insertion,
367 * the value is the number of atoms of the molecule to insert
368 * Only the energy difference due to the addition of the last molecule
369 * should be calculated.
373 /* Neighbor searching stuff */
380 /* QM-MM neighborlists */
383 /* Limit for printing large forces, negative is don't print */
386 /* coarse load balancing time measurement */
391 /* parameter needed for AdResS simulation */
393 gmx_bool badress_tf_full_box;
394 real adress_const_wf;
395 real adress_ex_width;
396 real adress_hy_width;
400 int n_adress_tf_grps;
401 int * adress_tf_table_index;
402 int *adress_group_explicit;
403 t_forcetable * atf_tabs;
404 real adress_ex_forcecap;
405 gmx_bool adress_do_hybridpairs;
407 /* User determined parameters, copied from the inputrec */
417 /* Thread local force and energy data */
418 /* FIXME move to bonded_thread_data_t */
424 /* Exclusion load distribution over the threads */
428 #define C6(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))]
429 #define C12(nbfp,ntp,ai,aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
430 #define BHAMC(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))]
431 #define BHAMA(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
432 #define BHAMB(nbfp,ntp,ai,aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]