3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GRoups of Organic Molecules in ACtion for Science
46 * The t_state struct should contain all the (possibly) non-static
47 * information required to define the state of the system.
48 * Currently the random seeds for SD and BD are missing.
51 /* for now, define the length of the NH chains here */
52 #define NHCHAINLENGTH 10
53 #define MAXLAMBDAS 1024
55 /* These enums are used in flags as (1<<est...).
56 * The order of these enums should not be changed,
57 * since that affects the checkpoint (.cpt) file format.
61 estBOX, estBOX_REL, estBOXV, estPRES_PREV, estNH_XI, estTC_INT,
62 estX, estV, estSDX, estCGP, estLD_RNG, estLD_RNGI,
63 estDISRE_INITF, estDISRE_RM3TAV,
64 estORIRE_INITF, estORIRE_DTAV,
65 estSVIR_PREV, estNH_VXI, estVETA, estVOL0, estNHPRES_XI, estNHPRES_VXI, estFVIR_PREV,
66 estFEPSTATE, estMC_RNG, estMC_RNGI,
70 #define EST_DISTR(e) (!(((e) >= estLAMBDA && (e) <= estTC_INT) || ((e) >= estSVIR_PREV && (e) <= estMC_RNGI)))
72 /* The names of the state entries, defined in src/gmxlib/checkpoint.c */
73 extern const char *est_names[estNR];
77 real disre_initf; /* The scaling factor for initializing the time av. */
78 int ndisrepairs; /* The number of distance restraints */
79 real *disre_rm3tav; /* The r^-3 time averaged pair distances */
80 real orire_initf; /* The scaling factor for initializing the time av. */
81 int norire_Dtav; /* The number of matrix element in dtav (npair*5) */
82 real *orire_Dtav; /* The time averaged orientation tensors */
85 /* Struct used for checkpointing only.
86 * This struct would not be required with unlimited precision.
87 * But because of limited precision, the COM motion removal implementation
88 * can cause the kinetic energy in the MD loop to differ by a few bits from
89 * the kinetic energy one would determine from state.v.
99 double *ekinscalef_nhc;
100 double *ekinscaleh_nhc;
106 /* energy history for delta_h histograms */
109 int nndh; /* the number of energy difference lists */
110 int *ndh; /* the number in each energy difference list */
111 real **dh; /* the energy difference lists */
113 double start_time; /* the start time of these energy diff blocks */
114 double start_lambda; /* lambda at start time */
116 gmx_bool start_lambda_set; /* whether the lambda value is set. Here
117 For backward-compatibility. */
122 int nlambda; /* total number of lambda states - for history*/
124 gmx_bool bEquil; /* Have we reached equilibration */
125 int *n_at_lam; /* number of points observed at each lambda */
126 real *wl_histo; /* histogram for WL flatness determination */
127 real wl_delta; /* current wang-landau delta */
129 real *sum_weights; /* weights of the states */
130 real *sum_dg; /* free energies of the states -- not actually used for weighting, but informational */
131 real *sum_minvar; /* corrections to weights for minimum variance */
132 real *sum_variance; /* variances of the states */
134 real **accum_p; /* accumulated bennett weights for n+1 */
135 real **accum_m; /* accumulated bennett weights for n-1 */
136 real **accum_p2; /* accumulated squared bennett weights for n+1 */
137 real **accum_m2; /* accumulated squared bennett weights for n-1 */
139 real **Tij; /* transition matrix */
140 real **Tij_empirical; /* Empirical transition matrix */
146 gmx_large_int_t nsteps; /* The number of steps in the history */
147 gmx_large_int_t nsum; /* The nr. of steps in the ener_ave and ener_sum */
148 double * ener_ave; /* Energy term history sum to get fluctuations */
149 double * ener_sum; /* Energy term history sum to get fluctuations */
150 int nener; /* Number of energy terms in two previous arrays */
151 gmx_large_int_t nsteps_sim; /* The number of steps in ener_sum_sim */
152 gmx_large_int_t nsum_sim; /* The number of frames in ener_sum_sim */
153 double * ener_sum_sim; /* Energy term history sum of the whole sim */
155 delta_h_history_t *dht; /* The BAR energy differences */
161 /* If one uses essential dynamics or flooding on a group of atoms from
162 * more than one molecule, we cannot make this group whole with
163 * do_pbc_first_mtop(). We assume that the ED group has the correct PBC
164 * representation at the beginning of the simulation and keep track
165 * of the shifts to always get it into that representation.
166 * For proper restarts from a checkpoint we store the positions of the
167 * reference group at the time of checkpoint writing */
168 gmx_bool bFromCpt; /* Did we start from a checkpoint file? */
169 int nED; /* No. of ED/Flooding data sets, if <1 no ED */
170 int *nref; /* No. of atoms in i'th reference structure */
171 int *nav; /* Same for average structure */
172 rvec **old_sref; /* Positions of the reference atoms
173 at the last time step (with correct PBC
175 rvec **old_sref_p; /* Pointer to these positions */
176 rvec **old_sav; /* Same for the average positions */
186 int nhchainlength; /* number of nose-hoover chains */
189 int flags; /* Flags telling which entries are present */
190 int fep_state; /* indicates which of the alchemical states we are in */
191 real *lambda; /* lambda vector */
192 matrix box; /* box vector coordinates */
193 matrix box_rel; /* Relitaive box vectors to preserve shape */
194 matrix boxv; /* box velocitites for Parrinello-Rahman pcoupl */
195 matrix pres_prev; /* Pressure of the previous step for pcoupl */
196 matrix svir_prev; /* Shake virial for previous step for pcoupl */
197 matrix fvir_prev; /* Force virial of the previous step for pcoupl */
198 double *nosehoover_xi; /* for Nose-Hoover tcoupl (ngtc) */
199 double *nosehoover_vxi; /* for N-H tcoupl (ngtc) */
200 double *nhpres_xi; /* for Nose-Hoover pcoupl for barostat */
201 double *nhpres_vxi; /* for Nose-Hoover pcoupl for barostat */
202 double *therm_integral; /* for N-H/V-rescale tcoupl (ngtc) */
203 real veta; /* trotter based isotropic P-coupling */
204 real vol0; /* initial volume,required for computing NPT conserverd quantity */
205 int nalloc; /* Allocation size for x, v and sd_x when !=NULL*/
206 rvec *x; /* the coordinates (natoms) */
207 rvec *v; /* the velocities (natoms) */
208 rvec *sd_X; /* random part of the x update for stoch. dyn. */
209 rvec *cg_p; /* p vector for conjugate gradient minimization */
211 unsigned int *ld_rng; /* RNG random state */
212 int *ld_rngi; /* RNG index */
214 int nmcrng; /* number of RNG states */
215 unsigned int *mc_rng; /* lambda MC RNG random state */
216 int *mc_rngi; /* lambda MC RNG index */
218 history_t hist; /* Time history for restraints */
220 ekinstate_t ekinstate; /* The state of the kinetic energy data */
222 energyhistory_t enerhist; /* Energy history for statistics */
223 df_history_t dfhist; /*Free energy history for free energy analysis */
224 edsamstate_t edsamstate; /* Essential dynamics / flooding history */
226 int ddp_count; /* The DD partitioning count for this state */
227 int ddp_count_cg_gl; /* The DD part. count for index_gl */
228 int ncg_gl; /* The number of local charge groups */
229 int *cg_gl; /* The global cg number of the local cgs */
230 int cg_gl_nalloc; /* Allocation size of cg_gl; */
235 double *Qinv; /* inverse mass of thermostat -- computed from inputs, but a good place to store */
236 double *QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
237 double Winv; /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
238 tensor Winvm; /* inverse pressure mass tensor, computed */
257 #endif /* _state_h_ */