Further improve getDDGridSetup
[gromacs.git] / src / gromacs / mdlib / mdebin_bar.h
blobe7ce32b6f3f104a084e89c0f98d87f1683ecae72
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,2018,2019, 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/energyoutput.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;
47 struct t_enxframe;
49 /* Data for one foreign lambda, or derivative. */
50 typedef struct
52 real *dh; /* the raw energy data. */
53 float *dhf; /* raw difference data -- in floats, for storage. */
54 unsigned int ndh; /* number of data points */
55 unsigned int ndhmax; /* the maximum number of points */
57 int nhist; /* the number of histograms. There can either be
58 0 (for no histograms)
59 1 (for 'foreign lambda' histograms)
60 2 (for derivative histograms: there's
61 a 'forward' and 'backward' histogram
62 containing the minimum and maximum
63 values, respectively). */
64 int *bin[2]; /* the histogram(s) */
65 double dx; /* the histogram spacing in kJ/mol. This is the
66 same for the two histograms? */
67 unsigned int nbins; /* the number of bins in the histograms*/
68 int64_t x0[2]; /* the starting point in units of spacing
69 of the histogram */
70 unsigned int maxbin[2]; /* highest bin number with data */
72 int type; /* the block type according to dhbtDH, etc. */
73 int derivative; /* The derivative direction (as an index in the lambda
74 vector) if this delta_h contains derivatives */
75 double *lambda; /* lambda vector (or NULL if not applicable) */
76 int nlambda; /* length of the lambda vector */
77 gmx_bool written; /* whether this data has already been written out */
79 int64_t subblock_meta_l[5]; /* metadata for an mdebin subblock for
80 I/O: for histogram counts, etc.*/
81 double *subblock_meta_d; /* metadata subblock for I/O, used for
82 communicating doubles (i.e. the lambda
83 vector) */
84 int subblock_meta_i[4]; /* metadata subblock for I/O, used for
85 communicating ints (i.e. derivative indices,
86 etc.) */
87 } t_mde_delta_h;
89 /* the type definition is in mdebin_bar.h */
90 struct t_mde_delta_h_coll
92 t_mde_delta_h *dh; /* the delta h data */
93 int ndh; /* the number of delta_h structures */
95 int nlambda; /* number of bar dU delta_h structures */
96 t_mde_delta_h *dh_du; /* the delta h data (pointer into dh) */
98 int ndhdl; /* number of bar dU delta_h structures */
99 t_mde_delta_h *dh_dhdl; /* the dhdl data (pointer into dh) */
101 t_mde_delta_h *dh_energy; /* energy output block (pointer into dh) */
102 t_mde_delta_h *dh_pv; /* pV output block (pointer into dh) */
103 t_mde_delta_h *dh_expanded; /* expanded ensemble output block (pointer
104 into dh) */
106 double start_time; /* start time of the current dh collection */
107 double delta_time; /* time difference between samples */
108 gmx_bool start_time_set; /* whether the start time has been set */
109 double start_lambda; /* starting lambda for continuous motion of state*/
110 double delta_lambda; /* delta lambda, for continuous motion of state */
111 double temperature; /* the temperature of the samples*/
113 double *native_lambda_vec; /* The lambda vector describing the current
114 lambda state if it is set (NULL otherwise) */
115 int n_lambda_vec; /* the size of the native lambda vector */
116 int *native_lambda_components; /* the native lambda (and by extension,
117 foreign lambda) components in terms
118 of efptFEP, efptMASS, etc. */
119 int lambda_index; /* the lambda_fep_state */
121 double *subblock_d; /* for writing a metadata mdebin subblock for I/O */
122 int *subblock_i; /* for writing a metadata mdebin subblock for I/O */
124 double *lambda_vec_subblock; /* native lambda vector data subblock for
125 I/O */
126 int *lambda_index_subblock; /* lambda vector index data subblock for I/O */
131 /* initialize a collection of delta h histograms/sets
132 dhc = the collection
133 ir = the input record */
135 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc,
136 const t_inputrec *ir);
138 void done_mde_delta_h_coll(t_mde_delta_h_coll *dhc);
140 /* add a bunch of samples to the delta_h collection
141 dhc = the collection
142 dhdl = the hamiltonian derivatives
143 U = the array with energies: from enerd->enerpart_lambda.
144 time = the current simulation time.
145 current_lambda = current lambda values : primarily useful for continuous processes
146 fep_state = current fep_state
149 /* add a bunch of samples - note fep_state is double to allow for better data storage */
150 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
151 double fep_state,
152 double energy,
153 double pV,
154 double *dhdl,
155 double *foreign_dU,
156 double time);
158 /* write the data associated with the du blocks collection as a collection
159 of mdebin blocks.
160 dhc = the collection
161 fr = the enxio frame
162 nblock = the current number of blocks */
163 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
164 t_enxframe *fr, int nblock);
167 /* reset the collection of delta_h buffers for a new round of
168 data gathering */
169 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc);
172 /* set the energyhistory variables to save state */
173 void mde_delta_h_coll_update_energyhistory(const t_mde_delta_h_coll *dhc,
174 energyhistory_t *enerhist);
176 /* restore the variables from an energyhistory */
177 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
178 const delta_h_history_t *deltaH);
180 #endif /* _mdebin_bar_h */