Decouple task assignment from task execution
[gromacs.git] / src / gromacs / mdtypes / forcerec.h
blob62f5030711ef961d34d5014ed27ea4c696dc2ee3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
37 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
38 #define GMX_MDTYPES_TYPES_FORCEREC_H
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"
50 struct ForceProviders;
52 /* Abstract type for PME that is defined only in the routine that use them. */
53 struct gmx_genborn_t;
54 struct gmx_ns_t;
55 struct gmx_pme_t;
56 struct nonbonded_verlet_t;
57 struct bonded_threading_t;
58 struct t_forcetable;
59 struct t_nblist;
60 struct t_nblists;
61 struct t_QMMMrec;
63 #ifdef __cplusplus
64 extern "C" {
65 #endif
67 /* macros for the cginfo data in forcerec
69 * Since the tpx format support max 256 energy groups, we do the same here.
70 * Note that we thus have bits 8-14 still unused.
72 * The maximum cg size in cginfo is 63
73 * because we only have space for 6 bits in cginfo,
74 * this cg size entry is actually only read with domain decomposition.
75 * But there is a smaller limit due to the t_excl data structure
76 * which is defined in nblist.h.
78 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
79 #define GET_CGINFO_GID(cgi) ( (cgi) & 255)
80 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1<<15))
81 #define GET_CGINFO_FEP(cgi) ( (cgi) & (1<<15))
82 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
83 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
84 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
85 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
86 #define SET_CGINFO_SOLOPT(cgi, opt) (cgi) = (((cgi) & ~(3<<18)) | ((opt)<<18))
87 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 3)
88 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1<<20))
89 #define GET_CGINFO_CONSTR(cgi) ( (cgi) & (1<<20))
90 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1<<21))
91 #define GET_CGINFO_SETTLE(cgi) ( (cgi) & (1<<21))
92 /* This bit is only used with bBondComm in the domain decomposition */
93 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1<<22))
94 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi) & (1<<22))
95 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1<<23))
96 #define GET_CGINFO_HAS_VDW(cgi) ( (cgi) & (1<<23))
97 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1<<24))
98 #define GET_CGINFO_HAS_Q(cgi) ( (cgi) & (1<<24))
99 #define SET_CGINFO_NATOMS(cgi, opt) (cgi) = (((cgi) & ~(63<<25)) | ((opt)<<25))
100 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>25) & 63)
103 /* Value to be used in mdrun for an infinite cut-off.
104 * Since we need to compare with the cut-off squared,
105 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
107 #define GMX_CUTOFF_INF 1E+18
109 /* enums for the neighborlist type */
110 enum {
111 enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
113 /* OOR is "one over r" -- standard coul */
114 enum {
115 enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
118 enum {
119 egCOULSR, egLJSR, egBHAMSR,
120 egCOUL14, egLJ14, egGB, egNR
122 extern const char *egrp_nm[egNR+1];
124 struct gmx_grppairener_t
126 int nener; /* The number of energy group pairs */
127 real *ener[egNR]; /* Energy terms for each pair of groups */
130 struct gmx_enerdata_t
132 real term[F_NRE]; /* The energies for all different interaction types */
133 struct gmx_grppairener_t grpp;
134 double dvdl_lin[efptNR]; /* Contributions to dvdl with linear lam-dependence */
135 double dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence */
136 int n_lambda;
137 int fep_state; /*current fep state -- just for printing */
138 double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
139 real foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
140 struct gmx_grppairener_t foreign_grpp; /* alternate array for storing foreign lambda energies */
142 /* The idea is that dvdl terms with linear lambda dependence will be added
143 * automatically to enerpart_lambda. Terms with non-linear lambda dependence
144 * should explicitly determine the energies at foreign lambda points
145 * when n_lambda > 0.
148 struct cginfo_mb_t
150 int cg_start;
151 int cg_end;
152 int cg_mod;
153 int *cginfo;
157 /* Forward declaration of type for managing Ewald tables */
158 struct gmx_ewald_tab_t;
160 struct ewald_corr_thread_t;
162 struct t_forcerec {
163 struct interaction_const_t *ic;
165 /* Domain Decomposition */
166 gmx_bool bDomDec;
168 /* PBC stuff */
169 int ePBC;
170 gmx_bool bMolPBC;
171 int rc_scaling;
172 rvec posres_com;
173 rvec posres_comB;
175 gmx_bool use_simd_kernels;
177 /* Interaction for calculated in kernels. In many cases this is similar to
178 * the electrostatics settings in the inputrecord, but the difference is that
179 * these variables always specify the actual interaction in the kernel - if
180 * we are tabulating reaction-field the inputrec will say reaction-field, but
181 * the kernel interaction will say cubic-spline-table. To be safe we also
182 * have a kernel-specific setting for the modifiers - if the interaction is
183 * tabulated we already included the inputrec modification there, so the kernel
184 * modification setting will say 'none' in that case.
186 int nbkernel_elec_interaction;
187 int nbkernel_vdw_interaction;
188 int nbkernel_elec_modifier;
189 int nbkernel_vdw_modifier;
191 /* Use special N*N kernels? */
192 gmx_bool bAllvsAll;
193 /* Private work data */
194 void *AllvsAll_work;
195 void *AllvsAll_workgb;
197 /* Cut-Off stuff.
198 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
200 real rlist;
202 /* Dielectric constant resp. multiplication factor for charges */
203 real zsquare, temp;
204 real epsilon_r, epsilon_rf, epsfac;
206 /* Constants for reaction fields */
207 real kappa, k_rf, c_rf;
209 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
210 double qsum[2];
211 double q2sum[2];
212 double c6sum[2];
213 rvec mu_tot[2];
215 /* Dispersion correction stuff */
216 int eDispCorr;
217 int numAtomsForDispersionCorrection;
218 struct t_forcetable *dispersionCorrectionTable;
220 /* The shift of the shift or user potentials */
221 real enershiftsix;
222 real enershifttwelve;
223 /* Integrated differces for energy and virial with cut-off functions */
224 real enerdiffsix;
225 real enerdifftwelve;
226 real virdiffsix;
227 real virdifftwelve;
228 /* Constant for long range dispersion correction (average dispersion)
229 * for topology A/B ([0]/[1]) */
230 real avcsix[2];
231 /* Constant for long range repulsion term. Relative difference of about
232 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
234 real avctwelve[2];
236 /* Fudge factors */
237 real fudgeQQ;
239 /* Table stuff */
240 gmx_bool bcoultab;
241 gmx_bool bvdwtab;
242 /* The normal tables are in the nblists struct(s) below */
244 struct t_forcetable *pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
246 /* PPPM & Shifting stuff */
247 int coulomb_modifier;
248 real rcoulomb_switch, rcoulomb;
249 real *phi;
251 /* VdW stuff */
252 int vdw_modifier;
253 double reppow;
254 real rvdw_switch, rvdw;
255 real bham_b_max;
257 /* Free energy */
258 int efep;
259 real sc_alphavdw;
260 real sc_alphacoul;
261 int sc_power;
262 real sc_r_power;
263 real sc_sigma6_def;
264 real sc_sigma6_min;
266 /* NS Stuff */
267 int eeltype;
268 int vdwtype;
269 int cg0, hcg;
270 /* solvent_opt contains the enum for the most common solvent
271 * in the system, which will be optimized.
272 * It can be set to esolNO to disable all water optimization */
273 int solvent_opt;
274 int nWatMol;
275 gmx_bool bGrid;
276 gmx_bool bExcl_IntraCGAll_InterCGNone;
277 struct cginfo_mb_t *cginfo_mb;
278 int *cginfo;
279 rvec *cg_cm;
280 int cg_nalloc;
281 rvec *shift_vec;
283 /* The neighborlists including tables */
284 int nnblists;
285 int *gid2nblists;
286 struct t_nblists *nblists;
288 int cutoff_scheme; /* group- or Verlet-style cutoff */
289 gmx_bool bNonbonded; /* true if nonbonded calculations are *not* turned off */
290 struct nonbonded_verlet_t *nbv;
292 /* The wall tables (if used) */
293 int nwall;
294 struct t_forcetable ***wall_tab;
296 /* The number of charge groups participating in do_force_lowlevel */
297 int ncg_force;
298 /* The number of atoms participating in do_force_lowlevel */
299 int natoms_force;
300 /* The number of atoms participating in force and constraints */
301 int natoms_force_constr;
302 /* The allocation size of vectors of size natoms_force */
303 int nalloc_force;
305 /* Forces that should not enter into the virial summation:
306 * PPPM/PME/Ewald/posres
307 * If such forces are present in the system, bF_NoVirSum=TRUE.
309 gmx_bool bF_NoVirSum;
310 #ifdef __cplusplus
311 /* TODO: Replace the pointer by an object once we got rid of C */
312 PaddedRVecVector *forceBufferNoVirialSummation;
313 #else
314 void *forceBufferNoVirialSummation_dummy;
315 #endif
316 /* Pointer that points to forceNoVirialSummation when virial is calcaluted,
317 * points to the normal force vector when the virial is not requested
318 * or when bF_NoVirSum == FALSE.
320 #ifdef __cplusplus
321 PaddedRVecVector *f_novirsum;
322 #else
323 void *f_novirsum_xdummy;
324 #endif
326 /* Long-range forces and virial for PPPM/PME/Ewald */
327 struct gmx_pme_t *pmedata;
328 int ljpme_combination_rule;
329 tensor vir_el_recip;
330 tensor vir_lj_recip;
332 /* PME/Ewald stuff */
333 gmx_bool bEwald;
334 real ewaldcoeff_q;
335 real ewaldcoeff_lj;
336 struct gmx_ewald_tab_t *ewald_table;
338 /* Virial Stuff */
339 rvec *fshift;
340 rvec vir_diag_posres;
341 dvec vir_wall_z;
343 /* Non bonded Parameter lists */
344 int ntype; /* Number of atom types */
345 gmx_bool bBHAM;
346 real *nbfp;
347 real *ljpme_c6grid; /* C6-values used on grid in LJPME */
349 /* Energy group pair flags */
350 int *egp_flags;
352 /* Shell molecular dynamics flexible constraints */
353 real fc_stepsize;
355 /* Generalized born implicit solvent */
356 gmx_bool bGB;
357 /* Generalized born stuff */
358 real gb_epsilon_solvent;
359 /* Table data for GB */
360 struct t_forcetable *gbtab;
361 /* VdW radius for each atomtype (dim is thus ntype) */
362 real *atype_radius;
363 /* Effective radius (derived from effective volume) for each type */
364 real *atype_vol;
365 /* Implicit solvent - surface tension for each atomtype */
366 real *atype_surftens;
367 /* Implicit solvent - radius for GB calculation */
368 real *atype_gb_radius;
369 /* Implicit solvent - overlap for HCT model */
370 real *atype_S_hct;
371 /* Generalized born interaction data */
372 struct gmx_genborn_t *born;
374 /* Table scale for GB */
375 real gbtabscale;
376 /* Table range for GB */
377 real gbtabr;
378 /* GB neighborlists (the sr list will contain for each atom all other atoms
379 * (for use in the SA calculation) and the lr list will contain
380 * for each atom all atoms 1-4 or greater (for use in the GB calculation)
382 struct t_nblist *gblist_sr;
383 struct t_nblist *gblist_lr;
384 struct t_nblist *gblist;
386 /* Inverse square root of the Born radii for implicit solvent */
387 real *invsqrta;
388 /* Derivatives of the potential with respect to the Born radii */
389 real *dvda;
390 /* Derivatives of the Born radii with respect to coordinates */
391 real *dadx;
392 real *dadx_rawptr;
393 int nalloc_dadx; /* Allocated size of dadx */
395 /* If > 0 signals Test Particle Insertion,
396 * the value is the number of atoms of the molecule to insert
397 * Only the energy difference due to the addition of the last molecule
398 * should be calculated.
400 gmx_bool n_tpi;
402 /* Neighbor searching stuff */
403 struct gmx_ns_t *ns;
405 /* QMMM stuff */
406 gmx_bool bQMMM;
407 struct t_QMMMrec *qr;
409 /* QM-MM neighborlists */
410 struct t_nblist *QMMMlist;
412 /* Limit for printing large forces, negative is don't print */
413 real print_force;
415 /* coarse load balancing time measurement */
416 double t_fnbf;
417 double t_wait;
418 int timesteps;
420 /* User determined parameters, copied from the inputrec */
421 int userint1;
422 int userint2;
423 int userint3;
424 int userint4;
425 real userreal1;
426 real userreal2;
427 real userreal3;
428 real userreal4;
430 /* Pointer to struct for managing threading of bonded force calculation */
431 struct bonded_threading_t *bonded_threading;
433 /* Ewald correction thread local virial and energy data */
434 int nthread_ewc;
435 struct ewald_corr_thread_t *ewc_t;
437 struct ForceProviders *forceProviders;
440 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
441 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
442 * in the code, but beware if you are using these macros externally.
444 #define C6(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))]
445 #define C12(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
446 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))]
447 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
448 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]
450 #ifdef __cplusplus
452 #endif
454 #endif