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_STATE_H
38 #define GMX_MDTYPES_STATE_H
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/mdtypes/md_enums.h"
42 #include "gromacs/utility/basedefinitions.h"
43 #include "gromacs/utility/real.h"
45 struct energyhistory_t
;
48 * The t_state struct should contain all the (possibly) non-static
49 * information required to define the state of the system.
50 * Currently the random seeds for SD and BD are missing.
53 /* These enums are used in flags as (1<<est...).
54 * The order of these enums should not be changed,
55 * since that affects the checkpoint (.cpt) file format.
59 estBOX
, estBOX_REL
, estBOXV
, estPRES_PREV
, estNH_XI
, estTC_INT
,
60 estX
, estV
, est_SDX_NOTSUPPORTED
, estCGP
,
61 estLD_RNG
, estLD_RNGI
,
62 estDISRE_INITF
, estDISRE_RM3TAV
,
63 estORIRE_INITF
, estORIRE_DTAV
,
64 estSVIR_PREV
, estNH_VXI
, estVETA
, estVOL0
, estNHPRES_XI
, estNHPRES_VXI
, estFVIR_PREV
,
65 estFEPSTATE
, estMC_RNG
, estMC_RNGI
,
69 #define EST_DISTR(e) (!(((e) >= estLAMBDA && (e) <= estTC_INT) || ((e) >= estSVIR_PREV && (e) <= estMC_RNGI)))
71 /* The names of the state entries, defined in src/gmxlib/checkpoint.c */
72 extern const char *est_names
[estNR
];
74 typedef struct history_t
76 real disre_initf
; /* The scaling factor for initializing the time av. */
77 int ndisrepairs
; /* The number of distance restraints */
78 real
*disre_rm3tav
; /* The r^-3 time averaged pair distances */
79 real orire_initf
; /* The scaling factor for initializing the time av. */
80 int norire_Dtav
; /* The number of matrix element in dtav (npair*5) */
81 real
*orire_Dtav
; /* The time averaged orientation tensors */
84 /* Struct used for checkpointing only.
85 * This struct would not be required with unlimited precision.
86 * But because of limited precision, the COM motion removal implementation
87 * can cause the kinetic energy in the MD loop to differ by a few bits from
88 * the kinetic energy one would determine from state.v.
90 typedef struct ekinstate_t
98 double *ekinscalef_nhc
;
99 double *ekinscaleh_nhc
;
105 typedef struct df_history_t
107 int nlambda
; /* total number of lambda states - for history*/
109 gmx_bool bEquil
; /* Have we reached equilibration */
110 int *n_at_lam
; /* number of points observed at each lambda */
111 real
*wl_histo
; /* histogram for WL flatness determination */
112 real wl_delta
; /* current wang-landau delta */
114 real
*sum_weights
; /* weights of the states */
115 real
*sum_dg
; /* free energies of the states -- not actually used for weighting, but informational */
116 real
*sum_minvar
; /* corrections to weights for minimum variance */
117 real
*sum_variance
; /* variances of the states */
119 real
**accum_p
; /* accumulated bennett weights for n+1 */
120 real
**accum_m
; /* accumulated bennett weights for n-1 */
121 real
**accum_p2
; /* accumulated squared bennett weights for n+1 */
122 real
**accum_m2
; /* accumulated squared bennett weights for n-1 */
124 real
**Tij
; /* transition matrix */
125 real
**Tij_empirical
; /* Empirical transition matrix */
129 typedef struct edsamstate_t
131 /* If one uses essential dynamics or flooding on a group of atoms from
132 * more than one molecule, we cannot make this group whole with
133 * do_pbc_first_mtop(). We assume that the ED group has the correct PBC
134 * representation at the beginning of the simulation and keep track
135 * of the shifts to always get it into that representation.
136 * For proper restarts from a checkpoint we store the positions of the
137 * reference group at the time of checkpoint writing */
138 gmx_bool bFromCpt
; /* Did we start from a checkpoint file? */
139 int nED
; /* No. of ED/Flooding data sets, if <1 no ED */
140 int *nref
; /* No. of atoms in i'th reference structure */
141 int *nav
; /* Same for average structure */
142 rvec
**old_sref
; /* Positions of the reference atoms
143 at the last time step (with correct PBC
145 rvec
**old_sref_p
; /* Pointer to these positions */
146 rvec
**old_sav
; /* Same for the average positions */
151 typedef struct swapstateIons_t
153 int nMolReq
[eCompNR
]; /* Requested # of molecules per compartment */
154 int *nMolReq_p
[eCompNR
]; /* Pointer to this data (for .cpt writing) */
155 int inflow_net
[eCompNR
]; /* Flux determined from the # of swaps */
156 int *inflow_net_p
[eCompNR
]; /* Pointer to this data */
157 int *nMolPast
[eCompNR
]; /* Array with nAverage entries for history */
158 int *nMolPast_p
[eCompNR
]; /* Pointer points to the first entry only */
160 /* Channel flux detection, this is counting only and has no influence on whether swaps */
161 /* are performed or not: */
162 int fluxfromAtoB
[eCompNR
]; /* Flux determined from the split cylinders */
163 int *fluxfromAtoB_p
[eCompNR
]; /* Pointer to this data */
164 int nMol
; /* # of molecules, size of the following arrays */
165 unsigned char *comp_from
; /* Ion came from which compartment? */
166 unsigned char *channel_label
; /* Through which channel did this ion pass? */
169 typedef struct swapstate_t
171 int eSwapCoords
; /* Swapping along x, y, or z-direction? */
172 int nIonTypes
; /* Number of ion types, this is the size of */
173 /* the following arrays */
174 int nAverage
; /* Use average over this many swap attempt */
175 /* steps when determining the ion counts */
176 int fluxleak
; /* Ions not going through any channel (bad!) */
177 int *fluxleak_p
; /* Pointer to this data */
179 /* To also make multimeric channel proteins whole, we save the last whole configuration */
180 /* of the channels in the checkpoint file. If we have no checkpoint file, we assume */
181 /* that the starting configuration has the correct PBC representation after making the */
182 /* individual molecules whole */
183 gmx_bool bFromCpt
; /* Did we start from a checkpoint file? */
184 int nat
[eChanNR
]; /* Size of xc_old_whole, i.e. the number of */
185 /* atoms in each channel */
186 rvec
*xc_old_whole
[eChanNR
]; /* Last known whole positions of the two */
187 /* channels (important for multimeric ch.!) */
188 rvec
**xc_old_whole_p
[eChanNR
]; /* Pointer to these positions */
189 swapstateIons_t
*ionType
;
194 typedef struct t_state
199 int nhchainlength
; /* number of nose-hoover chains */
200 int flags
; /* Flags telling which entries are present */
201 int fep_state
; /* indicates which of the alchemical states we are in */
202 real
*lambda
; /* lambda vector */
203 matrix box
; /* box vector coordinates */
204 matrix box_rel
; /* Relitaive box vectors to preserve shape */
205 matrix boxv
; /* box velocitites for Parrinello-Rahman pcoupl */
206 matrix pres_prev
; /* Pressure of the previous step for pcoupl */
207 matrix svir_prev
; /* Shake virial for previous step for pcoupl */
208 matrix fvir_prev
; /* Force virial of the previous step for pcoupl */
209 double *nosehoover_xi
; /* for Nose-Hoover tcoupl (ngtc) */
210 double *nosehoover_vxi
; /* for N-H tcoupl (ngtc) */
211 double *nhpres_xi
; /* for Nose-Hoover pcoupl for barostat */
212 double *nhpres_vxi
; /* for Nose-Hoover pcoupl for barostat */
213 double *therm_integral
; /* for N-H/V-rescale tcoupl (ngtc) */
214 real veta
; /* trotter based isotropic P-coupling */
215 real vol0
; /* initial volume,required for computing NPT conserverd quantity */
216 int nalloc
; /* Allocation size for x and v when !=NULL*/
217 rvec
*x
; /* the coordinates (natoms) */
218 rvec
*v
; /* the velocities (natoms) */
219 rvec
*cg_p
; /* p vector for conjugate gradient minimization */
221 history_t hist
; /* Time history for restraints */
223 ekinstate_t ekinstate
; /* The state of the kinetic energy data */
225 struct energyhistory_t
*enerhist
; /* Energy history for statistics */
226 swapstate_t swapstate
; /* Position swapping */
227 df_history_t dfhist
; /*Free energy history for free energy analysis */
228 edsamstate_t edsamstate
; /* Essential dynamics / flooding history */
230 int ddp_count
; /* The DD partitioning count for this state */
231 int ddp_count_cg_gl
; /* The DD part. count for index_gl */
232 int ncg_gl
; /* The number of local charge groups */
233 int *cg_gl
; /* The global cg number of the local cgs */
234 int cg_gl_nalloc
; /* Allocation size of cg_gl; */
237 typedef struct t_extmass
239 double *Qinv
; /* inverse mass of thermostat -- computed from inputs, but a good place to store */
240 double *QPinv
; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
241 double Winv
; /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
242 tensor Winvm
; /* inverse pressure mass tensor, computed */
256 void init_gtc_state(t_state
*state
, int ngtc
, int nnhpres
, int nhchainlength
);
258 void init_state(t_state
*state
, int natoms
, int ngtc
, int nnhpres
, int nhchainlength
, int nlambda
);
260 t_state
*serial_init_local_state(t_state
*state_global
);
262 void done_state(t_state
*state
);
264 void comp_state(const t_state
*st1
, const t_state
*st2
, gmx_bool bRMSD
, real ftol
, real abstol
);