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
54 /* These enums are used in flags as (1<<est...).
55 * The order of these enums should not be changed,
56 * since that affects the checkpoint (.cpt) file format.
59 estBOX, estBOX_REL, estBOXV, estPRES_PREV, estNH_XI, estTC_INT,
60 estX, estV, estSDX, estCGP, estLD_RNG, estLD_RNGI,
61 estDISRE_INITF, estDISRE_RM3TAV,
62 estORIRE_INITF, estORIRE_DTAV,
63 estSVIR_PREV, estNH_VXI, estVETA, estVOL0, estNHPRES_XI, estNHPRES_VXI,estFVIR_PREV,
66 #define EST_DISTR(e) (!(((e) >= estLAMBDA && (e) <= estTC_INT) || ((e) >= estSVIR_PREV && (e) <= estFVIR_PREV)))
68 /* The names of the state entries, defined in src/gmxib/checkpoint.c */
69 extern const char *est_names[estNR];
73 real disre_initf; /* The scaling factor for initializing the time av. */
74 int ndisrepairs; /* The number of distance restraints */
75 real *disre_rm3tav; /* The r^-3 time averaged pair distances */
76 real orire_initf; /* The scaling factor for initializing the time av. */
77 int norire_Dtav; /* The number of matrix element in dtav (npair*5) */
78 real *orire_Dtav; /* The time averaged orientation tensors */
81 /* Struct used for checkpointing only.
82 * This struct would not be required with unlimited precision.
83 * But because of limited precision, the COM motion removal implementation
84 * can cause the kinetic energy in the MD loop to differ by a few bits from
85 * the kinetic energy one would determine from state.v.
95 double *ekinscalef_nhc;
96 double *ekinscaleh_nhc;
102 /* energy history for delta_h histograms */
105 int nndh; /* the number of energy difference lists */
106 int *ndh; /* the number in each energy difference list */
107 real **dh; /* the energy difference lists */
109 double start_time; /* the start time of these energy diff blocks */
110 double start_lambda; /* lambda at start time */
112 gmx_bool start_lambda_set; /* whether the lambda value is set. Here
113 For backward-compatibility. */
118 gmx_large_int_t nsteps; /* The number of steps in the history */
119 gmx_large_int_t nsum; /* The nr. of steps in the ener_ave and ener_sum */
120 double * ener_ave; /* Energy term history sum to get fluctuations */
121 double * ener_sum; /* Energy term history sum to get fluctuations */
122 int nener; /* Number of energy terms in two previous arrays */
123 gmx_large_int_t nsteps_sim; /* The number of steps in ener_sum_sim */
124 gmx_large_int_t nsum_sim; /* The number of frames in ener_sum_sim */
125 double * ener_sum_sim; /* Energy term history sum of the whole sim */
127 delta_h_history_t *dht; /* The BAR energy differences */
136 int nhchainlength; /* length of each nose-hoover chain */
139 int flags; /* Flags telling which entries are present */
140 real lambda; /* the free energy switching parameter */
141 matrix box; /* box vector coordinates */
142 matrix box_rel; /* Relitaive box vectors to preserve shape */
143 matrix boxv; /* box velocitites for Parrinello-Rahman pcoupl */
144 matrix pres_prev; /* Pressure of the previous step for pcoupl */
145 matrix svir_prev; /* Shake virial for previous step for pcoupl */
146 matrix fvir_prev; /* Force virial of the previous step for pcoupl */
147 double *nosehoover_xi; /* for Nose-Hoover tcoupl (ngtc) */
148 double *nosehoover_vxi; /* for N-H tcoupl (ngtc) */
149 double *nhpres_xi; /* for Nose-Hoover pcoupl for barostat */
150 double *nhpres_vxi; /* for Nose-Hoover pcoupl for barostat */
151 double *therm_integral; /* for N-H/V-rescale tcoupl (ngtc) */
152 real veta; /* trotter based isotropic P-coupling */
153 real vol0; /* initial volume,required for computing NPT conserverd quantity */
154 int nalloc; /* Allocation size for x, v and sd_x when !=NULL*/
155 rvec *x; /* the coordinates (natoms) */
156 rvec *v; /* the velocities (natoms) */
157 rvec *sd_X; /* random part of the x update for stoch. dyn. */
158 rvec *cg_p; /* p vector for conjugate gradient minimization */
160 unsigned int *ld_rng; /* RNG random state */
161 int *ld_rngi; /* RNG index */
163 history_t hist; /* Time history for restraints */
165 ekinstate_t ekinstate; /* The state of the kinetic energy data */
167 energyhistory_t enerhist; /* Energy history for statistics */
169 int ddp_count; /* The DD partitioning count for this state */
170 int ddp_count_cg_gl; /* The DD part. count for index_gl */
171 int ncg_gl; /* The number of local charge groups */
172 int *cg_gl; /* The global cg number of the local cgs */
173 int cg_gl_nalloc; /* Allocation size of cg_gl; */
178 double *Qinv; /* inverse mass of thermostat -- computed from inputs, but a good place to store */
179 double *QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
180 double Winv; /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
181 tensor Winvm; /* inverse pressure mass tensor, computed */
200 #endif /* _state_h_ */