Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / awh / correlationtensor.h
blob6d2716e5d618be56aebfa85316c63d30c1bd9283
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \internal \file
38 * \brief
39 * Declares the CorrelationTensor class for correlation tensor statistics.
41 * \author Viveca Lindahl
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_awh
46 #ifndef GMX_AWH_CORRELATIONTENSOR_H
47 #define GMX_AWH_CORRELATIONTENSOR_H
49 #include <cstddef>
51 #include <vector>
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
55 #include "gromacs/utility/gmxassert.h"
57 namespace gmx
60 struct CorrelationBlockDataHistory;
62 /*! \internal \brief Correlation block averaging weight-only data.
64 class CorrelationBlockData
66 public:
67 /*! \internal \brief Correlation block averaging data.
69 struct CoordData
71 /*! \brief Constructor.
73 CoordData() :
74 blockSumWeightX(0),
75 sumOverBlocksBlockWeightBlockWeightX(0)
79 double blockSumWeightX; /**< Weighted sum of x for current block. */
80 double sumOverBlocksBlockWeightBlockWeightX; /**< Sum over all blocks in the simulation of block weight times sum_wx. */
83 /*! \brief Constructor.
85 * \param[in] numDim The dimensionality.
86 * \param[in] blockLengthInit The initial block length.
88 CorrelationBlockData(int numDim,
89 double blockLengthInit) :
90 blockSumWeight_(0),
91 blockSumSquareWeight_(0),
92 sumOverBlocksSquareBlockWeight_(0),
93 sumOverBlocksBlockSquareWeight_(0),
94 blockLength_(blockLengthInit),
95 previousBlockIndex_(-1),
96 coordData_(numDim),
97 correlationIntegral_(numDim*(numDim + 1)/2)
101 /*! \brief Restore the state from history.
103 * \param[in] blockHistory The block data history containing the weight sums.
104 * \param[in] coordData The coordinate data.
105 * \param[in] correlationIntegral The correlation integral for all tensor elements.
107 void restoreFromHistory(const CorrelationBlockDataHistory &blockHistory,
108 const std::vector<CoordData> &coordData,
109 const std::vector<double> &correlationIntegral);
111 /*! \brief Adds a weighted data vector to one point in the correlation grid.
113 * \param[in] weight The weight of the data.
114 * \param[in] data One data point for each grid dimension.
116 void addData(double weight,
117 gmx::ArrayRef<const double> data)
119 GMX_ASSERT(data.size() == coordData_.size(), "Size of data should match the size of coordData");
121 blockSumWeight_ += weight;
122 blockSumSquareWeight_ += weight*weight;
124 for (size_t d = 0; d < coordData_.size(); d++)
126 coordData_[d].blockSumWeightX += weight*data[d];
130 /*! \brief Adds a filled data block to correlation time integral.
132 void addBlockToCorrelationIntegral();
134 /*! \brief Returns the sum weights for current block. */
135 double blockSumWeight() const
137 return blockSumWeight_;
140 /*! \brief Returns the sum weights^2 for current block. */
141 double blockSumSquareWeight() const
143 return blockSumSquareWeight_;
146 /*! \brief Returns the sum over blocks of block weight^2. */
147 double sumOverBlocksSquareBlockWeight() const
149 return sumOverBlocksSquareBlockWeight_;
152 /*! \brief Returns the sum over blocks of weight^2. */
153 double sumOverBlocksBlockSquareWeight() const
155 return sumOverBlocksBlockSquareWeight_;
158 /*! \brief Returns the length of each block used for block averaging. */
159 double blockLength() const
161 return blockLength_;
164 /*! \brief Double the length of each block used for block averaging. */
165 void doubleBlockLength()
167 blockLength_ *= 2;
170 /*! \brief Return the last block index data was added to (needed only for block length in terms of time). */
171 int previousBlockIndex() const
173 return previousBlockIndex_;
176 /*! \brief Set the last block index data was added to.
178 * \param[in] blockIndex The block index.
180 void setPreviousBlockIndex(int blockIndex)
182 previousBlockIndex_ = blockIndex;
185 /*! \brief Return sums for each coordinate dimension. */
186 const std::vector<CoordData> &coordData() const
188 return coordData_;
191 /*! \brief Return the correlation integral tensor. */
192 const std::vector<double> &correlationIntegral() const
194 return correlationIntegral_;
197 private:
198 /* Weight sum data, indentical for all dimensions */
199 double blockSumWeight_; /**< Sum weights for current block. */
200 double blockSumSquareWeight_; /**< Sum weights^2 for current block. */
201 double sumOverBlocksSquareBlockWeight_; /**< Sum over all blocks in the simulation of block weight^2. */
202 double sumOverBlocksBlockSquareWeight_; /**< Sum over all blocks in the simulation of weight^2. */
203 double blockLength_; /**< The length of each block used for block averaging */
204 int previousBlockIndex_; /**< The last block index data was added to (needed only for block length in terms of time). */
206 /* Sums for each coordinate dimension. */
207 std::vector<CoordData> coordData_; /**< Array with sums for each coordinate dimension. */
209 /* Correlation tensor. */
210 std::vector<double> correlationIntegral_; /**< Array with the correlation elements corr(x, y) in the tensor, where x, y are vector components. */
213 /*! \internal
214 * \brief Correlation data for computing the correlation tensor of one grid point.
216 * The time integrated autocorrelation of the desired quantity is computed using
217 * block averages, which is a computationally efficient and low memory method.
218 * Most of the work here goes into computing the block averages for weights
219 * and the coordinate quantity. This is done for a number of blocks in
220 * the range of \p c_numCorrelationBlocks/2 + 1 to \p c_numCorrelationBlocks,
221 * depending on the current simulation length. As the simulation time
222 * progresses, the blocks get longer. This is implemented in an efficient
223 * way by keeping track of log2(\p c_numCorrelationBlocks) \p BlockData
224 * data blocks with block length increasing progressively by a factor of 2.
225 * Once \p c_numCorrelationBlocks are reached, all block lengths are doubled.
227 class CorrelationTensor
229 public:
230 /*! \brief 64 blocks is a good trade-off between signal and noise */
231 static constexpr int c_numCorrelationBlocks = 64;
233 /*! \brief Constructor.
235 * \param[in] numDim The dimensionality.
236 * \param[in] numBlockData The number of data block data structs.
237 * \param[in] blockLengthInit The initial block length.
239 CorrelationTensor(int numDim,
240 int numBlockData,
241 double blockLengthInit);
243 /*! \brief Get a const reference to the list of block data.
245 const std::vector<CorrelationBlockData> &blockDataList() const
247 return blockDataList_;
250 /*! \brief
251 * Get the total weight of the data in the correlation matrix.
253 * \returns the weight of the added data.
255 double getWeight() const
257 /* The last blockdata has only 1 block containing all data */
258 return blockDataList().back().blockSumWeight();
261 /*! \brief Restore a correlation element from history.
263 * \param[in] blockDataBuffer The linear correlation grid data history buffer.
264 * \param[in,out] bufferIndex The index in \p blockDataBuffer to start reading, is increased with the number of blocks read.
266 void restoreFromHistory(const std::vector<CorrelationBlockDataHistory> &blockDataBuffer,
267 size_t *bufferIndex);
269 private:
270 /*! \brief Updates the block length by doubling.
272 * The length of all blocks is doubled. This is achieved by removing
273 * the shortest block, moving all other blocks and duplicating
274 * the data of longest block to a nw block of double length (but
275 * currenly only half filled with data).
277 void doubleBlockLengths();
279 /*! \brief Updates the block length such that data fits.
281 * \param[in] samplingLength Sampling length of all data, in time or weight.
283 void updateBlockLengths(double samplingLength);
285 public:
286 /*! \brief Adds a weighted data vector to one point in the correlation grid.
288 * \note To avoid rounding noise, data with weight smaller than 1e-6
289 * is ignored.
291 * \param[in] weight The weight of the data.
292 * \param[in] data One data point for each grid dimension.
293 * \param[in] blockLengthInWeight If true, a block is measured in probability weight, otherwise in time.
294 * \param[in] t The simulation time.
296 void addData(double weight,
297 gmx::ArrayRef<const double> data,
298 bool blockLengthInWeight,
299 double t);
301 /*! \brief Returns an element of the time integrated correlation tensor at a given point in the grid.
303 * The units of the integral are time*(units of data)^2. This will be friction units time/length^2
304 * if the data unit is 1/length.
306 * The correlation index lists the elements of the upper-triangular
307 * correlation matrix row-wise, so e.g. in 3D:
308 * 0 (0,0), 1 (1,0), 2 (1,1), 3 (2,0), 4 (2,1), 5 (2,2).
309 * (TODO: this should ideally not have to be known by the caller.)
311 * \param[in] tensorIndex Index in the tensor.
312 * \param[in] dtSample The sampling interval length.
313 * \returns the integral.
315 double getTimeIntegral(int tensorIndex,
316 double dtSample) const;
318 /*! \brief
319 * Returns the volume element of the correlation metric.
321 * The matrix of the metric equals the time-integrated correlation matrix. The volume element of
322 * the metric therefore equals the square-root of the absolute value of its determinant according
323 * to the standard formula for a volume element in a metric space.
325 * Since the units of the metric matrix elements are time*(units of data)^2, the volume element has
326 * units of (sqrt(time)*(units of data))^(ndim of data).
328 * \param[in] dtSample The sampling interval length.
329 * \returns the volume element.
331 double getVolumeElement(double dtSample) const;
333 private:
334 std::vector<CorrelationBlockData> blockDataList_; /**< The block data for different, consecutively doubling block lengths. */
337 } // namespace gmx
339 #endif /* GMX_AWH_CORRELATIONTENSOR_H */