Clean up some mathematical constants
[gromacs.git] / src / mdlib / mdebin_bar.h
blobbb4ffd21ed607e4fac913ce3cd36d513a371ba97
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
5 * G R O M A C S
7 * GROningen MAchine for Chemical Simulations
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.
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
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _mdebin_bar_h
37 #define _mdebin_bar_h
39 #ifdef __cplusplus
40 extern "C" {
41 #endif
43 /* The functions & data structures here describe writing
44 energy differences (or their histogram )for use with g_bar */
46 /* Data for one foreign lambda, or derivative. */
47 typedef struct
49 real *dh; /* the raw energy difference data -- actually, store more in here. */
50 float *dhf; /* raw difference data -- in floats, for storage. */
51 unsigned int ndh; /* number of data points */
52 unsigned int ndhmax; /* the maximum number of points */
54 int nhist; /* the number of histograms. There can either be
55 0 (for no histograms)
56 1 (for 'foreign lambda' histograms)
57 2 (for derivative histograms: there's
58 a 'forward' and 'backward' histogram
59 containing the minimum and maximum
60 values, respectively). */
61 int *bin[2]; /* the histogram(s) */
62 double dx; /* the histogram spacing in kJ/mol. This is the
63 same for the two histograms? */
64 unsigned int nbins; /* the number of bins in the histograms*/
65 gmx_large_int_t x0[2]; /* the starting point in units of spacing
66 of the histogram */
67 unsigned int maxbin[2]; /* highest bin number with data */
69 gmx_bool derivative; /* whether this delta_h contains derivatives */
70 double lambda; /* current lambda */
71 gmx_bool written; /* whether this data has already been written out */
73 double subblock_d[4]; /* data for an mdebin subblock for I/O. */
74 gmx_large_int_t subblock_l[4]; /* data for an mdebin subblock for I/O. */
75 int subblock_i[4]; /* data for an mdebin subblock for I/O. */
76 } t_mde_delta_h;
78 /* the type definition is in mdebin_bar.h */
79 struct t_mde_delta_h_coll
81 t_mde_delta_h *dh; /* the delta h data */
82 int ndh; /* the number of delta_h structures */
83 double start_time; /* start time of the current dh collection */
84 double delta_time; /* time difference between samples */
85 gmx_bool start_time_set; /* whether the start time has been set */
86 double start_lambda; /* starting lambda for continuous motion of state*/
87 double delta_lambda; /* delta lambda, for continuous motion of state */
88 double temperature; /* the temperature of the samples*/
89 double subblock_d[5]; /* data for writing an mdebin subblock for I/O */
94 /* initialize a collection of delta h histograms/sets
95 dhc = the collection
96 ir = the input record */
98 void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc,
99 const t_inputrec *ir);
101 /* add a bunch of samples to the delta_h collection
102 dhc = the collection
103 dhdl = the hamiltonian derivatives
104 U = the array with energies: from enerd->enerpart_lambda.
105 time = the current simulation time.
106 current_lambda = current lambda values : primarily useful for continuous processes
107 fep_state = current fep_state
110 /* add a bunch of samples - note fep_state is double to allow for better data storage */
111 void mde_delta_h_coll_add_dh(t_mde_delta_h_coll *dhc,
112 double fep_state,
113 double energy,
114 double pV,
115 int bExpanded,
116 int bPrintEnergy,
117 int bPressure,
118 int ndhdl,
119 int nlambda,
120 double *dhdl,
121 double *foreign_dU,
122 double time);
124 /* write the data associated with the du blocks collection as a collection
125 of mdebin blocks.
126 dhc = the collection
127 fr = the enxio frame
128 nblock = the current number of blocks */
129 void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
130 t_enxframe *fr, int nblock);
133 /* reset the collection of delta_h buffers for a new round of
134 data gathering */
135 void mde_delta_h_coll_reset(t_mde_delta_h_coll *dhc);
138 /* set the energyhistory variables to save state */
139 void mde_delta_h_coll_update_energyhistory(t_mde_delta_h_coll *dhc,
140 energyhistory_t *enerhist);
142 /* restore the variables from an energyhistory */
143 void mde_delta_h_coll_restore_energyhistory(t_mde_delta_h_coll *dhc,
144 energyhistory_t *enerhist);
147 #ifdef __cplusplus
149 #endif
151 #endif /* _mdebin_bar_h */