Properly finalize MPI on mdrun -version. Fixes #1313
[gromacs.git] / include / types / forcerec.h
blob2e78634a43718cd815d38a402a82ff97bd01bdcf
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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "ns.h"
40 #include "genborn.h"
41 #include "qmmmrec.h"
42 #include "idef.h"
43 #include "nb_verlet.h"
44 #include "interaction_const.h"
45 #include "hw_info.h"
47 #ifdef __cplusplus
48 extern "C" {
49 #endif
50 #if 0
51 } /* fixes auto-indentation problems */
52 #endif
54 /* Abstract type for PME that is defined only in the routine that use them. */
55 typedef struct gmx_pme *gmx_pme_t;
59 /* Structure describing the data in a single table */
60 typedef struct
62 enum gmx_table_interaction interaction; /* Types of interactions stored in this table */
63 enum gmx_table_format format; /* Interpolation type and data format */
65 real r; /* range of the table */
66 int n; /* n+1 is the number of table points */
67 real scale; /* distance (nm) between two table points */
68 real scale_exp; /* distance for exponential part of VdW table, not always used */
69 real * data; /* the actual table data */
71 /* Some information about the table layout. This can also be derived from the interpolation
72 * type and the table interactions, but it is convenient to have here for sanity checks, and it makes it
73 * much easier to access the tables in the nonbonded kernels when we can set the data from variables.
74 * It is always true that stride = formatsize*ninteractions
76 int formatsize; /* Number of fp variables for each table point (1 for F, 2 for VF, 4 for YFGH, etc.) */
77 int ninteractions; /* Number of interactions in table, 1 for coul-only, 3 for coul+rep+disp. */
78 int stride; /* Distance to next table point (number of fp variables per table point in total) */
79 } t_forcetable;
81 typedef struct
83 t_forcetable table_elec;
84 t_forcetable table_vdw;
85 t_forcetable table_elec_vdw;
87 /* The actual neighbor lists, short and long range, see enum above
88 * for definition of neighborlist indices.
90 t_nblist nlist_sr[eNL_NR];
91 t_nblist nlist_lr[eNL_NR];
92 } t_nblists;
94 /* macros for the cginfo data in forcerec */
95 /* The maximum cg size in cginfo is 63
96 * because we only have space for 6 bits in cginfo,
97 * this cg size entry is actually only read with domain decomposition.
98 * But there is a smaller limit due to the t_excl data structure
99 * which is defined in nblist.h.
101 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~65535) | (gid) )
102 #define GET_CGINFO_GID(cgi) ( (cgi) & 65535)
103 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
104 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
105 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
106 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
107 #define SET_CGINFO_SOLOPT(cgi, opt) (cgi) = (((cgi) & ~(3<<18)) | ((opt)<<18))
108 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 3)
109 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1<<20))
110 #define GET_CGINFO_CONSTR(cgi) ( (cgi) & (1<<20))
111 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1<<21))
112 #define GET_CGINFO_SETTLE(cgi) ( (cgi) & (1<<21))
113 /* This bit is only used with bBondComm in the domain decomposition */
114 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1<<22))
115 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi) & (1<<22))
116 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1<<23))
117 #define GET_CGINFO_HAS_VDW(cgi) ( (cgi) & (1<<23))
118 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1<<24))
119 #define GET_CGINFO_HAS_Q(cgi) ( (cgi) & (1<<24))
120 #define SET_CGINFO_NATOMS(cgi, opt) (cgi) = (((cgi) & ~(63<<25)) | ((opt)<<25))
121 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>25) & 63)
124 /* Value to be used in mdrun for an infinite cut-off.
125 * Since we need to compare with the cut-off squared,
126 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
128 #define GMX_CUTOFF_INF 1E+18
130 /* enums for the neighborlist type */
131 enum {
132 enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
134 /* OOR is "one over r" -- standard coul */
135 enum {
136 enbcoulNONE, enbcoulOOR, enbcoulRF, enbcoulTAB, enbcoulGB, enbcoulFEWALD, enbcoulNR
139 enum {
140 egCOULSR, egLJSR, egBHAMSR, egCOULLR, egLJLR, egBHAMLR,
141 egCOUL14, egLJ14, egGB, egNR
144 typedef struct {
145 int nener; /* The number of energy group pairs */
146 real *ener[egNR]; /* Energy terms for each pair of groups */
147 } gmx_grppairener_t;
149 typedef struct {
150 real term[F_NRE]; /* The energies for all different interaction types */
151 gmx_grppairener_t grpp;
152 double dvdl_lin[efptNR]; /* Contributions to dvdl with linear lam-dependence */
153 double dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence */
154 int n_lambda;
155 int fep_state; /*current fep state -- just for printing */
156 double *enerpart_lambda; /* Partial energy for lambda and flambda[] */
157 real foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
158 gmx_grppairener_t foreign_grpp; /* alternate array for storing foreign lambda energies */
159 } gmx_enerdata_t;
160 /* The idea is that dvdl terms with linear lambda dependence will be added
161 * automatically to enerpart_lambda. Terms with non-linear lambda dependence
162 * should explicitly determine the energies at foreign lambda points
163 * when n_lambda > 0.
166 typedef struct {
167 int cg_start;
168 int cg_end;
169 int cg_mod;
170 int *cginfo;
171 } cginfo_mb_t;
174 /* ewald table type */
175 typedef struct ewald_tab *ewald_tab_t;
177 typedef struct {
178 rvec *f;
179 int f_nalloc;
180 unsigned red_mask; /* Mask for marking which parts of f are filled */
181 rvec *fshift;
182 real ener[F_NRE];
183 gmx_grppairener_t grpp;
184 real Vcorr;
185 real dvdl[efptNR];
186 tensor vir;
187 } f_thread_t;
189 typedef struct {
190 interaction_const_t *ic;
192 /* Domain Decomposition */
193 gmx_bool bDomDec;
195 /* PBC stuff */
196 int ePBC;
197 gmx_bool bMolPBC;
198 int rc_scaling;
199 rvec posres_com;
200 rvec posres_comB;
202 const gmx_hw_info_t *hwinfo;
203 gmx_bool use_cpu_acceleration;
205 /* Interaction for calculated in kernels. In many cases this is similar to
206 * the electrostatics settings in the inputrecord, but the difference is that
207 * these variables always specify the actual interaction in the kernel - if
208 * we are tabulating reaction-field the inputrec will say reaction-field, but
209 * the kernel interaction will say cubic-spline-table. To be safe we also
210 * have a kernel-specific setting for the modifiers - if the interaction is
211 * tabulated we already included the inputrec modification there, so the kernel
212 * modification setting will say 'none' in that case.
214 int nbkernel_elec_interaction;
215 int nbkernel_vdw_interaction;
216 int nbkernel_elec_modifier;
217 int nbkernel_vdw_modifier;
219 /* Use special N*N kernels? */
220 gmx_bool bAllvsAll;
221 /* Private work data */
222 void *AllvsAll_work;
223 void *AllvsAll_workgb;
225 /* Cut-Off stuff.
226 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
228 real rlist, rlistlong;
230 /* Dielectric constant resp. multiplication factor for charges */
231 real zsquare, temp;
232 real epsilon_r, epsilon_rf, epsfac;
234 /* Constants for reaction fields */
235 real kappa, k_rf, c_rf;
237 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
238 double qsum[2];
239 double q2sum[2];
240 rvec mu_tot[2];
242 /* Dispersion correction stuff */
243 int eDispCorr;
245 /* The shift of the shift or user potentials */
246 real enershiftsix;
247 real enershifttwelve;
248 /* Integrated differces for energy and virial with cut-off functions */
249 real enerdiffsix;
250 real enerdifftwelve;
251 real virdiffsix;
252 real virdifftwelve;
253 /* Constant for long range dispersion correction (average dispersion)
254 * for topology A/B ([0]/[1]) */
255 real avcsix[2];
256 /* Constant for long range repulsion term. Relative difference of about
257 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
259 real avctwelve[2];
261 /* Fudge factors */
262 real fudgeQQ;
264 /* Table stuff */
265 gmx_bool bcoultab;
266 gmx_bool bvdwtab;
267 /* The normal tables are in the nblists struct(s) below */
268 t_forcetable tab14; /* for 1-4 interactions only */
270 /* PPPM & Shifting stuff */
271 int coulomb_modifier;
272 real rcoulomb_switch, rcoulomb;
273 real *phi;
275 /* VdW stuff */
276 int vdw_modifier;
277 double reppow;
278 real rvdw_switch, rvdw;
279 real bham_b_max;
281 /* Free energy */
282 int efep;
283 real sc_alphavdw;
284 real sc_alphacoul;
285 int sc_power;
286 real sc_r_power;
287 real sc_sigma6_def;
288 real sc_sigma6_min;
289 gmx_bool bSepDVDL;
291 /* NS Stuff */
292 int eeltype;
293 int vdwtype;
294 int cg0, hcg;
295 /* solvent_opt contains the enum for the most common solvent
296 * in the system, which will be optimized.
297 * It can be set to esolNO to disable all water optimization */
298 int solvent_opt;
299 int nWatMol;
300 gmx_bool bGrid;
301 gmx_bool bExcl_IntraCGAll_InterCGNone;
302 cginfo_mb_t *cginfo_mb;
303 int *cginfo;
304 rvec *cg_cm;
305 int cg_nalloc;
306 rvec *shift_vec;
308 /* The neighborlists including tables */
309 int nnblists;
310 int *gid2nblists;
311 t_nblists *nblists;
313 int cutoff_scheme; /* group- or Verlet-style cutoff */
314 gmx_bool bNonbonded; /* true if nonbonded calculations are *not* turned off */
315 nonbonded_verlet_t *nbv;
317 /* The wall tables (if used) */
318 int nwall;
319 t_forcetable **wall_tab;
321 /* The number of charge groups participating in do_force_lowlevel */
322 int ncg_force;
323 /* The number of atoms participating in do_force_lowlevel */
324 int natoms_force;
325 /* The number of atoms participating in force and constraints */
326 int natoms_force_constr;
327 /* The allocation size of vectors of size natoms_force */
328 int nalloc_force;
330 /* Twin Range stuff, f_twin has size natoms_force */
331 gmx_bool bTwinRange;
332 int nlr;
333 rvec *f_twin;
335 /* Forces that should not enter into the virial summation:
336 * PPPM/PME/Ewald/posres
338 gmx_bool bF_NoVirSum;
339 int f_novirsum_n;
340 int f_novirsum_nalloc;
341 rvec *f_novirsum_alloc;
342 /* Pointer that points to f_novirsum_alloc when pressure is calcaluted,
343 * points to the normal force vectors wen pressure is not requested.
345 rvec *f_novirsum;
347 /* Long-range forces and virial for PPPM/PME/Ewald */
348 gmx_pme_t pmedata;
349 tensor vir_el_recip;
351 /* PME/Ewald stuff */
352 gmx_bool bEwald;
353 real ewaldcoeff;
354 ewald_tab_t ewald_table;
356 /* Virial Stuff */
357 rvec *fshift;
358 rvec vir_diag_posres;
359 dvec vir_wall_z;
361 /* Non bonded Parameter lists */
362 int ntype; /* Number of atom types */
363 gmx_bool bBHAM;
364 real *nbfp;
366 /* Energy group pair flags */
367 int *egp_flags;
369 /* xmdrun flexible constraints */
370 real fc_stepsize;
372 /* Generalized born implicit solvent */
373 gmx_bool bGB;
374 /* Generalized born stuff */
375 real gb_epsilon_solvent;
376 /* Table data for GB */
377 t_forcetable gbtab;
378 /* VdW radius for each atomtype (dim is thus ntype) */
379 real *atype_radius;
380 /* Effective radius (derived from effective volume) for each type */
381 real *atype_vol;
382 /* Implicit solvent - surface tension for each atomtype */
383 real *atype_surftens;
384 /* Implicit solvent - radius for GB calculation */
385 real *atype_gb_radius;
386 /* Implicit solvent - overlap for HCT model */
387 real *atype_S_hct;
388 /* Generalized born interaction data */
389 gmx_genborn_t *born;
391 /* Table scale for GB */
392 real gbtabscale;
393 /* Table range for GB */
394 real gbtabr;
395 /* GB neighborlists (the sr list will contain for each atom all other atoms
396 * (for use in the SA calculation) and the lr list will contain
397 * for each atom all atoms 1-4 or greater (for use in the GB calculation)
399 t_nblist gblist_sr;
400 t_nblist gblist_lr;
401 t_nblist gblist;
403 /* Inverse square root of the Born radii for implicit solvent */
404 real *invsqrta;
405 /* Derivatives of the potential with respect to the Born radii */
406 real *dvda;
407 /* Derivatives of the Born radii with respect to coordinates */
408 real *dadx;
409 real *dadx_rawptr;
410 int nalloc_dadx; /* Allocated size of dadx */
412 /* If > 0 signals Test Particle Insertion,
413 * the value is the number of atoms of the molecule to insert
414 * Only the energy difference due to the addition of the last molecule
415 * should be calculated.
417 gmx_bool n_tpi;
419 /* Neighbor searching stuff */
420 gmx_ns_t ns;
422 /* QMMM stuff */
423 gmx_bool bQMMM;
424 t_QMMMrec *qr;
426 /* QM-MM neighborlists */
427 t_nblist QMMMlist;
429 /* Limit for printing large forces, negative is don't print */
430 real print_force;
432 /* coarse load balancing time measurement */
433 double t_fnbf;
434 double t_wait;
435 int timesteps;
437 /* parameter needed for AdResS simulation */
438 int adress_type;
439 gmx_bool badress_tf_full_box;
440 real adress_const_wf;
441 real adress_ex_width;
442 real adress_hy_width;
443 int adress_icor;
444 int adress_site;
445 rvec adress_refs;
446 int n_adress_tf_grps;
447 int * adress_tf_table_index;
448 int *adress_group_explicit;
449 t_forcetable * atf_tabs;
450 real adress_ex_forcecap;
451 gmx_bool adress_do_hybridpairs;
453 /* User determined parameters, copied from the inputrec */
454 int userint1;
455 int userint2;
456 int userint3;
457 int userint4;
458 real userreal1;
459 real userreal2;
460 real userreal3;
461 real userreal4;
463 /* Thread local force and energy data */
464 /* FIXME move to bonded_thread_data_t */
465 int nthreads;
466 int red_ashift;
467 int red_nblock;
468 f_thread_t *f_t;
470 /* Exclusion load distribution over the threads */
471 int *excl_load;
472 } t_forcerec;
474 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
475 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
476 * in the code, but beware if you are using these macros externally.
478 #define C6(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))]
479 #define C12(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
480 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))]
481 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
482 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]
484 #ifdef __cplusplus
486 #endif