Remove support for implicit solvation
[alexxy/gromacs.git] / src / gromacs / mdtypes / 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,2016,2017,2018, 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 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
38 #define GMX_MDTYPES_TYPES_FORCEREC_H
39
40 #include "gromacs/math/vectypes.h"
41 #ifdef __cplusplus
42 #include "gromacs/math/paddedvector.h"
43 #endif
44 #include "gromacs/mdtypes/interaction_const.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/topology/idef.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
49
50 struct ForceProviders;
51
52 /* Abstract type for PME that is defined only in the routine that use them. */
53 struct gmx_ns_t;
54 struct gmx_pme_t;
55 struct nonbonded_verlet_t;
56 struct bonded_threading_t;
57 struct t_forcetable;
58 struct t_nblist;
59 struct t_nblists;
60 struct t_QMMMrec;
61
62 #ifdef __cplusplus
63 extern "C" {
64 #endif
65
66 /* macros for the cginfo data in forcerec
67  *
68  * Since the tpx format support max 256 energy groups, we do the same here.
69  * Note that we thus have bits 8-14 still unused.
70  *
71  * The maximum cg size in cginfo is 63
72  * because we only have space for 6 bits in cginfo,
73  * this cg size entry is actually only read with domain decomposition.
74  * But there is a smaller limit due to the t_excl data structure
75  * which is defined in nblist.h.
76  */
77 #define SET_CGINFO_GID(cgi, gid)     (cgi) = (((cgi)  &  ~255) | (gid))
78 #define GET_CGINFO_GID(cgi)        ( (cgi)            &   255)
79 #define SET_CGINFO_FEP(cgi)          (cgi) =  ((cgi)  |  (1<<15))
80 #define GET_CGINFO_FEP(cgi)        ( (cgi)            &  (1<<15))
81 #define SET_CGINFO_EXCL_INTRA(cgi)   (cgi) =  ((cgi)  |  (1<<16))
82 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi)            &  (1<<16))
83 #define SET_CGINFO_EXCL_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<17))
84 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi)            &  (1<<17))
85 #define SET_CGINFO_SOLOPT(cgi, opt)  (cgi) = (((cgi)  & ~(3<<18)) | ((opt)<<18))
86 #define GET_CGINFO_SOLOPT(cgi)     (((cgi)>>18)       &   3)
87 #define SET_CGINFO_CONSTR(cgi)       (cgi) =  ((cgi)  |  (1<<20))
88 #define GET_CGINFO_CONSTR(cgi)     ( (cgi)            &  (1<<20))
89 #define SET_CGINFO_SETTLE(cgi)       (cgi) =  ((cgi)  |  (1<<21))
90 #define GET_CGINFO_SETTLE(cgi)     ( (cgi)            &  (1<<21))
91 /* This bit is only used with bBondComm in the domain decomposition */
92 #define SET_CGINFO_BOND_INTER(cgi)   (cgi) =  ((cgi)  |  (1<<22))
93 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi)            &  (1<<22))
94 #define SET_CGINFO_HAS_VDW(cgi)      (cgi) =  ((cgi)  |  (1<<23))
95 #define GET_CGINFO_HAS_VDW(cgi)    ( (cgi)            &  (1<<23))
96 #define SET_CGINFO_HAS_Q(cgi)        (cgi) =  ((cgi)  |  (1<<24))
97 #define GET_CGINFO_HAS_Q(cgi)      ( (cgi)            &  (1<<24))
98 #define SET_CGINFO_NATOMS(cgi, opt)  (cgi) = (((cgi)  & ~(63<<25)) | ((opt)<<25))
99 #define GET_CGINFO_NATOMS(cgi)     (((cgi)>>25)       &   63)
100
101
102 /* Value to be used in mdrun for an infinite cut-off.
103  * Since we need to compare with the cut-off squared,
104  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
105  */
106 #define GMX_CUTOFF_INF 1E+18
107
108 /* enums for the neighborlist type */
109 enum {
110     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
111 };
112
113 enum {
114     egCOULSR, egLJSR, egBHAMSR,
115     egCOUL14, egLJ14, egNR
116 };
117 extern const char *egrp_nm[egNR+1];
118
119 struct gmx_grppairener_t
120 {
121     int   nener;      /* The number of energy group pairs     */
122     real *ener[egNR]; /* Energy terms for each pair of groups */
123 };
124
125 struct gmx_enerdata_t
126 {
127     real                     term[F_NRE];         /* The energies for all different interaction types */
128     struct gmx_grppairener_t grpp;
129     double                   dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
130     double                   dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
131     int                      n_lambda;
132     int                      fep_state;           /*current fep state -- just for printing */
133     double                  *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
134     real                     foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
135     struct gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
136 };
137 /* The idea is that dvdl terms with linear lambda dependence will be added
138  * automatically to enerpart_lambda. Terms with non-linear lambda dependence
139  * should explicitly determine the energies at foreign lambda points
140  * when n_lambda > 0.
141  */
142
143 struct cginfo_mb_t
144 {
145     int  cg_start;
146     int  cg_end;
147     int  cg_mod;
148     int *cginfo;
149 };
150
151
152 /* Forward declaration of type for managing Ewald tables */
153 struct gmx_ewald_tab_t;
154
155 struct ewald_corr_thread_t;
156
157 struct t_forcerec {
158     struct interaction_const_t *ic;
159
160     /* Domain Decomposition */
161     gmx_bool bDomDec;
162
163     /* PBC stuff */
164     int                         ePBC;
165     gmx_bool                    bMolPBC;
166     int                         rc_scaling;
167     rvec                        posres_com;
168     rvec                        posres_comB;
169
170     gmx_bool                    use_simd_kernels;
171
172     /* Interaction for calculated in kernels. In many cases this is similar to
173      * the electrostatics settings in the inputrecord, but the difference is that
174      * these variables always specify the actual interaction in the kernel - if
175      * we are tabulating reaction-field the inputrec will say reaction-field, but
176      * the kernel interaction will say cubic-spline-table. To be safe we also
177      * have a kernel-specific setting for the modifiers - if the interaction is
178      * tabulated we already included the inputrec modification there, so the kernel
179      * modification setting will say 'none' in that case.
180      */
181     int nbkernel_elec_interaction;
182     int nbkernel_vdw_interaction;
183     int nbkernel_elec_modifier;
184     int nbkernel_vdw_modifier;
185
186     /* Use special N*N kernels? */
187     gmx_bool bAllvsAll;
188     /* Private work data */
189     void    *AllvsAll_work;
190
191     /* Cut-Off stuff.
192      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
193      */
194     real rlist;
195
196     /* Parameters for generalized reaction field */
197     real zsquare, temp;
198
199     /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
200     double qsum[2];
201     double q2sum[2];
202     double c6sum[2];
203     rvec   mu_tot[2];
204
205     /* Dispersion correction stuff */
206     int                  eDispCorr;
207     int                  numAtomsForDispersionCorrection;
208     struct t_forcetable *dispersionCorrectionTable;
209
210     /* The shift of the shift or user potentials */
211     real enershiftsix;
212     real enershifttwelve;
213     /* Integrated differces for energy and virial with cut-off functions */
214     real enerdiffsix;
215     real enerdifftwelve;
216     real virdiffsix;
217     real virdifftwelve;
218     /* Constant for long range dispersion correction (average dispersion)
219      * for topology A/B ([0]/[1]) */
220     real avcsix[2];
221     /* Constant for long range repulsion term. Relative difference of about
222      * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
223      */
224     real avctwelve[2];
225
226     /* Fudge factors */
227     real fudgeQQ;
228
229     /* Table stuff */
230     gmx_bool             bcoultab;
231     gmx_bool             bvdwtab;
232     /* The normal tables are in the nblists struct(s) below */
233
234     struct t_forcetable *pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
235
236     /* Free energy */
237     int      efep;
238     real     sc_alphavdw;
239     real     sc_alphacoul;
240     int      sc_power;
241     real     sc_r_power;
242     real     sc_sigma6_def;
243     real     sc_sigma6_min;
244
245     /* NS Stuff */
246     int  cg0, hcg;
247     /* solvent_opt contains the enum for the most common solvent
248      * in the system, which will be optimized.
249      * It can be set to esolNO to disable all water optimization */
250     int                 solvent_opt;
251     int                 nWatMol;
252     gmx_bool            bGrid;
253     gmx_bool            bExcl_IntraCGAll_InterCGNone;
254     struct cginfo_mb_t *cginfo_mb;
255     int                *cginfo;
256     rvec               *cg_cm;
257     int                 cg_nalloc;
258     rvec               *shift_vec;
259
260     /* The neighborlists including tables */
261     int                        nnblists;
262     int                       *gid2nblists;
263     struct t_nblists          *nblists;
264
265     int                        cutoff_scheme; /* group- or Verlet-style cutoff */
266     gmx_bool                   bNonbonded;    /* true if nonbonded calculations are *not* turned off */
267     struct nonbonded_verlet_t *nbv;
268
269     /* The wall tables (if used) */
270     int                    nwall;
271     struct t_forcetable ***wall_tab;
272
273     /* The number of charge groups participating in do_force_lowlevel */
274     int ncg_force;
275     /* The number of atoms participating in do_force_lowlevel */
276     int natoms_force;
277     /* The number of atoms participating in force and constraints */
278     int natoms_force_constr;
279     /* The allocation size of vectors of size natoms_force */
280     int nalloc_force;
281
282     /* Forces that should not enter into the coord x force virial summation:
283      * PPPM/PME/Ewald/posres/ForceProviders
284      */
285     /* True when we have contributions that are directly added to the virial */
286     gmx_bool          haveDirectVirialContributions;
287 #ifdef __cplusplus
288     /* TODO: Replace the pointer by an object once we got rid of C */
289     std::vector<gmx::RVec>  *forceBufferForDirectVirialContributions;
290 #else
291     void                    *forceBufferForDirectVirialContributions_dummy;
292 #endif
293
294     /* Data for PPPM/PME/Ewald */
295     struct gmx_pme_t *pmedata;
296     int               ljpme_combination_rule;
297
298     /* PME/Ewald stuff */
299     struct gmx_ewald_tab_t *ewald_table;
300
301     /* Virial Stuff */
302     rvec *fshift;
303     dvec  vir_wall_z;
304
305     /* Non bonded Parameter lists */
306     int      ntype; /* Number of atom types */
307     gmx_bool bBHAM;
308     real    *nbfp;
309     real    *ljpme_c6grid; /* C6-values used on grid in LJPME */
310
311     /* Energy group pair flags */
312     int *egp_flags;
313
314     /* Shell molecular dynamics flexible constraints */
315     real fc_stepsize;
316
317     /* If > 0 signals Test Particle Insertion,
318      * the value is the number of atoms of the molecule to insert
319      * Only the energy difference due to the addition of the last molecule
320      * should be calculated.
321      */
322     gmx_bool n_tpi;
323
324     /* Neighbor searching stuff */
325     struct gmx_ns_t *ns;
326
327     /* QMMM stuff */
328     gmx_bool          bQMMM;
329     struct t_QMMMrec *qr;
330
331     /* QM-MM neighborlists */
332     struct t_nblist        *QMMMlist;
333
334     /* Limit for printing large forces, negative is don't print */
335     real print_force;
336
337     /* coarse load balancing time measurement */
338     double t_fnbf;
339     double t_wait;
340     int    timesteps;
341
342     /* User determined parameters, copied from the inputrec */
343     int  userint1;
344     int  userint2;
345     int  userint3;
346     int  userint4;
347     real userreal1;
348     real userreal2;
349     real userreal3;
350     real userreal4;
351
352     /* Pointer to struct for managing threading of bonded force calculation */
353     struct bonded_threading_t *bonded_threading;
354
355     /* Ewald correction thread local virial and energy data */
356     int                         nthread_ewc;
357     struct ewald_corr_thread_t *ewc_t;
358
359     struct ForceProviders      *forceProviders;
360 };
361
362 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
363  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
364  * in the code, but beware if you are using these macros externally.
365  */
366 #define C6(nbfp, ntp, ai, aj)     (nbfp)[2*((ntp)*(ai)+(aj))]
367 #define C12(nbfp, ntp, ai, aj)    (nbfp)[2*((ntp)*(ai)+(aj))+1]
368 #define BHAMC(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))]
369 #define BHAMA(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+1]
370 #define BHAMB(nbfp, ntp, ai, aj)  (nbfp)[3*((ntp)*(ai)+(aj))+2]
371
372 #ifdef __cplusplus
373 }
374 #endif
375
376 #endif