Make virial reproducible
[alexxy/gromacs.git] / src / gromacs / mdtypes / forcerec.h
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-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016, 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 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
38 #define GMX_MDTYPES_TYPES_FORCEREC_H
39
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/mdtypes/interaction_const.h"
42 #include "gromacs/mdtypes/md_enums.h"
43 #include "gromacs/topology/idef.h"
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/real.h"
46
47 #ifdef __cplusplus
48 extern "C" {
49 #endif
50 #if 0
51 } /* fixes auto-indentation problems */
52 #endif
53
54 /* Abstract type for PME that is defined only in the routine that use them. */
55 struct gmx_genborn_t;
56 struct gmx_ns_t;
57 struct gmx_pme_t;
58 struct nonbonded_verlet_t;
59 struct bonded_threading_t;
60 struct t_forcetable;
61 struct t_nblist;
62 struct t_nblists;
63 struct t_QMMMrec;
64 struct gmx_hw_info_t;
65 struct gmx_gpu_opt_t;
66
67 /* macros for the cginfo data in forcerec
68  *
69  * Since the tpx format support max 256 energy groups, we do the same here.
70  * Note that we thus have bits 8-14 still unused.
71  *
72  * The maximum cg size in cginfo is 63
73  * because we only have space for 6 bits in cginfo,
74  * this cg size entry is actually only read with domain decomposition.
75  * But there is a smaller limit due to the t_excl data structure
76  * which is defined in nblist.h.
77  */
78 #define SET_CGINFO_GID(cgi, gid)     (cgi) = (((cgi)  &  ~255) | (gid))
79 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   255)
80 #define SET_CGINFO_FEP(cgi)          (cgi) =  ((cgi)  |  (1<<15))
81 #define GET_CGINFO_FEP(cgi)        ( (cgi)            &  (1<<15))
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)
101
102
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).
106  */
107 #define GMX_CUTOFF_INF 1E+18
108
109 /* enums for the neighborlist type */
110 enum {
111     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
112 };
113 /* OOR is "one over r" -- standard coul */
114 enum {
115     enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
116 };
117
118 enum {
119     egCOULSR, egLJSR, egBHAMSR,
120     egCOUL14, egLJ14, egGB, egNR
121 };
122 extern const char *egrp_nm[egNR+1];
123
124 typedef struct gmx_grppairener_t {
125     int   nener;      /* The number of energy group pairs     */
126     real *ener[egNR]; /* Energy terms for each pair of groups */
127 } gmx_grppairener_t;
128
129 typedef struct gmx_enerdata_t {
130     real              term[F_NRE];         /* The energies for all different interaction types */
131     gmx_grppairener_t grpp;
132     double            dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
133     double            dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
134     int               n_lambda;
135     int               fep_state;           /*current fep state -- just for printing */
136     double           *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
137     real              foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
138     gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
139 } gmx_enerdata_t;
140 /* The idea is that dvdl terms with linear lambda dependence will be added
141  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
142  * should explicitly determine the energies at foreign lambda points
143  * when n_lambda > 0.
144  */
145
146 typedef struct {
147     int  cg_start;
148     int  cg_end;
149     int  cg_mod;
150     int *cginfo;
151 } cginfo_mb_t;
152
153
154 /* Forward declaration of type for managing Ewald tables */
155 struct gmx_ewald_tab_t;
156
157 typedef struct ewald_corr_thread_t ewald_corr_thread_t;
158
159 typedef struct t_forcerec {
160     interaction_const_t *ic;
161
162     /* Domain Decomposition */
163     gmx_bool bDomDec;
164
165     /* PBC stuff */
166     int                         ePBC;
167     gmx_bool                    bMolPBC;
168     int                         rc_scaling;
169     rvec                        posres_com;
170     rvec                        posres_comB;
171
172     const struct gmx_hw_info_t *hwinfo;
173     const struct gmx_gpu_opt_t *gpu_opt;
174     gmx_bool                    use_simd_kernels;
175
176     /* Interaction for calculated in kernels. In many cases this is similar to
177      * the electrostatics settings in the inputrecord, but the difference is that
178      * these variables always specify the actual interaction in the kernel - if
179      * we are tabulating reaction-field the inputrec will say reaction-field, but
180      * the kernel interaction will say cubic-spline-table. To be safe we also
181      * have a kernel-specific setting for the modifiers - if the interaction is
182      * tabulated we already included the inputrec modification there, so the kernel
183      * modification setting will say 'none' in that case.
184      */
185     int nbkernel_elec_interaction;
186     int nbkernel_vdw_interaction;
187     int nbkernel_elec_modifier;
188     int nbkernel_vdw_modifier;
189
190     /* Use special N*N kernels? */
191     gmx_bool bAllvsAll;
192     /* Private work data */
193     void    *AllvsAll_work;
194     void    *AllvsAll_workgb;
195
196     /* Cut-Off stuff.
197      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
198      */
199     real rlist;
200
201     /* Dielectric constant resp. multiplication factor for charges */
202     real zsquare, temp;
203     real epsilon_r, epsilon_rf, epsfac;
204
205     /* Constants for reaction fields */
206     real kappa, k_rf, c_rf;
207
208     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
209     double qsum[2];
210     double q2sum[2];
211     double c6sum[2];
212     rvec   mu_tot[2];
213
214     /* Dispersion correction stuff */
215     int                  eDispCorr;
216     int                  numAtomsForDispersionCorrection;
217     struct t_forcetable *dispersionCorrectionTable;
218
219     /* The shift of the shift or user potentials */
220     real enershiftsix;
221     real enershifttwelve;
222     /* Integrated differces for energy and virial with cut-off functions */
223     real enerdiffsix;
224     real enerdifftwelve;
225     real virdiffsix;
226     real virdifftwelve;
227     /* Constant for long range dispersion correction (average dispersion)
228      * for topology A/B ([0]/[1]) */
229     real avcsix[2];
230     /* Constant for long range repulsion term. Relative difference of about
231      * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
232      */
233     real avctwelve[2];
234
235     /* Fudge factors */
236     real fudgeQQ;
237
238     /* Table stuff */
239     gmx_bool             bcoultab;
240     gmx_bool             bvdwtab;
241     /* The normal tables are in the nblists struct(s) below */
242
243     struct t_forcetable *pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
244
245     /* PPPM & Shifting stuff */
246     int   coulomb_modifier;
247     real  rcoulomb_switch, rcoulomb;
248     real *phi;
249
250     /* VdW stuff */
251     int    vdw_modifier;
252     double reppow;
253     real   rvdw_switch, rvdw;
254     real   bham_b_max;
255
256     /* Free energy */
257     int      efep;
258     real     sc_alphavdw;
259     real     sc_alphacoul;
260     int      sc_power;
261     real     sc_r_power;
262     real     sc_sigma6_def;
263     real     sc_sigma6_min;
264
265     /* NS Stuff */
266     int  eeltype;
267     int  vdwtype;
268     int  cg0, hcg;
269     /* solvent_opt contains the enum for the most common solvent
270      * in the system, which will be optimized.
271      * It can be set to esolNO to disable all water optimization */
272     int          solvent_opt;
273     int          nWatMol;
274     gmx_bool     bGrid;
275     gmx_bool     bExcl_IntraCGAll_InterCGNone;
276     cginfo_mb_t *cginfo_mb;
277     int         *cginfo;
278     rvec        *cg_cm;
279     int          cg_nalloc;
280     rvec        *shift_vec;
281
282     /* The neighborlists including tables */
283     int                        nnblists;
284     int                       *gid2nblists;
285     struct t_nblists          *nblists;
286
287     int                        cutoff_scheme; /* group- or Verlet-style cutoff */
288     gmx_bool                   bNonbonded;    /* true if nonbonded calculations are *not* turned off */
289     struct nonbonded_verlet_t *nbv;
290
291     /* The wall tables (if used) */
292     int                    nwall;
293     struct t_forcetable ***wall_tab;
294
295     /* The number of charge groups participating in do_force_lowlevel */
296     int ncg_force;
297     /* The number of atoms participating in do_force_lowlevel */
298     int natoms_force;
299     /* The number of atoms participating in force and constraints */
300     int natoms_force_constr;
301     /* The allocation size of vectors of size natoms_force */
302     int nalloc_force;
303
304     /* Forces that should not enter into the virial summation:
305      * PPPM/PME/Ewald/posres
306      */
307     gmx_bool bF_NoVirSum;
308     int      f_novirsum_n;
309     int      f_novirsum_nalloc;
310     rvec    *f_novirsum_alloc;
311     /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
312      * points to the normal force vectors wen pressure is not requested.
313      */
314     rvec *f_novirsum;
315
316     /* Long-range forces and virial for PPPM/PME/Ewald */
317     struct gmx_pme_t *pmedata;
318     int               ljpme_combination_rule;
319     tensor            vir_el_recip;
320     tensor            vir_lj_recip;
321
322     /* PME/Ewald stuff */
323     gmx_bool                bEwald;
324     real                    ewaldcoeff_q;
325     real                    ewaldcoeff_lj;
326     struct gmx_ewald_tab_t *ewald_table;
327
328     /* Virial Stuff */
329     rvec *fshift;
330     rvec  vir_diag_posres;
331     dvec  vir_wall_z;
332
333     /* Non bonded Parameter lists */
334     int      ntype; /* Number of atom types */
335     gmx_bool bBHAM;
336     real    *nbfp;
337     real    *ljpme_c6grid; /* C6-values used on grid in LJPME */
338
339     /* Energy group pair flags */
340     int *egp_flags;
341
342     /* Shell molecular dynamics flexible constraints */
343     real fc_stepsize;
344
345     /* Generalized born implicit solvent */
346     gmx_bool              bGB;
347     /* Generalized born stuff */
348     real                  gb_epsilon_solvent;
349     /* Table data for GB */
350     struct t_forcetable  *gbtab;
351     /* VdW radius for each atomtype (dim is thus ntype) */
352     real                 *atype_radius;
353     /* Effective radius (derived from effective volume) for each type */
354     real                 *atype_vol;
355     /* Implicit solvent - surface tension for each atomtype */
356     real                 *atype_surftens;
357     /* Implicit solvent - radius for GB calculation */
358     real                 *atype_gb_radius;
359     /* Implicit solvent - overlap for HCT model */
360     real                 *atype_S_hct;
361     /* Generalized born interaction data */
362     struct gmx_genborn_t *born;
363
364     /* Table scale for GB */
365     real gbtabscale;
366     /* Table range for GB */
367     real gbtabr;
368     /* GB neighborlists (the sr list will contain for each atom all other atoms
369      * (for use in the SA calculation) and the lr list will contain
370      * for each atom all atoms 1-4 or greater (for use in the GB calculation)
371      */
372     struct t_nblist *gblist_sr;
373     struct t_nblist *gblist_lr;
374     struct t_nblist *gblist;
375
376     /* Inverse square root of the Born radii for implicit solvent */
377     real *invsqrta;
378     /* Derivatives of the potential with respect to the Born radii */
379     real *dvda;
380     /* Derivatives of the Born radii with respect to coordinates */
381     real *dadx;
382     real *dadx_rawptr;
383     int   nalloc_dadx; /* Allocated size of dadx */
384
385     /* If > 0 signals Test Particle Insertion,
386      * the value is the number of atoms of the molecule to insert
387      * Only the energy difference due to the addition of the last molecule
388      * should be calculated.
389      */
390     gmx_bool n_tpi;
391
392     /* Neighbor searching stuff */
393     struct gmx_ns_t *ns;
394
395     /* QMMM stuff */
396     gmx_bool          bQMMM;
397     struct t_QMMMrec *qr;
398
399     /* QM-MM neighborlists */
400     struct t_nblist        *QMMMlist;
401
402     /* Limit for printing large forces, negative is don't print */
403     real print_force;
404
405     /* coarse load balancing time measurement */
406     double t_fnbf;
407     double t_wait;
408     int    timesteps;
409
410     /* User determined parameters, copied from the inputrec */
411     int  userint1;
412     int  userint2;
413     int  userint3;
414     int  userint4;
415     real userreal1;
416     real userreal2;
417     real userreal3;
418     real userreal4;
419
420     /* Pointer to struct for managing threading of bonded force calculation */
421     struct bonded_threading_t *bonded_threading;
422
423     /* Ewald correction thread local virial and energy data */
424     int                  nthread_ewc;
425     ewald_corr_thread_t *ewc_t;
426     /* Ewald charge correction load distribution over the threads */
427     int                 *excl_load;
428 } t_forcerec;
429
430 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
431  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
432  * in the code, but beware if you are using these macros externally.
433  */
434 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
435 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
436 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
437 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
438 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
439
440 #ifdef __cplusplus
441 }
442 #endif
443 #endif