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