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