Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / 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, 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
38 #include "ns.h"
39 #include "genborn.h"
40 #include "qmmmrec.h"
41 #include "../../topology/idef.h"
42 #include "nb_verlet.h"
43 #include "interaction_const.h"
44 #include "hw_info.h"
45
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
49 #if 0
50 } /* fixes auto-indentation problems */
51 #endif
52
53 /* Abstract type for PME that is defined only in the routine that use them. */
54 typedef struct gmx_pme *gmx_pme_t;
55
56
57
58 /* Structure describing the data in a single table */
59 typedef struct
60 {
61     enum gmx_table_interaction  interaction; /* Types of interactions stored in this table */
62     enum gmx_table_format       format;      /* Interpolation type and data format */
63
64     real                        r;           /* range of the table */
65     int                         n;           /* n+1 is the number of table points */
66     real                        scale;       /* distance (nm) between two table points */
67     real                        scale_exp;   /* distance for exponential part of VdW table, not always used */
68     real *                      data;        /* the actual table data */
69
70     /* Some information about the table layout. This can also be derived from the interpolation
71      * type and the table interactions, but it is convenient to have here for sanity checks, and it makes it
72      * much easier to access the tables in the nonbonded kernels when we can set the data from variables.
73      * It is always true that stride = formatsize*ninteractions
74      */
75     int                         formatsize;    /* Number of fp variables for each table point (1 for F, 2 for VF, 4 for YFGH, etc.) */
76     int                         ninteractions; /* Number of interactions in table, 1 for coul-only, 3 for coul+rep+disp. */
77     int                         stride;        /* Distance to next table point (number of fp variables per table point in total) */
78 } t_forcetable;
79
80 typedef struct
81 {
82     t_forcetable   table_elec;
83     t_forcetable   table_vdw;
84     t_forcetable   table_elec_vdw;
85
86     /* The actual neighbor lists, short and long range, see enum above
87      * for definition of neighborlist indices.
88      */
89     t_nblist nlist_sr[eNL_NR];
90     t_nblist nlist_lr[eNL_NR];
91 } t_nblists;
92
93 /* macros for the cginfo data in forcerec
94  *
95  * Since the tpx format support max 256 energy groups, we do the same here.
96  * Note that we thus have bits 8-14 still unused.
97  *
98  * The maximum cg size in cginfo is 63
99  * because we only have space for 6 bits in cginfo,
100  * this cg size entry is actually only read with domain decomposition.
101  * But there is a smaller limit due to the t_excl data structure
102  * which is defined in nblist.h.
103  */
104 #define SET_CGINFO_GID(cgi, gid)     (cgi) = (((cgi)  &  ~255) | (gid))
105 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   255)
106 #define SET_CGINFO_FEP(cgi)          (cgi) =  ((cgi)  |  (1<<15))
107 #define GET_CGINFO_FEP(cgi)        ( (cgi)            &  (1<<15))
108 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
109 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
110 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
111 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
112 #define SET_CGINFO_SOLOPT(cgi, opt)  (cgi) = (((cgi)  & ~(3<<18)) | ((opt)<<18))
113 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   3)
114 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
115 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
116 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
117 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
118 /* This bit is only used with bBondComm in the domain decomposition */
119 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
120 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
121 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
122 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
123 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
124 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
125 #define SET_CGINFO_NATOMS(cgi, opt)  (cgi) = (((cgi)  & ~(63<<25)) | ((opt)<<25))
126 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>25)       &   63)
127
128
129 /* Value to be used in mdrun for an infinite cut-off.
130  * Since we need to compare with the cut-off squared,
131  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
132  */
133 #define GMX_CUTOFF_INF 1E+18
134
135 /* enums for the neighborlist type */
136 enum {
137     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
138 };
139 /* OOR is "one over r" -- standard coul */
140 enum {
141     enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
142 };
143
144 enum {
145     egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
146     egCOUL14, egLJ14, egGB, egNR
147 };
148
149 typedef struct {
150     int   nener;      /* The number of energy group pairs     */
151     real *ener[egNR]; /* Energy terms for each pair of groups */
152 } gmx_grppairener_t;
153
154 typedef struct {
155     real              term[F_NRE];         /* The energies for all different interaction types */
156     gmx_grppairener_t grpp;
157     double            dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
158     double            dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
159     int               n_lambda;
160     int               fep_state;           /*current fep state -- just for printing */
161     double           *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
162     real              foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
163     gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
164 } gmx_enerdata_t;
165 /* The idea is that dvdl terms with linear lambda dependence will be added
166  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
167  * should explicitly determine the energies at foreign lambda points
168  * when n_lambda > 0.
169  */
170
171 typedef struct {
172     int  cg_start;
173     int  cg_end;
174     int  cg_mod;
175     int *cginfo;
176 } cginfo_mb_t;
177
178
179 /* ewald table type */
180 typedef struct ewald_tab *ewald_tab_t;
181
182 typedef struct {
183     rvec             *f;
184     int               f_nalloc;
185     unsigned          red_mask; /* Mask for marking which parts of f are filled */
186     rvec             *fshift;
187     real              ener[F_NRE];
188     gmx_grppairener_t grpp;
189     real              Vcorr_q;
190     real              Vcorr_lj;
191     real              dvdl[efptNR];
192     tensor            vir_q;
193     tensor            vir_lj;
194 } f_thread_t;
195
196 typedef struct {
197     interaction_const_t *ic;
198
199     /* Domain Decomposition */
200     gmx_bool bDomDec;
201
202     /* PBC stuff */
203     int                  ePBC;
204     gmx_bool             bMolPBC;
205     int                  rc_scaling;
206     rvec                 posres_com;
207     rvec                 posres_comB;
208
209     const gmx_hw_info_t *hwinfo;
210     const gmx_gpu_opt_t *gpu_opt;
211     gmx_bool             use_simd_kernels;
212
213     /* Interaction for calculated in kernels. In many cases this is similar to
214      * the electrostatics settings in the inputrecord, but the difference is that
215      * these variables always specify the actual interaction in the kernel - if
216      * we are tabulating reaction-field the inputrec will say reaction-field, but
217      * the kernel interaction will say cubic-spline-table. To be safe we also
218      * have a kernel-specific setting for the modifiers - if the interaction is
219      * tabulated we already included the inputrec modification there, so the kernel
220      * modification setting will say 'none' in that case.
221      */
222     int nbkernel_elec_interaction;
223     int nbkernel_vdw_interaction;
224     int nbkernel_elec_modifier;
225     int nbkernel_vdw_modifier;
226
227     /* Use special N*N kernels? */
228     gmx_bool bAllvsAll;
229     /* Private work data */
230     void    *AllvsAll_work;
231     void    *AllvsAll_workgb;
232
233     /* Cut-Off stuff.
234      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
235      */
236     real rlist, rlistlong;
237
238     /* Dielectric constant resp. multiplication factor for charges */
239     real zsquare, temp;
240     real epsilon_r, epsilon_rf, epsfac;
241
242     /* Constants for reaction fields */
243     real kappa, k_rf, c_rf;
244
245     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
246     double qsum[2];
247     double q2sum[2];
248     double c6sum[2];
249     rvec   mu_tot[2];
250
251     /* Dispersion correction stuff */
252     int  eDispCorr;
253
254     /* The shift of the shift or user potentials */
255     real enershiftsix;
256     real enershifttwelve;
257     /* Integrated differces for energy and virial with cut-off functions */
258     real enerdiffsix;
259     real enerdifftwelve;
260     real virdiffsix;
261     real virdifftwelve;
262     /* Constant for long range dispersion correction (average dispersion)
263      * for topology A/B ([0]/[1]) */
264     real avcsix[2];
265     /* Constant for long range repulsion term. Relative difference of about
266      * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
267      */
268     real avctwelve[2];
269
270     /* Fudge factors */
271     real fudgeQQ;
272
273     /* Table stuff */
274     gmx_bool     bcoultab;
275     gmx_bool     bvdwtab;
276     /* The normal tables are in the nblists struct(s) below */
277     t_forcetable tab14; /* for 1-4 interactions only */
278
279     /* PPPM & Shifting stuff */
280     int   coulomb_modifier;
281     real  rcoulomb_switch, rcoulomb;
282     real *phi;
283
284     /* VdW stuff */
285     int    vdw_modifier;
286     double reppow;
287     real   rvdw_switch, rvdw;
288     real   bham_b_max;
289
290     /* Free energy */
291     int      efep;
292     real     sc_alphavdw;
293     real     sc_alphacoul;
294     int      sc_power;
295     real     sc_r_power;
296     real     sc_sigma6_def;
297     real     sc_sigma6_min;
298     gmx_bool bSepDVDL;
299
300     /* NS Stuff */
301     int  eeltype;
302     int  vdwtype;
303     int  cg0, hcg;
304     /* solvent_opt contains the enum for the most common solvent
305      * in the system, which will be optimized.
306      * It can be set to esolNO to disable all water optimization */
307     int          solvent_opt;
308     int          nWatMol;
309     gmx_bool     bGrid;
310     gmx_bool     bExcl_IntraCGAll_InterCGNone;
311     cginfo_mb_t *cginfo_mb;
312     int         *cginfo;
313     rvec        *cg_cm;
314     int          cg_nalloc;
315     rvec        *shift_vec;
316
317     /* The neighborlists including tables */
318     int                 nnblists;
319     int                *gid2nblists;
320     t_nblists          *nblists;
321
322     int                 cutoff_scheme; /* group- or Verlet-style cutoff */
323     gmx_bool            bNonbonded;    /* true if nonbonded calculations are *not* turned off */
324     nonbonded_verlet_t *nbv;
325
326     /* The wall tables (if used) */
327     int            nwall;
328     t_forcetable **wall_tab;
329
330     /* The number of charge groups participating in do_force_lowlevel */
331     int ncg_force;
332     /* The number of atoms participating in do_force_lowlevel */
333     int natoms_force;
334     /* The number of atoms participating in force and constraints */
335     int natoms_force_constr;
336     /* The allocation size of vectors of size natoms_force */
337     int nalloc_force;
338
339     /* Twin Range stuff, f_twin has size natoms_force */
340     gmx_bool bTwinRange;
341     int      nlr;
342     rvec    *f_twin;
343     /* Constraint virial correction for multiple time stepping */
344     tensor   vir_twin_constr;
345
346     /* Forces that should not enter into the virial summation:
347      * PPPM/PME/Ewald/posres
348      */
349     gmx_bool bF_NoVirSum;
350     int      f_novirsum_n;
351     int      f_novirsum_nalloc;
352     rvec    *f_novirsum_alloc;
353     /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
354      * points to the normal force vectors wen pressure is not requested.
355      */
356     rvec *f_novirsum;
357
358     /* Long-range forces and virial for PPPM/PME/Ewald */
359     gmx_pme_t pmedata;
360     int       ljpme_combination_rule;
361     tensor    vir_el_recip;
362     tensor    vir_lj_recip;
363
364     /* PME/Ewald stuff */
365     gmx_bool    bEwald;
366     real        ewaldcoeff_q;
367     real        ewaldcoeff_lj;
368     ewald_tab_t ewald_table;
369
370     /* Virial Stuff */
371     rvec *fshift;
372     rvec  vir_diag_posres;
373     dvec  vir_wall_z;
374
375     /* Non bonded Parameter lists */
376     int      ntype; /* Number of atom types */
377     gmx_bool bBHAM;
378     real    *nbfp;
379     real    *ljpme_c6grid; /* C6-values used on grid in LJPME */
380
381     /* Energy group pair flags */
382     int *egp_flags;
383
384     /* Shell molecular dynamics flexible constraints */
385     real fc_stepsize;
386
387     /* Generalized born implicit solvent */
388     gmx_bool       bGB;
389     /* Generalized born stuff */
390     real           gb_epsilon_solvent;
391     /* Table data for GB */
392     t_forcetable   gbtab;
393     /* VdW radius for each atomtype (dim is thus ntype) */
394     real          *atype_radius;
395     /* Effective radius (derived from effective volume) for each type */
396     real          *atype_vol;
397     /* Implicit solvent - surface tension for each atomtype */
398     real          *atype_surftens;
399     /* Implicit solvent - radius for GB calculation */
400     real          *atype_gb_radius;
401     /* Implicit solvent - overlap for HCT model */
402     real          *atype_S_hct;
403     /* Generalized born interaction data */
404     gmx_genborn_t *born;
405
406     /* Table scale for GB */
407     real gbtabscale;
408     /* Table range for GB */
409     real gbtabr;
410     /* GB neighborlists (the sr list will contain for each atom all other atoms
411      * (for use in the SA calculation) and the lr list will contain
412      * for each atom all atoms 1-4 or greater (for use in the GB calculation)
413      */
414     t_nblist gblist_sr;
415     t_nblist gblist_lr;
416     t_nblist gblist;
417
418     /* Inverse square root of the Born radii for implicit solvent */
419     real *invsqrta;
420     /* Derivatives of the potential with respect to the Born radii */
421     real *dvda;
422     /* Derivatives of the Born radii with respect to coordinates */
423     real *dadx;
424     real *dadx_rawptr;
425     int   nalloc_dadx; /* Allocated size of dadx */
426
427     /* If > 0 signals Test Particle Insertion,
428      * the value is the number of atoms of the molecule to insert
429      * Only the energy difference due to the addition of the last molecule
430      * should be calculated.
431      */
432     gmx_bool n_tpi;
433
434     /* Neighbor searching stuff */
435     gmx_ns_t ns;
436
437     /* QMMM stuff */
438     gmx_bool         bQMMM;
439     t_QMMMrec       *qr;
440
441     /* QM-MM neighborlists */
442     t_nblist QMMMlist;
443
444     /* Limit for printing large forces, negative is don't print */
445     real print_force;
446
447     /* coarse load balancing time measurement */
448     double t_fnbf;
449     double t_wait;
450     int    timesteps;
451
452     /* parameter needed for AdResS simulation */
453     int             adress_type;
454     gmx_bool        badress_tf_full_box;
455     real            adress_const_wf;
456     real            adress_ex_width;
457     real            adress_hy_width;
458     int             adress_icor;
459     int             adress_site;
460     rvec            adress_refs;
461     int             n_adress_tf_grps;
462     int           * adress_tf_table_index;
463     int            *adress_group_explicit;
464     t_forcetable *  atf_tabs;
465     real            adress_ex_forcecap;
466     gmx_bool        adress_do_hybridpairs;
467
468     /* User determined parameters, copied from the inputrec */
469     int  userint1;
470     int  userint2;
471     int  userint3;
472     int  userint4;
473     real userreal1;
474     real userreal2;
475     real userreal3;
476     real userreal4;
477
478     /* Thread local force and energy data */
479     /* FIXME move to bonded_thread_data_t */
480     int         nthreads;
481     int         red_ashift;
482     int         red_nblock;
483     f_thread_t *f_t;
484
485     /* Exclusion load distribution over the threads */
486     int  *excl_load;
487 } t_forcerec;
488
489 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
490  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
491  * in the code, but beware if you are using these macros externally.
492  */
493 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
494 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
495 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
496 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
497 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
498
499 #ifdef __cplusplus
500 }
501 #endif