Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / mdebin_bar.h
bloba506d7e0d23154511e2ef30c3f8d9d808ae6e66f
1 /*
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,2017, 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.
38 #ifndef _mdebin_bar_h
39 #define _mdebin_bar_h
41 #include "gromacs/mdlib/mdebin.h"
43 /* The functions & data structures here describe writing
44 energy differences (or their histogram )for use with g_bar */
46 class delta_h_history_t;
48 /* Data for one foreign lambda, or derivative. */
49 typedef struct
51 real *dh; /* the raw energy data. */
52 float *dhf; /* raw difference data -- in floats, for storage. */
53 unsigned int ndh; /* number of data points */
54 unsigned int ndhmax; /* the maximum number of points */
56 int nhist; /* the number of histograms. There can either be
57 0 (for no histograms)
58 1 (for 'foreign lambda' histograms)
59 2 (for derivative histograms: there's
60 a 'forward' and 'backward' histogram
61 containing the minimum and maximum
62 values, respectively). */
63 int *bin[2]; /* the histogram(s) */
64 double dx; /* the histogram spacing in kJ/mol. This is the
65 same for the two histograms? */
66 unsigned int nbins; /* the number of bins in the histograms*/
67 gmx_int64_t x0[2]; /* the starting point in units of spacing
68 of the histogram */
69 unsigned int maxbin[2]; /* highest bin number with data */
71 int type; /* the block type according to dhbtDH, etc. */
72 int derivative; /* The derivative direction (as an index in the lambda
73 vector) if this delta_h contains derivatives */
74 double *lambda; /* lambda vector (or NULL if not applicable) */
75 int nlambda; /* length of the lambda vector */
76 gmx_bool written; /* whether this data has already been written out */
78 gmx_int64_t subblock_meta_l[5]; /* metadata for an mdebin subblock for
79 I/O: for histogram counts, etc.*/
80 double *subblock_meta_d; /* metadata subblock for I/O, used for
81 communicating doubles (i.e. the lambda
82 vector) */
83 int subblock_meta_i[4]; /* metadata subblock for I/O, used for
84 communicating ints (i.e. derivative indices,
85 etc.) */
86 } t_mde_delta_h;
88 /* the type definition is in mdebin_bar.h */
89 struct t_mde_delta_h_coll
91 t_mde_delta_h *dh; /* the delta h data */
92 int ndh; /* the number of delta_h structures */
94 int nlambda; /* number of bar dU delta_h structures */
95 t_mde_delta_h *dh_du; /* the delta h data (pointer into dh) */
97 int ndhdl; /* number of bar dU delta_h structures */
98 t_mde_delta_h *dh_dhdl; /* the dhdl data (pointer into dh) */
100 t_mde_delta_h *dh_energy; /* energy output block (pointer into dh) */
101 t_mde_delta_h *dh_pv; /* pV output block (pointer into dh) */
102 t_mde_delta_h *dh_expanded; /* expanded ensemble output block (pointer
103 into dh) */
105 double start_time; /* start time of the current dh collection */
106 double delta_time; /* time difference between samples */
107 gmx_bool start_time_set; /* whether the start time has been set */
108 double start_lambda; /* starting lambda for continuous motion of state*/
109 double delta_lambda; /* delta lambda, for continuous motion of state */
110 double temperature; /* the temperature of the samples*/
112 double *native_lambda_vec; /* The lambda vector describing the current
113 lambda state if it is set (NULL otherwise) */
114 int n_lambda_vec; /* the size of the native lambda vector */
115 int *native_lambda_components; /* the native lambda (and by extension,
116 foreign lambda) components in terms
117 of efptFEP, efptMASS, etc. */
118 int lambda_index; /* the lambda_fep_state */
120 double *subblock_d; /* for writing a metadata mdebin subblock for I/O */
121 int *subblock_i; /* for writing a metadata mdebin subblock for I/O */
123 double *lambda_vec_subblock; /* native lambda vector data subblock for
124 I/O */
125 int *lambda_index_subblock; /* lambda vector index data subblock for I/O */
130 /* initialize a collection of delta h histograms/sets
131 dhc = the collection
132 ir = the input record */
134 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc,
135 const t_inputrec *ir);
137 /* add a bunch of samples to the delta_h collection
138 dhc = the collection
139 dhdl = the hamiltonian derivatives
140 U = the array with energies: from enerd->enerpart_lambda.
141 time = the current simulation time.
142 current_lambda = current lambda values : primarily useful for continuous processes
143 fep_state = current fep_state
146 /* add a bunch of samples - note fep_state is double to allow for better data storage */
147 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
148 double fep_state,
149 double energy,
150 double pV,
151 double *dhdl,
152 double *foreign_dU,
153 double time);
155 /* write the data associated with the du blocks collection as a collection
156 of mdebin blocks.
157 dhc = the collection
158 fr = the enxio frame
159 nblock = the current number of blocks */
160 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
161 t_enxframe *fr, int nblock);
164 /* reset the collection of delta_h buffers for a new round of
165 data gathering */
166 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc);
169 /* set the energyhistory variables to save state */
170 void mde_delta_h_coll_update_energyhistory(const t_mde_delta_h_coll *dhc,
171 energyhistory_t *enerhist);
173 /* restore the variables from an energyhistory */
174 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
175 const delta_h_history_t *deltaH);
177 #endif /* _mdebin_bar_h */