Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / mdebin_bar.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gromacs Runs On Most of All Computer Systems
34  */
35
36 #ifndef _mdebin_bar_h
37 #define _mdebin_bar_h
38
39 #ifdef __cplusplus
40 extern "C" {
41 #endif
42
43
44 /* The functions & data structures here describe writing
45    energy differences (or their histogram )for use with g_bar */
46
47 /* Data for one foreign lambda, or derivative. */
48 typedef struct
49 {
50     real        *dh;                    /* the raw energy data. */
51     float       *dhf;                   /* raw difference data -- in floats, for storage. */
52     unsigned int ndh;                   /* number of data points */
53     unsigned int ndhmax;                /* the maximum number of points */
54
55     int          nhist;                 /* the number of histograms. There can either be
56                                            0 (for no histograms)
57                                            1 (for 'foreign lambda' histograms)
58                                            2 (for derivative histograms: there's
59                                               a 'forward' and 'backward' histogram
60                                               containing the minimum and maximum
61                                               values, respectively). */
62     int            *bin[2];             /* the histogram(s) */
63     double          dx;                 /* the histogram spacing in kJ/mol. This is the
64                                            same for the two histograms? */
65     unsigned int    nbins;              /* the number of bins in the histograms*/
66     gmx_large_int_t x0[2];              /* the starting point in units of spacing
67                                            of the histogram */
68     unsigned int    maxbin[2];          /* highest bin number with data */
69
70     int             type;               /* the block type according to dhbtDH, etc. */
71     int             derivative;         /* The derivative direction (as an index in the lambda
72                                            vector) if this delta_h contains derivatives */
73     double         *lambda;             /* lambda vector (or NULL if not applicable) */
74     int             nlambda;            /* length of the lambda vector */
75     gmx_bool        written;            /* whether this data has already been written out */
76
77     gmx_large_int_t subblock_meta_l[5]; /* metadata for an mdebin subblock for
78                                            I/O: for histogram counts, etc.*/
79     double         *subblock_meta_d;    /* metadata subblock for I/O, used for
80                                            communicating doubles (i.e. the lambda
81                                            vector) */
82     int subblock_meta_i[4];             /* metadata subblock for I/O, used for
83                                            communicating ints (i.e. derivative indices,
84                                            etc.) */
85 } t_mde_delta_h;
86
87 /* the type definition is in mdebin_bar.h */
88 struct t_mde_delta_h_coll
89 {
90     t_mde_delta_h *dh;                 /* the delta h data */
91     int            ndh;                /* the number of delta_h structures */
92
93     int            nlambda;            /* number of bar dU delta_h structures */
94     t_mde_delta_h *dh_du;              /* the delta h data (pointer into dh) */
95
96     int            ndhdl;              /* number of bar dU delta_h structures */
97     t_mde_delta_h *dh_dhdl;            /* the dhdl data (pointer into dh) */
98
99     t_mde_delta_h *dh_energy;          /* energy output block (pointer into dh) */
100     t_mde_delta_h *dh_pv;              /* pV output block (pointer into dh) */
101     t_mde_delta_h *dh_expanded;        /* expanded ensemble output block (pointer
102                                           into dh) */
103
104     double   start_time;               /* start time of the current dh collection */
105     double   delta_time;               /* time difference between samples */
106     gmx_bool start_time_set;           /* whether the start time has been set */
107     double   start_lambda;             /* starting lambda for continuous motion of state*/
108     double   delta_lambda;             /* delta lambda, for continuous motion of state */
109     double   temperature;              /* the temperature of the samples*/
110
111     double  *native_lambda_vec;        /* The lambda vector describing the current
112                                           lambda state if it is set (NULL otherwise) */
113     int      n_lambda_vec;             /* the size of the native lambda vector */
114     int     *native_lambda_components; /* the native lambda (and by extension,
115                                           foreign lambda) components in terms
116                                           of efptFEP, efptMASS, etc. */
117     int     lambda_index;              /* the lambda_fep_state */
118
119     double *subblock_d;                /* for writing a metadata mdebin subblock for I/O */
120     int    *subblock_i;                /* for writing a metadata mdebin subblock for I/O */
121
122     double *lambda_vec_subblock;       /* native lambda vector data subblock for
123                                           I/O */
124     int    *lambda_index_subblock;     /* lambda vector index data subblock for I/O */
125 };
126
127
128
129 /* initialize a collection of delta h histograms/sets
130     dhc = the collection
131     ir = the input record */
132
133 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc,
134                            const t_inputrec   *ir);
135
136 /* add a bunch of samples to the delta_h collection
137     dhc = the collection
138     dhdl = the hamiltonian derivatives
139     U = the array with energies: from enerd->enerpart_lambda.
140     time = the current simulation time.
141     current_lambda = current lambda values : primarily useful for continuous processes
142     fep_state = current fep_state
143  */
144
145 /* add a bunch of samples - note fep_state is double to allow for better data storage */
146 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
147                              double              fep_state,
148                              double              energy,
149                              double              pV,
150                              double             *dhdl,
151                              double             *foreign_dU,
152                              double              time);
153
154 /* write the data associated with the du blocks collection as a collection
155     of mdebin blocks.
156     dhc = the collection
157     fr = the enxio frame
158     nblock = the current number of blocks */
159 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
160                                    t_enxframe *fr, int nblock);
161
162
163 /* reset the collection of delta_h buffers for a new round of
164    data gathering */
165 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc);
166
167
168 /* set the energyhistory variables to save state */
169 void mde_delta_h_coll_update_energyhistory(t_mde_delta_h_coll *dhc,
170                                            energyhistory_t    *enerhist);
171
172 /* restore the variables from an energyhistory */
173 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
174                                             energyhistory_t    *enerhist);
175
176
177 #ifdef __cplusplus
178 }
179 #endif
180
181 #endif  /* _mdebin_bar_h */