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