59cd9f4279d61196597b5d97cccad70202711039
[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/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 typedef struct gmx_pme *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 {
189     rvec             *f;
190     int               f_nalloc;
191     unsigned          red_mask; /* Mask for marking which parts of f are filled */
192     rvec             *fshift;
193     real              ener[F_NRE];
194     gmx_grppairener_t grpp;
195     real              Vcorr_q;
196     real              Vcorr_lj;
197     real              dvdl[efptNR];
198     tensor            vir_q;
199     tensor            vir_lj;
200 } f_thread_t;
201
202 typedef struct {
203     interaction_const_t *ic;
204
205     /* Domain Decomposition */
206     gmx_bool bDomDec;
207
208     /* PBC stuff */
209     int                  ePBC;
210     gmx_bool             bMolPBC;
211     int                  rc_scaling;
212     rvec                 posres_com;
213     rvec                 posres_comB;
214
215     const gmx_hw_info_t *hwinfo;
216     const gmx_gpu_opt_t *gpu_opt;
217     gmx_bool             use_simd_kernels;
218
219     /* Interaction for calculated in kernels. In many cases this is similar to
220      * the electrostatics settings in the inputrecord, but the difference is that
221      * these variables always specify the actual interaction in the kernel - if
222      * we are tabulating reaction-field the inputrec will say reaction-field, but
223      * the kernel interaction will say cubic-spline-table. To be safe we also
224      * have a kernel-specific setting for the modifiers - if the interaction is
225      * tabulated we already included the inputrec modification there, so the kernel
226      * modification setting will say 'none' in that case.
227      */
228     int nbkernel_elec_interaction;
229     int nbkernel_vdw_interaction;
230     int nbkernel_elec_modifier;
231     int nbkernel_vdw_modifier;
232
233     /* Use special N*N kernels? */
234     gmx_bool bAllvsAll;
235     /* Private work data */
236     void    *AllvsAll_work;
237     void    *AllvsAll_workgb;
238
239     /* Cut-Off stuff.
240      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
241      */
242     real rlist, rlistlong;
243
244     /* Dielectric constant resp. multiplication factor for charges */
245     real zsquare, temp;
246     real epsilon_r, epsilon_rf, epsfac;
247
248     /* Constants for reaction fields */
249     real kappa, k_rf, c_rf;
250
251     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
252     double qsum[2];
253     double q2sum[2];
254     double c6sum[2];
255     rvec   mu_tot[2];
256
257     /* Dispersion correction stuff */
258     int  eDispCorr;
259
260     /* The shift of the shift or user potentials */
261     real enershiftsix;
262     real enershifttwelve;
263     /* Integrated differces for energy and virial with cut-off functions */
264     real enerdiffsix;
265     real enerdifftwelve;
266     real virdiffsix;
267     real virdifftwelve;
268     /* Constant for long range dispersion correction (average dispersion)
269      * for topology A/B ([0]/[1]) */
270     real avcsix[2];
271     /* Constant for long range repulsion term. Relative difference of about
272      * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
273      */
274     real avctwelve[2];
275
276     /* Fudge factors */
277     real fudgeQQ;
278
279     /* Table stuff */
280     gmx_bool     bcoultab;
281     gmx_bool     bvdwtab;
282     /* The normal tables are in the nblists struct(s) below */
283     t_forcetable tab14; /* for 1-4 interactions only */
284
285     /* PPPM & Shifting stuff */
286     int   coulomb_modifier;
287     real  rcoulomb_switch, rcoulomb;
288     real *phi;
289
290     /* VdW stuff */
291     int    vdw_modifier;
292     double reppow;
293     real   rvdw_switch, rvdw;
294     real   bham_b_max;
295
296     /* Free energy */
297     int      efep;
298     real     sc_alphavdw;
299     real     sc_alphacoul;
300     int      sc_power;
301     real     sc_r_power;
302     real     sc_sigma6_def;
303     real     sc_sigma6_min;
304
305     /* NS Stuff */
306     int  eeltype;
307     int  vdwtype;
308     int  cg0, hcg;
309     /* solvent_opt contains the enum for the most common solvent
310      * in the system, which will be optimized.
311      * It can be set to esolNO to disable all water optimization */
312     int          solvent_opt;
313     int          nWatMol;
314     gmx_bool     bGrid;
315     gmx_bool     bExcl_IntraCGAll_InterCGNone;
316     cginfo_mb_t *cginfo_mb;
317     int         *cginfo;
318     rvec        *cg_cm;
319     int          cg_nalloc;
320     rvec        *shift_vec;
321
322     /* The neighborlists including tables */
323     int                        nnblists;
324     int                       *gid2nblists;
325     t_nblists                 *nblists;
326
327     int                        cutoff_scheme; /* group- or Verlet-style cutoff */
328     gmx_bool                   bNonbonded;    /* true if nonbonded calculations are *not* turned off */
329     struct nonbonded_verlet_t *nbv;
330
331     /* The wall tables (if used) */
332     int            nwall;
333     t_forcetable **wall_tab;
334
335     /* The number of charge groups participating in do_force_lowlevel */
336     int ncg_force;
337     /* The number of atoms participating in do_force_lowlevel */
338     int natoms_force;
339     /* The number of atoms participating in force and constraints */
340     int natoms_force_constr;
341     /* The allocation size of vectors of size natoms_force */
342     int nalloc_force;
343
344     /* Twin Range stuff, f_twin has size natoms_force */
345     gmx_bool bTwinRange;
346     int      nlr;
347     rvec    *f_twin;
348     /* Constraint virial correction for multiple time stepping */
349     tensor   vir_twin_constr;
350
351     /* Forces that should not enter into the virial summation:
352      * PPPM/PME/Ewald/posres
353      */
354     gmx_bool bF_NoVirSum;
355     int      f_novirsum_n;
356     int      f_novirsum_nalloc;
357     rvec    *f_novirsum_alloc;
358     /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
359      * points to the normal force vectors wen pressure is not requested.
360      */
361     rvec *f_novirsum;
362
363     /* Long-range forces and virial for PPPM/PME/Ewald */
364     gmx_pme_t pmedata;
365     int       ljpme_combination_rule;
366     tensor    vir_el_recip;
367     tensor    vir_lj_recip;
368
369     /* PME/Ewald stuff */
370     gmx_bool                bEwald;
371     real                    ewaldcoeff_q;
372     real                    ewaldcoeff_lj;
373     struct gmx_ewald_tab_t *ewald_table;
374
375     /* Virial Stuff */
376     rvec *fshift;
377     rvec  vir_diag_posres;
378     dvec  vir_wall_z;
379
380     /* Non bonded Parameter lists */
381     int      ntype; /* Number of atom types */
382     gmx_bool bBHAM;
383     real    *nbfp;
384     real    *ljpme_c6grid; /* C6-values used on grid in LJPME */
385
386     /* Energy group pair flags */
387     int *egp_flags;
388
389     /* Shell molecular dynamics flexible constraints */
390     real fc_stepsize;
391
392     /* Generalized born implicit solvent */
393     gmx_bool       bGB;
394     /* Generalized born stuff */
395     real           gb_epsilon_solvent;
396     /* Table data for GB */
397     t_forcetable   gbtab;
398     /* VdW radius for each atomtype (dim is thus ntype) */
399     real          *atype_radius;
400     /* Effective radius (derived from effective volume) for each type */
401     real          *atype_vol;
402     /* Implicit solvent - surface tension for each atomtype */
403     real          *atype_surftens;
404     /* Implicit solvent - radius for GB calculation */
405     real          *atype_gb_radius;
406     /* Implicit solvent - overlap for HCT model */
407     real          *atype_S_hct;
408     /* Generalized born interaction data */
409     gmx_genborn_t *born;
410
411     /* Table scale for GB */
412     real gbtabscale;
413     /* Table range for GB */
414     real gbtabr;
415     /* GB neighborlists (the sr list will contain for each atom all other atoms
416      * (for use in the SA calculation) and the lr list will contain
417      * for each atom all atoms 1-4 or greater (for use in the GB calculation)
418      */
419     t_nblist gblist_sr;
420     t_nblist gblist_lr;
421     t_nblist gblist;
422
423     /* Inverse square root of the Born radii for implicit solvent */
424     real *invsqrta;
425     /* Derivatives of the potential with respect to the Born radii */
426     real *dvda;
427     /* Derivatives of the Born radii with respect to coordinates */
428     real *dadx;
429     real *dadx_rawptr;
430     int   nalloc_dadx; /* Allocated size of dadx */
431
432     /* If > 0 signals Test Particle Insertion,
433      * the value is the number of atoms of the molecule to insert
434      * Only the energy difference due to the addition of the last molecule
435      * should be calculated.
436      */
437     gmx_bool n_tpi;
438
439     /* Neighbor searching stuff */
440     gmx_ns_t ns;
441
442     /* QMMM stuff */
443     gmx_bool         bQMMM;
444     t_QMMMrec       *qr;
445
446     /* QM-MM neighborlists */
447     t_nblist QMMMlist;
448
449     /* Limit for printing large forces, negative is don't print */
450     real print_force;
451
452     /* coarse load balancing time measurement */
453     double t_fnbf;
454     double t_wait;
455     int    timesteps;
456
457     /* parameter needed for AdResS simulation */
458     int             adress_type;
459     gmx_bool        badress_tf_full_box;
460     real            adress_const_wf;
461     real            adress_ex_width;
462     real            adress_hy_width;
463     int             adress_icor;
464     int             adress_site;
465     rvec            adress_refs;
466     int             n_adress_tf_grps;
467     int           * adress_tf_table_index;
468     int            *adress_group_explicit;
469     t_forcetable *  atf_tabs;
470     real            adress_ex_forcecap;
471     gmx_bool        adress_do_hybridpairs;
472
473     /* User determined parameters, copied from the inputrec */
474     int  userint1;
475     int  userint2;
476     int  userint3;
477     int  userint4;
478     real userreal1;
479     real userreal2;
480     real userreal3;
481     real userreal4;
482
483     /* Thread local force and energy data */
484     /* FIXME move to bonded_thread_data_t */
485     int         nthreads;
486     int         red_ashift;
487     int         red_nblock;
488     f_thread_t *f_t;
489
490     /* Exclusion load distribution over the threads */
491     int  *excl_load;
492 } t_forcerec;
493
494 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
495  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
496  * in the code, but beware if you are using these macros externally.
497  */
498 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
499 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
500 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
501 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
502 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
503
504 #ifdef __cplusplus
505 }
506 #endif
507 #endif