Calculating LJ-interactions using Particle mesh Ewald
[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, 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 "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 /* The maximum cg size in cginfo is 63
95  * because we only have space for 6 bits in cginfo,
96  * this cg size entry is actually only read with domain decomposition.
97  * But there is a smaller limit due to the t_excl data structure
98  * which is defined in nblist.h.
99  */
100 #define SET_CGINFO_GID(cgi, gid)      (cgi) = (((cgi)  &  ~65535)  |  (gid)   )
101 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   65535)
102 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
103 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
104 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
105 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
106 #define SET_CGINFO_SOLOPT(cgi, opt)   (cgi) = (((cgi)  & ~(3<<18)) | ((opt)<<18))
107 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   3)
108 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
109 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
110 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
111 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
112 /* This bit is only used with bBondComm in the domain decomposition */
113 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
114 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
115 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
116 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
117 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
118 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
119 #define SET_CGINFO_NATOMS(cgi, opt)   (cgi) = (((cgi)  & ~(63<<25)) | ((opt)<<25))
120 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>25)       &   63)
121
122
123 /* Value to be used in mdrun for an infinite cut-off.
124  * Since we need to compare with the cut-off squared,
125  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
126  */
127 #define GMX_CUTOFF_INF 1E+18
128
129 /* enums for the neighborlist type */
130 enum {
131     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
132 };
133 /* OOR is "one over r" -- standard coul */
134 enum {
135     enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
136 };
137
138 enum {
139     egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
140     egCOUL14, egLJ14, egGB, egNR
141 };
142
143 typedef struct {
144     int   nener;      /* The number of energy group pairs     */
145     real *ener[egNR]; /* Energy terms for each pair of groups */
146 } gmx_grppairener_t;
147
148 typedef struct {
149     real              term[F_NRE];         /* The energies for all different interaction types */
150     gmx_grppairener_t grpp;
151     double            dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
152     double            dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
153     int               n_lambda;
154     int               fep_state;           /*current fep state -- just for printing */
155     double           *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
156     real              foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
157     gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
158 } gmx_enerdata_t;
159 /* The idea is that dvdl terms with linear lambda dependence will be added
160  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
161  * should explicitly determine the energies at foreign lambda points
162  * when n_lambda > 0.
163  */
164
165 typedef struct {
166     int  cg_start;
167     int  cg_end;
168     int  cg_mod;
169     int *cginfo;
170 } cginfo_mb_t;
171
172
173 /* ewald table type */
174 typedef struct ewald_tab *ewald_tab_t;
175
176 typedef struct {
177     rvec             *f;
178     int               f_nalloc;
179     unsigned          red_mask; /* Mask for marking which parts of f are filled */
180     rvec             *fshift;
181     real              ener[F_NRE];
182     gmx_grppairener_t grpp;
183     real              Vcorr_q;
184     real              Vcorr_lj;
185     real              dvdl[efptNR];
186     tensor            vir_q;
187     tensor            vir_lj;
188 } 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_cpu_acceleration;
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     gmx_bool bSepDVDL;
293
294     /* NS Stuff */
295     int  eeltype;
296     int  vdwtype;
297     int  cg0, hcg;
298     /* solvent_opt contains the enum for the most common solvent
299      * in the system, which will be optimized.
300      * It can be set to esolNO to disable all water optimization */
301     int          solvent_opt;
302     int          nWatMol;
303     gmx_bool     bGrid;
304     gmx_bool     bExcl_IntraCGAll_InterCGNone;
305     cginfo_mb_t *cginfo_mb;
306     int         *cginfo;
307     rvec        *cg_cm;
308     int          cg_nalloc;
309     rvec        *shift_vec;
310
311     /* The neighborlists including tables */
312     int                 nnblists;
313     int                *gid2nblists;
314     t_nblists          *nblists;
315
316     int                 cutoff_scheme; /* group- or Verlet-style cutoff */
317     gmx_bool            bNonbonded;    /* true if nonbonded calculations are *not* turned off */
318     nonbonded_verlet_t *nbv;
319
320     /* The wall tables (if used) */
321     int            nwall;
322     t_forcetable **wall_tab;
323
324     /* The number of charge groups participating in do_force_lowlevel */
325     int ncg_force;
326     /* The number of atoms participating in do_force_lowlevel */
327     int natoms_force;
328     /* The number of atoms participating in force and constraints */
329     int natoms_force_constr;
330     /* The allocation size of vectors of size natoms_force */
331     int nalloc_force;
332
333     /* Twin Range stuff, f_twin has size natoms_force */
334     gmx_bool bTwinRange;
335     int      nlr;
336     rvec    *f_twin;
337
338     /* Forces that should not enter into the virial summation:
339      * PPPM/PME/Ewald/posres
340      */
341     gmx_bool bF_NoVirSum;
342     int      f_novirsum_n;
343     int      f_novirsum_nalloc;
344     rvec    *f_novirsum_alloc;
345     /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
346      * points to the normal force vectors wen pressure is not requested.
347      */
348     rvec *f_novirsum;
349
350     /* Long-range forces and virial for PPPM/PME/Ewald */
351     gmx_pme_t pmedata;
352     int       ljpme_combination_rule;
353     tensor    vir_el_recip;
354     tensor    vir_lj_recip;
355
356     /* PME/Ewald stuff */
357     gmx_bool    bEwald;
358     real        ewaldcoeff_q;
359     real        ewaldcoeff_lj;
360     ewald_tab_t ewald_table;
361
362     /* Virial Stuff */
363     rvec *fshift;
364     rvec  vir_diag_posres;
365     dvec  vir_wall_z;
366
367     /* Non bonded Parameter lists */
368     int      ntype; /* Number of atom types */
369     gmx_bool bBHAM;
370     real    *nbfp;
371
372     /* Energy group pair flags */
373     int *egp_flags;
374
375     /* Shell molecular dynamics flexible constraints */
376     real fc_stepsize;
377
378     /* Generalized born implicit solvent */
379     gmx_bool       bGB;
380     /* Generalized born stuff */
381     real           gb_epsilon_solvent;
382     /* Table data for GB */
383     t_forcetable   gbtab;
384     /* VdW radius for each atomtype (dim is thus ntype) */
385     real          *atype_radius;
386     /* Effective radius (derived from effective volume) for each type */
387     real          *atype_vol;
388     /* Implicit solvent - surface tension for each atomtype */
389     real          *atype_surftens;
390     /* Implicit solvent - radius for GB calculation */
391     real          *atype_gb_radius;
392     /* Implicit solvent - overlap for HCT model */
393     real          *atype_S_hct;
394     /* Generalized born interaction data */
395     gmx_genborn_t *born;
396
397     /* Table scale for GB */
398     real gbtabscale;
399     /* Table range for GB */
400     real gbtabr;
401     /* GB neighborlists (the sr list will contain for each atom all other atoms
402      * (for use in the SA calculation) and the lr list will contain
403      * for each atom all atoms 1-4 or greater (for use in the GB calculation)
404      */
405     t_nblist gblist_sr;
406     t_nblist gblist_lr;
407     t_nblist gblist;
408
409     /* Inverse square root of the Born radii for implicit solvent */
410     real *invsqrta;
411     /* Derivatives of the potential with respect to the Born radii */
412     real *dvda;
413     /* Derivatives of the Born radii with respect to coordinates */
414     real *dadx;
415     real *dadx_rawptr;
416     int   nalloc_dadx; /* Allocated size of dadx */
417
418     /* If > 0 signals Test Particle Insertion,
419      * the value is the number of atoms of the molecule to insert
420      * Only the energy difference due to the addition of the last molecule
421      * should be calculated.
422      */
423     gmx_bool n_tpi;
424
425     /* Neighbor searching stuff */
426     gmx_ns_t ns;
427
428     /* QMMM stuff */
429     gmx_bool         bQMMM;
430     t_QMMMrec       *qr;
431
432     /* QM-MM neighborlists */
433     t_nblist QMMMlist;
434
435     /* Limit for printing large forces, negative is don't print */
436     real print_force;
437
438     /* coarse load balancing time measurement */
439     double t_fnbf;
440     double t_wait;
441     int    timesteps;
442
443     /* parameter needed for AdResS simulation */
444     int             adress_type;
445     gmx_bool        badress_tf_full_box;
446     real            adress_const_wf;
447     real            adress_ex_width;
448     real            adress_hy_width;
449     int             adress_icor;
450     int             adress_site;
451     rvec            adress_refs;
452     int             n_adress_tf_grps;
453     int           * adress_tf_table_index;
454     int            *adress_group_explicit;
455     t_forcetable *  atf_tabs;
456     real            adress_ex_forcecap;
457     gmx_bool        adress_do_hybridpairs;
458
459     /* User determined parameters, copied from the inputrec */
460     int  userint1;
461     int  userint2;
462     int  userint3;
463     int  userint4;
464     real userreal1;
465     real userreal2;
466     real userreal3;
467     real userreal4;
468
469     /* Thread local force and energy data */
470     /* FIXME move to bonded_thread_data_t */
471     int         nthreads;
472     int         red_ashift;
473     int         red_nblock;
474     f_thread_t *f_t;
475
476     /* Exclusion load distribution over the threads */
477     int  *excl_load;
478 } t_forcerec;
479
480 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
481  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
482  * in the code, but beware if you are using these macros externally.
483  */
484 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
485 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
486 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
487 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
488 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
489
490 #ifdef __cplusplus
491 }
492 #endif