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.
50 /* check kernel/toppush.c when you change these numbers */
52 #define MAXFORCEPARAM 12
56 typedef atom_id t_iatom
;
58 /* this MUST correspond to the
59 t_interaction_function[F_NRE] in gmxlib/ifunc.c */
136 F_VTEMP_NOLONGERUSED
,
146 F_DVDL_TEMPERATURE
, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
147 F_NRE
/* This number is for the total number of energies */
150 #define IS_RESTRAINT_TYPE(ifunc) (((ifunc == F_POSRES) || (ifunc == F_DISRES) || (ifunc == F_RESTRBONDS) || (ifunc == F_DISRESVIOL) || (ifunc == F_ORIRES) || (ifunc == F_ORIRESDEV) || (ifunc == F_ANGRES) || (ifunc == F_ANGRESZ) || (ifunc == F_DIHRES)))
152 /* A macro for checking if ftype is an explicit pair-listed LJ or COULOMB
154 * bonded LJ (usually 1-4), or special listed non-bonded for FEP.
156 #define IS_LISTED_LJ_C(ftype) ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB)
160 /* Some parameters have A and B values for free energy calculations.
161 * The B values are not used for regular simulations of course.
162 * Free Energy for nonbondeds can be computed by changing the atom type.
163 * The harmonic type is used for all harmonic potentials:
164 * bonds, angles and improper dihedrals
170 real rA
, krA
, rB
, krB
;
173 real klinA
, aA
, klinB
, aB
;
176 real lowA
, up1A
, up2A
, kA
, lowB
, up1B
, up2B
, kB
;
178 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
189 real r1e
, r2e
, r3e
, krt
;
192 real thetaA
, kthetaA
, r13A
, kUBA
, thetaB
, kthetaB
, r13B
, kUBB
;
201 real alpha
, drcut
, khyp
;
204 real al_x
, al_y
, al_z
, rOH
, rHH
, rOD
;
207 real a
, alpha1
, alpha2
, rfac
;
213 real c6A
, c12A
, c6B
, c12B
;
216 real fqq
, qi
, qj
, c6
, c12
;
219 real qi
, qj
, c6
, c12
;
221 /* Proper dihedrals can not have different multiplicity when
222 * doing free energy calculations, because the potential would not
223 * be periodic anymore.
226 real phiA
, cpA
; int mult
; real phiB
, cpB
;
231 /* Settle can not be used for Free energy calculations of water bond geometry.
232 * Use shake (or lincs) instead if you have to change the water bonds.
238 real b0A
, cbA
, betaA
, b0B
, cbB
, betaB
;
241 real pos0A
[DIM
], fcA
[DIM
], pos0B
[DIM
], fcB
[DIM
];
244 real rbcA
[NR_RBDIHS
], rbcB
[NR_RBDIHS
];
247 real a
, b
, c
, d
, e
, f
;
252 /* NOTE: npair is only set after reading the tpx file */
254 real low
, up1
, up2
, kfac
; int type
, label
, npair
;
257 real phiA
, dphiA
, kfacA
, phiB
, dphiB
, kfacB
;
260 int ex
, power
, label
; real c
, obs
, kfac
;
263 int table
; real kA
; real kB
;
266 real sar
, st
, pi
, gbr
, bmlt
;
272 real buf
[MAXFORCEPARAM
];
273 } generic
; /* Conversion */
276 typedef int t_functype
;
279 * The nonperturbed/perturbed interactions are now separated (sorted) in the
280 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
281 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
293 * The struct t_ilist defines a list of atoms with their interactions.
294 * General field description:
296 * the size (nr elements) of the interactions array (iatoms[]).
298 * specifies which atoms are involved in an interaction of a certain
299 * type. The layout of this array is as follows:
301 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
302 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
303 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
305 * So for interaction type type1 3 atoms are needed, and for type2 and
306 * type3 only 2. The type identifier is used to select the function to
307 * calculate the interaction and its actual parameters. This type
308 * identifier is an index in a params[] and functype[] array.
313 real
*cmap
; /* Has length 4*grid_spacing*grid_spacing, */
314 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
319 int ngrid
; /* Number of allocated cmap (cmapdata_t ) grids */
320 int grid_spacing
; /* Grid spacing */
321 cmapdata_t
*cmapdata
; /* Pointer to grid with actual, pre-interpolated data */
329 t_functype
*functype
;
331 double reppow
; /* The repulsion power for VdW: C12*r^-reppow */
332 real fudgeQQ
; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
333 gmx_cmap_t cmap_grid
; /* The dihedral correction maps */
337 ilsortUNKNOWN
, ilsortNO_FE
, ilsortFE_UNSORTED
, ilsortFE_SORTED
344 t_functype
*functype
;
347 gmx_cmap_t cmap_grid
;
348 t_iparams
*iparams_posres
;
349 int iparams_posres_nalloc
;
356 * The struct t_idef defines all the interactions for the complete
357 * simulation. The structure is setup in such a way that the multinode
358 * version of the program can use it as easy as the single node version.
359 * General field description:
361 * defines the number of elements in functype[] and param[].
363 * the node id (if parallel machines)
365 * the number of atomtypes
366 * t_functype *functype
367 * array of length ntypes, defines for every force type what type of
368 * function to use. Every "bond" with the same function but different
369 * force parameters is a different force type. The type identifier in the
370 * forceatoms[] array is an index in this array.
372 * array of length ntypes, defines the parameters for every interaction
373 * type. The type identifier in the actual interaction list
374 * (ilist[ftype].iatoms[]) is an index in this array.
375 * gmx_cmap_t cmap_grid
376 * the grid for the dihedral pair correction maps.
377 * t_iparams *iparams_posres
378 * defines the parameters for position restraints only.
379 * Position restraints are the only interactions that have different
380 * parameters (reference positions) for different molecules
381 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
383 * The list of interactions for each type. Note that some,
384 * such as LJ and COUL will have 0 entries.
388 int n
; /* n+1 is the number of points */
389 real scale
; /* distance between two points */
390 real
*data
; /* the actual table data, per point there are 4 numbers */