Fix remaining copyright headers
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / state.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef _state_h_
38 #define _state_h_
39
40
41 #include "simple.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47 /*
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.
51  */
52
53 /* for now, define the length of the NH chains here */
54 #define NHCHAINLENGTH 10
55 #define MAXLAMBDAS 1024
56
57 /* These enums are used in flags as (1<<est...).
58  * The order of these enums should not be changed,
59  * since that affects the checkpoint (.cpt) file format.
60  */
61 enum {
62     estLAMBDA,
63     estBOX, estBOX_REL, estBOXV, estPRES_PREV, estNH_XI,  estTC_INT,
64     estX,   estV,       estSDX,  estCGP,       estLD_RNG, estLD_RNGI,
65     estDISRE_INITF, estDISRE_RM3TAV,
66     estORIRE_INITF, estORIRE_DTAV,
67     estSVIR_PREV, estNH_VXI, estVETA, estVOL0, estNHPRES_XI, estNHPRES_VXI, estFVIR_PREV,
68     estFEPSTATE, estMC_RNG, estMC_RNGI,
69     estNR
70 };
71
72 #define EST_DISTR(e) (!(((e) >= estLAMBDA && (e) <= estTC_INT) || ((e) >= estSVIR_PREV && (e) <= estMC_RNGI)))
73
74 /* The names of the state entries, defined in src/gmxlib/checkpoint.c */
75 extern const char *est_names[estNR];
76
77 typedef struct
78 {
79     real  disre_initf;  /* The scaling factor for initializing the time av. */
80     int   ndisrepairs;  /* The number of distance restraints                */
81     real *disre_rm3tav; /* The r^-3 time averaged pair distances            */
82     real  orire_initf;  /* The scaling factor for initializing the time av. */
83     int   norire_Dtav;  /* The number of matrix element in dtav (npair*5)   */
84     real *orire_Dtav;   /* The time averaged orientation tensors            */
85 } history_t;
86
87 /* Struct used for checkpointing only.
88  * This struct would not be required with unlimited precision.
89  * But because of limited precision, the COM motion removal implementation
90  * can cause the kinetic energy in the MD loop to differ by a few bits from
91  * the kinetic energy one would determine from state.v.
92  */
93 typedef struct
94 {
95     gmx_bool     bUpToDate;
96     int          ekin_n;
97     tensor      *ekinh;
98     tensor      *ekinf;
99     tensor      *ekinh_old;
100     tensor       ekin_total;
101     double      *ekinscalef_nhc;
102     double      *ekinscaleh_nhc;
103     double      *vscale_nhc;
104     real         dekindl;
105     real         mvcos;
106 } ekinstate_t;
107
108 /* energy history for delta_h histograms */
109 typedef struct
110 {
111     int      nndh;             /* the number of energy difference lists */
112     int     *ndh;              /* the number in each energy difference list */
113     real   **dh;               /* the energy difference lists */
114
115     double   start_time;       /* the start time of these energy diff blocks */
116     double   start_lambda;     /* lambda at start time */
117
118     gmx_bool start_lambda_set; /* whether the lambda value is set. Here
119                                   For backward-compatibility. */
120 } delta_h_history_t;
121
122 typedef struct
123 {
124     int      nlambda;        /* total number of lambda states - for history*/
125
126     gmx_bool bEquil;         /* Have we reached equilibration */
127     int     *n_at_lam;       /* number of points observed at each lambda */
128     real    *wl_histo;       /* histogram for WL flatness determination */
129     real     wl_delta;       /* current wang-landau delta */
130
131     real    *sum_weights;    /* weights of the states */
132     real    *sum_dg;         /* free energies of the states -- not actually used for weighting, but informational */
133     real    *sum_minvar;     /* corrections to weights for minimum variance */
134     real    *sum_variance;   /* variances of the states */
135
136     real   **accum_p;        /* accumulated bennett weights for n+1 */
137     real   **accum_m;        /* accumulated bennett weights for n-1 */
138     real   **accum_p2;       /* accumulated squared bennett weights for n+1 */
139     real   **accum_m2;       /* accumulated squared bennett weights for n-1 */
140
141     real   **Tij;            /* transition matrix */
142     real   **Tij_empirical;  /* Empirical transition matrix */
143
144 } df_history_t;
145
146 typedef struct
147 {
148     gmx_int64_t        nsteps;       /* The number of steps in the history            */
149     gmx_int64_t        nsum;         /* The nr. of steps in the ener_ave and ener_sum */
150     double         *   ener_ave;     /* Energy term history sum to get fluctuations   */
151     double         *   ener_sum;     /* Energy term history sum to get fluctuations   */
152     int                nener;        /* Number of energy terms in two previous arrays */
153     gmx_int64_t        nsteps_sim;   /* The number of steps in ener_sum_sim      */
154     gmx_int64_t        nsum_sim;     /* The number of frames in ener_sum_sim     */
155     double         *   ener_sum_sim; /* Energy term history sum of the whole sim      */
156
157     delta_h_history_t *dht;          /* The BAR energy differences */
158 }
159 energyhistory_t;
160
161 typedef struct
162 {
163     /* If one uses essential dynamics or flooding on a group of atoms from
164      * more than one molecule, we cannot make this group whole with
165      * do_pbc_first_mtop(). We assume that the ED group has the correct PBC
166      * representation at the beginning of the simulation and keep track
167      * of the shifts to always get it into that representation.
168      * For proper restarts from a checkpoint we store the positions of the
169      * reference group at the time of checkpoint writing */
170     gmx_bool      bFromCpt;     /* Did we start from a checkpoint file?       */
171     int           nED;          /* No. of ED/Flooding data sets, if <1 no ED  */
172     int          *nref;         /* No. of atoms in i'th reference structure   */
173     int          *nav;          /* Same for average structure                 */
174     rvec        **old_sref;     /* Positions of the reference atoms
175                                    at the last time step (with correct PBC
176                                    representation)                            */
177     rvec        **old_sref_p;   /* Pointer to these positions                 */
178     rvec        **old_sav;      /* Same for the average positions             */
179     rvec        **old_sav_p;
180 }
181 edsamstate_t;
182
183 typedef struct
184 {
185     int              natoms;
186     int              ngtc;
187     int              nnhpres;
188     int              nhchainlength; /* number of nose-hoover chains               */
189     int              nrng;
190     int              nrngi;
191     int              flags;           /* Flags telling which entries are present      */
192     int              fep_state;       /* indicates which of the alchemical states we are in                 */
193     real            *lambda;          /* lambda vector                               */
194     matrix           box;             /* box vector coordinates                         */
195     matrix           box_rel;         /* Relitaive box vectors to preserve shape        */
196     matrix           boxv;            /* box velocitites for Parrinello-Rahman pcoupl */
197     matrix           pres_prev;       /* Pressure of the previous step for pcoupl  */
198     matrix           svir_prev;       /* Shake virial for previous step for pcoupl */
199     matrix           fvir_prev;       /* Force virial of the previous step for pcoupl  */
200     double          *nosehoover_xi;   /* for Nose-Hoover tcoupl (ngtc)       */
201     double          *nosehoover_vxi;  /* for N-H tcoupl (ngtc)               */
202     double          *nhpres_xi;       /* for Nose-Hoover pcoupl for barostat     */
203     double          *nhpres_vxi;      /* for Nose-Hoover pcoupl for barostat     */
204     double          *therm_integral;  /* for N-H/V-rescale tcoupl (ngtc)     */
205     real             veta;            /* trotter based isotropic P-coupling             */
206     real             vol0;            /* initial volume,required for computing NPT conserverd quantity */
207     int              nalloc;          /* Allocation size for x, v and sd_x when !=NULL*/
208     rvec            *x;               /* the coordinates (natoms)                     */
209     rvec            *v;               /* the velocities (natoms)                      */
210     rvec            *sd_X;            /* random part of the x update for stoch. dyn.  */
211     rvec            *cg_p;            /* p vector for conjugate gradient minimization */
212
213     unsigned int    *ld_rng;          /* RNG random state                           */
214     int             *ld_rngi;         /* RNG index                                  */
215
216     int              nmcrng;          /* number of RNG states                       */
217     unsigned int    *mc_rng;          /* lambda MC RNG random state                 */
218     int             *mc_rngi;         /* lambda MC RNG index                        */
219
220     history_t        hist;            /* Time history for restraints                  */
221
222     ekinstate_t      ekinstate;       /* The state of the kinetic energy data      */
223
224     energyhistory_t  enerhist;        /* Energy history for statistics           */
225     df_history_t     dfhist;          /*Free energy history for free energy analysis  */
226     edsamstate_t     edsamstate;      /* Essential dynamics / flooding history */
227
228     int              ddp_count;       /* The DD partitioning count for this state  */
229     int              ddp_count_cg_gl; /* The DD part. count for index_gl     */
230     int              ncg_gl;          /* The number of local charge groups            */
231     int             *cg_gl;           /* The global cg number of the local cgs        */
232     int              cg_gl_nalloc;    /* Allocation size of cg_gl;              */
233 } t_state;
234
235 typedef struct
236 {
237     double *Qinv;  /* inverse mass of thermostat -- computed from inputs, but a good place to store */
238     double *QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
239     double  Winv;  /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
240     tensor  Winvm; /* inverse pressure mass tensor, computed       */
241 } t_extmass;
242
243
244 typedef struct
245 {
246     real    veta;
247     double  rscale;
248     double  vscale;
249     double  rvscale;
250     double  alpha;
251     double *vscale_nhc;
252 } t_vetavars;
253
254 #ifdef __cplusplus
255 }
256 #endif
257
258
259 #endif /* _state_h_ */