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