Created paddedvector.h
[gromacs.git] / src / gromacs / mdtypes / forcerec.h
blob7cabf08377e11c2f7eed9c21137adb6980b7b6ec
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, 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 #ifdef __cplusplus
52 extern "C" {
54 struct t_mdatoms;
55 struct t_commrec;
57 /*! \libinternal \brief
58 * Interface for a component that provides forces during MD.
60 * This is typically part of a larger structure/class managing its own
61 * data, such that it has the information on what to do stored locally.
62 * \todo Implement returning of energy and dH/dlambda.
63 * \inlibraryapi
65 struct IForceProvider
67 public:
68 /*! \brief Compute forces
70 * \todo This is specific for electric fields and needs to be generalized.
71 * \param[in] cr Communication record for parallel operations
72 * \param[in] mdatoms Atom information
73 * \param[inout] force The forces
74 * \param[in] t The actual time in the simulation (ps)
76 virtual void calculateForces(const t_commrec *cr,
77 const t_mdatoms *mdatoms,
78 PaddedRVecVector *force,
79 double t) = 0;
81 protected:
82 ~IForceProvider() {}
84 #endif
86 /* Abstract type for PME that is defined only in the routine that use them. */
87 struct gmx_genborn_t;
88 struct gmx_ns_t;
89 struct gmx_pme_t;
90 struct nonbonded_verlet_t;
91 struct bonded_threading_t;
92 struct t_forcetable;
93 struct t_nblist;
94 struct t_nblists;
95 struct t_QMMMrec;
96 struct gmx_hw_info_t;
97 struct gmx_gpu_opt_t;
99 /* macros for the cginfo data in forcerec
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.
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.
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)
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).
139 #define GMX_CUTOFF_INF 1E+18
141 /* enums for the neighborlist type */
142 enum {
143 enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
145 /* OOR is "one over r" -- standard coul */
146 enum {
147 enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
150 enum {
151 egCOULSR, egLJSR, egBHAMSR,
152 egCOUL14, egLJ14, egGB, egNR
154 extern const char *egrp_nm[egNR+1];
156 typedef struct gmx_grppairener_t {
157 int nener; /* The number of energy group pairs */
158 real *ener[egNR]; /* Energy terms for each pair of groups */
159 } gmx_grppairener_t;
161 typedef struct gmx_enerdata_t {
162 real term[F_NRE]; /* The energies for all different interaction types */
163 gmx_grppairener_t grpp;
164 double dvdl_lin[efptNR]; /* Contributions to dvdl with linear lam-dependence */
165 double dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence */
166 int n_lambda;
167 int fep_state; /*current fep state -- just for printing */
168 double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
169 real foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
170 gmx_grppairener_t foreign_grpp; /* alternate array for storing foreign lambda energies */
171 } gmx_enerdata_t;
172 /* The idea is that dvdl terms with linear lambda dependence will be added
173 * automatically to enerpart_lambda. Terms with non-linear lambda dependence
174 * should explicitly determine the energies at foreign lambda points
175 * when n_lambda > 0.
178 typedef struct {
179 int cg_start;
180 int cg_end;
181 int cg_mod;
182 int *cginfo;
183 } cginfo_mb_t;
186 /* Forward declaration of type for managing Ewald tables */
187 struct gmx_ewald_tab_t;
189 typedef struct ewald_corr_thread_t ewald_corr_thread_t;
191 typedef struct t_forcerec {
192 interaction_const_t *ic;
194 /* Domain Decomposition */
195 gmx_bool bDomDec;
197 /* PBC stuff */
198 int ePBC;
199 gmx_bool bMolPBC;
200 int rc_scaling;
201 rvec posres_com;
202 rvec posres_comB;
204 const struct gmx_hw_info_t *hwinfo;
205 const struct gmx_gpu_opt_t *gpu_opt;
206 gmx_bool use_simd_kernels;
208 /* Interaction for calculated in kernels. In many cases this is similar to
209 * the electrostatics settings in the inputrecord, but the difference is that
210 * these variables always specify the actual interaction in the kernel - if
211 * we are tabulating reaction-field the inputrec will say reaction-field, but
212 * the kernel interaction will say cubic-spline-table. To be safe we also
213 * have a kernel-specific setting for the modifiers - if the interaction is
214 * tabulated we already included the inputrec modification there, so the kernel
215 * modification setting will say 'none' in that case.
217 int nbkernel_elec_interaction;
218 int nbkernel_vdw_interaction;
219 int nbkernel_elec_modifier;
220 int nbkernel_vdw_modifier;
222 /* Use special N*N kernels? */
223 gmx_bool bAllvsAll;
224 /* Private work data */
225 void *AllvsAll_work;
226 void *AllvsAll_workgb;
228 /* Cut-Off stuff.
229 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
231 real rlist;
233 /* Dielectric constant resp. multiplication factor for charges */
234 real zsquare, temp;
235 real epsilon_r, epsilon_rf, epsfac;
237 /* Constants for reaction fields */
238 real kappa, k_rf, c_rf;
240 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
241 double qsum[2];
242 double q2sum[2];
243 double c6sum[2];
244 rvec mu_tot[2];
246 /* Dispersion correction stuff */
247 int eDispCorr;
248 int numAtomsForDispersionCorrection;
249 struct t_forcetable *dispersionCorrectionTable;
251 /* The shift of the shift or user potentials */
252 real enershiftsix;
253 real enershifttwelve;
254 /* Integrated differces for energy and virial with cut-off functions */
255 real enerdiffsix;
256 real enerdifftwelve;
257 real virdiffsix;
258 real virdifftwelve;
259 /* Constant for long range dispersion correction (average dispersion)
260 * for topology A/B ([0]/[1]) */
261 real avcsix[2];
262 /* Constant for long range repulsion term. Relative difference of about
263 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
265 real avctwelve[2];
267 /* Fudge factors */
268 real fudgeQQ;
270 /* Table stuff */
271 gmx_bool bcoultab;
272 gmx_bool bvdwtab;
273 /* The normal tables are in the nblists struct(s) below */
275 struct t_forcetable *pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
277 /* PPPM & Shifting stuff */
278 int coulomb_modifier;
279 real rcoulomb_switch, rcoulomb;
280 real *phi;
282 /* VdW stuff */
283 int vdw_modifier;
284 double reppow;
285 real rvdw_switch, rvdw;
286 real bham_b_max;
288 /* Free energy */
289 int efep;
290 real sc_alphavdw;
291 real sc_alphacoul;
292 int sc_power;
293 real sc_r_power;
294 real sc_sigma6_def;
295 real sc_sigma6_min;
297 /* NS Stuff */
298 int eeltype;
299 int vdwtype;
300 int cg0, hcg;
301 /* solvent_opt contains the enum for the most common solvent
302 * in the system, which will be optimized.
303 * It can be set to esolNO to disable all water optimization */
304 int solvent_opt;
305 int nWatMol;
306 gmx_bool bGrid;
307 gmx_bool bExcl_IntraCGAll_InterCGNone;
308 cginfo_mb_t *cginfo_mb;
309 int *cginfo;
310 rvec *cg_cm;
311 int cg_nalloc;
312 rvec *shift_vec;
314 /* The neighborlists including tables */
315 int nnblists;
316 int *gid2nblists;
317 struct t_nblists *nblists;
319 int cutoff_scheme; /* group- or Verlet-style cutoff */
320 gmx_bool bNonbonded; /* true if nonbonded calculations are *not* turned off */
321 struct nonbonded_verlet_t *nbv;
323 /* The wall tables (if used) */
324 int nwall;
325 struct t_forcetable ***wall_tab;
327 /* The number of charge groups participating in do_force_lowlevel */
328 int ncg_force;
329 /* The number of atoms participating in do_force_lowlevel */
330 int natoms_force;
331 /* The number of atoms participating in force and constraints */
332 int natoms_force_constr;
333 /* The allocation size of vectors of size natoms_force */
334 int nalloc_force;
336 /* Forces that should not enter into the virial summation:
337 * PPPM/PME/Ewald/posres
338 * If such forces are present in the system, bF_NoVirSum=TRUE.
340 gmx_bool bF_NoVirSum;
341 #ifdef __cplusplus
342 /* TODO: Replace the pointer by an object once we got rid of C */
343 PaddedRVecVector *forceBufferNoVirialSummation;
344 #else
345 void *forceBufferNoVirialSummation_dummy;
346 #endif
347 /* Pointer that points to forceNoVirialSummation when virial is calcaluted,
348 * points to the normal force vector when the virial is not requested
349 * or when bF_NoVirSum == FALSE.
351 #ifdef __cplusplus
352 PaddedRVecVector *f_novirsum;
353 #else
354 void *f_novirsum_xdummy;
355 #endif
357 /* Long-range forces and virial for PPPM/PME/Ewald */
358 struct gmx_pme_t *pmedata;
359 int ljpme_combination_rule;
360 tensor vir_el_recip;
361 tensor vir_lj_recip;
363 /* PME/Ewald stuff */
364 gmx_bool bEwald;
365 real ewaldcoeff_q;
366 real ewaldcoeff_lj;
367 struct gmx_ewald_tab_t *ewald_table;
369 /* Virial Stuff */
370 rvec *fshift;
371 rvec vir_diag_posres;
372 dvec vir_wall_z;
374 /* Non bonded Parameter lists */
375 int ntype; /* Number of atom types */
376 gmx_bool bBHAM;
377 real *nbfp;
378 real *ljpme_c6grid; /* C6-values used on grid in LJPME */
380 /* Energy group pair flags */
381 int *egp_flags;
383 /* Shell molecular dynamics flexible constraints */
384 real fc_stepsize;
386 /* Generalized born implicit solvent */
387 gmx_bool bGB;
388 /* Generalized born stuff */
389 real gb_epsilon_solvent;
390 /* Table data for GB */
391 struct t_forcetable *gbtab;
392 /* VdW radius for each atomtype (dim is thus ntype) */
393 real *atype_radius;
394 /* Effective radius (derived from effective volume) for each type */
395 real *atype_vol;
396 /* Implicit solvent - surface tension for each atomtype */
397 real *atype_surftens;
398 /* Implicit solvent - radius for GB calculation */
399 real *atype_gb_radius;
400 /* Implicit solvent - overlap for HCT model */
401 real *atype_S_hct;
402 /* Generalized born interaction data */
403 struct gmx_genborn_t *born;
405 /* Table scale for GB */
406 real gbtabscale;
407 /* Table range for GB */
408 real gbtabr;
409 /* GB neighborlists (the sr list will contain for each atom all other atoms
410 * (for use in the SA calculation) and the lr list will contain
411 * for each atom all atoms 1-4 or greater (for use in the GB calculation)
413 struct t_nblist *gblist_sr;
414 struct t_nblist *gblist_lr;
415 struct t_nblist *gblist;
417 /* Inverse square root of the Born radii for implicit solvent */
418 real *invsqrta;
419 /* Derivatives of the potential with respect to the Born radii */
420 real *dvda;
421 /* Derivatives of the Born radii with respect to coordinates */
422 real *dadx;
423 real *dadx_rawptr;
424 int nalloc_dadx; /* Allocated size of dadx */
426 /* If > 0 signals Test Particle Insertion,
427 * the value is the number of atoms of the molecule to insert
428 * Only the energy difference due to the addition of the last molecule
429 * should be calculated.
431 gmx_bool n_tpi;
433 /* Neighbor searching stuff */
434 struct gmx_ns_t *ns;
436 /* QMMM stuff */
437 gmx_bool bQMMM;
438 struct t_QMMMrec *qr;
440 /* QM-MM neighborlists */
441 struct t_nblist *QMMMlist;
443 /* Limit for printing large forces, negative is don't print */
444 real print_force;
446 /* coarse load balancing time measurement */
447 double t_fnbf;
448 double t_wait;
449 int timesteps;
451 /* User determined parameters, copied from the inputrec */
452 int userint1;
453 int userint2;
454 int userint3;
455 int userint4;
456 real userreal1;
457 real userreal2;
458 real userreal3;
459 real userreal4;
461 /* Pointer to struct for managing threading of bonded force calculation */
462 struct bonded_threading_t *bonded_threading;
464 /* Ewald correction thread local virial and energy data */
465 int nthread_ewc;
466 ewald_corr_thread_t *ewc_t;
468 struct IForceProvider *efield;
469 } t_forcerec;
471 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
472 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
473 * in the code, but beware if you are using these macros externally.
475 #define C6(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))]
476 #define C12(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
477 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))]
478 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
479 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]
481 #ifdef __cplusplus
483 #endif
484 #endif