2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
5 * Copyright (c) 2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
40 * Declares the CorrelationTensor class for correlation tensor statistics.
42 * \author Viveca Lindahl
43 * \author Berk Hess <hess@kth.se>
47 #ifndef GMX_AWH_CORRELATIONTENSOR_H
48 #define GMX_AWH_CORRELATIONTENSOR_H
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/gmxassert.h"
61 struct CorrelationBlockDataHistory;
63 /*! \internal \brief Correlation block averaging weight-only data.
65 class CorrelationBlockData
68 /*! \internal \brief Correlation block averaging data.
72 /*! \brief Constructor.
74 CoordData() : blockSumWeightX(0), sumOverBlocksBlockWeightBlockWeightX(0) {}
76 double blockSumWeightX; /**< Weighted sum of x for current block. */
77 double sumOverBlocksBlockWeightBlockWeightX; /**< Sum over all blocks in the simulation of block weight times sum_wx. */
80 /*! \brief Constructor.
82 * \param[in] numDim The dimensionality.
83 * \param[in] blockLengthInit The initial block length.
85 CorrelationBlockData(int numDim, double blockLengthInit) :
87 blockSumSquareWeight_(0),
88 sumOverBlocksSquareBlockWeight_(0),
89 sumOverBlocksBlockSquareWeight_(0),
90 blockLength_(blockLengthInit),
91 previousBlockIndex_(-1),
93 correlationIntegral_(numDim * (numDim + 1) / 2)
97 /*! \brief Restore the state from history.
99 * \param[in] blockHistory The block data history containing the weight sums.
100 * \param[in] coordData The coordinate data.
101 * \param[in] correlationIntegral The correlation integral for all tensor elements.
103 void restoreFromHistory(const CorrelationBlockDataHistory& blockHistory,
104 const std::vector<CoordData>& coordData,
105 const std::vector<double>& correlationIntegral);
107 /*! \brief Adds a weighted data vector to one point in the correlation grid.
109 * \param[in] weight The weight of the data.
110 * \param[in] data One data point for each grid dimension.
112 void addData(double weight, gmx::ArrayRef<const double> data)
114 GMX_ASSERT(data.size() == coordData_.size(),
115 "Size of data should match the size of coordData");
117 blockSumWeight_ += weight;
118 blockSumSquareWeight_ += weight * weight;
120 for (size_t d = 0; d < coordData_.size(); d++)
122 coordData_[d].blockSumWeightX += weight * data[d];
126 /*! \brief Adds a filled data block to correlation time integral.
128 void addBlockToCorrelationIntegral();
130 /*! \brief Returns the sum weights for current block. */
131 double blockSumWeight() const { return blockSumWeight_; }
133 /*! \brief Returns the sum weights^2 for current block. */
134 double blockSumSquareWeight() const { return blockSumSquareWeight_; }
136 /*! \brief Returns the sum over blocks of block weight^2. */
137 double sumOverBlocksSquareBlockWeight() const { return sumOverBlocksSquareBlockWeight_; }
139 /*! \brief Returns the sum over blocks of weight^2. */
140 double sumOverBlocksBlockSquareWeight() const { return sumOverBlocksBlockSquareWeight_; }
142 /*! \brief Returns the length of each block used for block averaging. */
143 double blockLength() const { return blockLength_; }
145 /*! \brief Double the length of each block used for block averaging. */
146 void doubleBlockLength() { blockLength_ *= 2; }
148 /*! \brief Return the last block index data was added to (needed only for block length in terms of time). */
149 int previousBlockIndex() const { return previousBlockIndex_; }
151 /*! \brief Set the last block index data was added to.
153 * \param[in] blockIndex The block index.
155 void setPreviousBlockIndex(int blockIndex) { previousBlockIndex_ = blockIndex; }
157 /*! \brief Return sums for each coordinate dimension. */
158 const std::vector<CoordData>& coordData() const { return coordData_; }
160 /*! \brief Return the correlation integral tensor. */
161 const std::vector<double>& correlationIntegral() const { return correlationIntegral_; }
164 /* Weight sum data, indentical for all dimensions */
165 double blockSumWeight_; /**< Sum weights for current block. */
166 double blockSumSquareWeight_; /**< Sum weights^2 for current block. */
167 double sumOverBlocksSquareBlockWeight_; /**< Sum over all blocks in the simulation of block weight^2. */
168 double sumOverBlocksBlockSquareWeight_; /**< Sum over all blocks in the simulation of weight^2. */
169 double blockLength_; /**< The length of each block used for block averaging */
170 int previousBlockIndex_; /**< The last block index data was added to (needed only for block length in terms of time). */
172 /* Sums for each coordinate dimension. */
173 std::vector<CoordData> coordData_; /**< Array with sums for each coordinate dimension. */
175 /* Correlation tensor. */
176 std::vector<double> correlationIntegral_; /**< Array with the correlation elements corr(x, y) in the tensor, where x, y are vector components. */
180 * \brief Correlation data for computing the correlation tensor of one grid point.
182 * The time integrated autocorrelation of the desired quantity is computed using
183 * block averages, which is a computationally efficient and low memory method.
184 * Most of the work here goes into computing the block averages for weights
185 * and the coordinate quantity. This is done for a number of blocks in
186 * the range of \p c_numCorrelationBlocks/2 + 1 to \p c_numCorrelationBlocks,
187 * depending on the current simulation length. As the simulation time
188 * progresses, the blocks get longer. This is implemented in an efficient
189 * way by keeping track of log2(\p c_numCorrelationBlocks) \p BlockData
190 * data blocks with block length increasing progressively by a factor of 2.
191 * Once \p c_numCorrelationBlocks are reached, all block lengths are doubled.
193 class CorrelationTensor
196 /*! \brief 64 blocks is a good trade-off between signal and noise */
197 static constexpr int c_numCorrelationBlocks = 64;
199 /*! \brief Constructor.
201 * \param[in] numDim The dimensionality.
202 * \param[in] numBlockData The number of data block data structs.
203 * \param[in] blockLengthInit The initial block length.
205 CorrelationTensor(int numDim, int numBlockData, double blockLengthInit);
207 /*! \brief Get a const reference to the list of block data.
209 const std::vector<CorrelationBlockData>& blockDataList() const { return blockDataList_; }
212 * Get the total weight of the data in the correlation matrix.
214 * \returns the weight of the added data.
216 double getWeight() const
218 /* The last blockdata has only 1 block containing all data */
219 return blockDataList().back().blockSumWeight();
222 /*! \brief Restore a correlation element from history.
224 * \param[in] blockDataBuffer The linear correlation grid data history buffer.
225 * \param[in,out] bufferIndex The index in \p blockDataBuffer to start reading, is increased with the number of blocks read.
227 void restoreFromHistory(const std::vector<CorrelationBlockDataHistory>& blockDataBuffer,
228 size_t* bufferIndex);
231 /*! \brief Updates the block length by doubling.
233 * The length of all blocks is doubled. This is achieved by removing
234 * the shortest block, moving all other blocks and duplicating
235 * the data of longest block to a nw block of double length (but
236 * currenly only half filled with data).
238 void doubleBlockLengths();
240 /*! \brief Updates the block length such that data fits.
242 * \param[in] samplingLength Sampling length of all data, in time or weight.
244 void updateBlockLengths(double samplingLength);
247 /*! \brief Adds a weighted data vector to one point in the correlation grid.
249 * \note To avoid rounding noise, data with weight smaller than 1e-6
252 * \param[in] weight The weight of the data.
253 * \param[in] data One data point for each grid dimension.
254 * \param[in] blockLengthInWeight If true, a block is measured in probability weight, otherwise
256 * \param[in] t The simulation time.
258 void addData(double weight, gmx::ArrayRef<const double> data, bool blockLengthInWeight, double t);
260 /*! \brief Returns an element of the time integrated correlation tensor at a given point in the grid.
262 * The units of the integral are time*(units of data)^2. This will be friction units time/length^2
263 * if the data unit is 1/length.
265 * The correlation index lists the elements of the upper-triangular
266 * correlation matrix row-wise, so e.g. in 3D:
267 * 0 (0,0), 1 (1,0), 2 (1,1), 3 (2,0), 4 (2,1), 5 (2,2).
268 * (TODO: this should ideally not have to be known by the caller.)
270 * \param[in] tensorIndex Index in the tensor.
271 * \param[in] dtSample The sampling interval length.
272 * \returns the integral.
274 double getTimeIntegral(int tensorIndex, double dtSample) const;
277 * Returns the volume element of the correlation metric.
279 * The matrix of the metric equals the time-integrated correlation matrix. The volume element of
280 * the metric therefore equals the square-root of the absolute value of its determinant
281 * according to the standard formula for a volume element in a metric space.
283 * Since the units of the metric matrix elements are time*(units of data)^2, the volume element
284 * has units of (sqrt(time)*(units of data))^(ndim of data).
286 * \param[in] dtSample The sampling interval length.
287 * \returns the volume element.
289 double getVolumeElement(double dtSample) const;
292 std::vector<CorrelationBlockData> blockDataList_; /**< The block data for different, consecutively doubling block lengths. */
297 #endif /* GMX_AWH_CORRELATIONTENSOR_H */