96f680428d80d288dee41e48795ffa80b9a71373
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / correlationtensor.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36
37 /*! \internal \file
38  *
39  * \brief
40  * Declares the CorrelationTensor class for correlation tensor statistics.
41  *
42  * \author Viveca Lindahl
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_awh
45  */
46
47 #ifndef GMX_AWH_CORRELATIONTENSOR_H
48 #define GMX_AWH_CORRELATIONTENSOR_H
49
50 #include <cstddef>
51
52 #include <vector>
53
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/gmxassert.h"
57
58 namespace gmx
59 {
60
61 struct CorrelationBlockDataHistory;
62
63 /*! \internal \brief Correlation block averaging weight-only data.
64  */
65 class CorrelationBlockData
66 {
67 public:
68     /*! \internal \brief Correlation block averaging data.
69      */
70     struct CoordData
71     {
72         /*! \brief Constructor.
73          */
74         CoordData() : blockSumWeightX(0), sumOverBlocksBlockWeightBlockWeightX(0) {}
75
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. */
78     };
79
80     /*! \brief Constructor.
81      *
82      * \param[in] numDim           The dimensionality.
83      * \param[in] blockLengthInit  The initial block length.
84      */
85     CorrelationBlockData(int numDim, double blockLengthInit) :
86         blockSumWeight_(0),
87         blockSumSquareWeight_(0),
88         sumOverBlocksSquareBlockWeight_(0),
89         sumOverBlocksBlockSquareWeight_(0),
90         blockLength_(blockLengthInit),
91         previousBlockIndex_(-1),
92         coordData_(numDim),
93         correlationIntegral_(numDim * (numDim + 1) / 2)
94     {
95     }
96
97     /*! \brief Restore the state from history.
98      *
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.
102      */
103     void restoreFromHistory(const CorrelationBlockDataHistory& blockHistory,
104                             const std::vector<CoordData>&      coordData,
105                             const std::vector<double>&         correlationIntegral);
106
107     /*! \brief Adds a weighted data vector to one point in the correlation grid.
108      *
109      * \param[in] weight  The weight of the data.
110      * \param[in] data    One data point for each grid dimension.
111      */
112     void addData(double weight, gmx::ArrayRef<const double> data)
113     {
114         GMX_ASSERT(data.size() == coordData_.size(),
115                    "Size of data should match the size of coordData");
116
117         blockSumWeight_ += weight;
118         blockSumSquareWeight_ += weight * weight;
119
120         for (size_t d = 0; d < coordData_.size(); d++)
121         {
122             coordData_[d].blockSumWeightX += weight * data[d];
123         }
124     }
125
126     /*! \brief Adds a filled data block to correlation time integral.
127      */
128     void addBlockToCorrelationIntegral();
129
130     /*! \brief Returns the sum weights for current block. */
131     double blockSumWeight() const { return blockSumWeight_; }
132
133     /*! \brief Returns the sum weights^2 for current block. */
134     double blockSumSquareWeight() const { return blockSumSquareWeight_; }
135
136     /*! \brief Returns the sum over blocks of block weight^2. */
137     double sumOverBlocksSquareBlockWeight() const { return sumOverBlocksSquareBlockWeight_; }
138
139     /*! \brief Returns the sum over blocks of weight^2. */
140     double sumOverBlocksBlockSquareWeight() const { return sumOverBlocksBlockSquareWeight_; }
141
142     /*! \brief Returns the length of each block used for block averaging. */
143     double blockLength() const { return blockLength_; }
144
145     /*! \brief Double the length of each block used for block averaging. */
146     void doubleBlockLength() { blockLength_ *= 2; }
147
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_; }
150
151     /*! \brief Set the last block index data was added to.
152      *
153      * \param[in] blockIndex  The block index.
154      */
155     void setPreviousBlockIndex(int blockIndex) { previousBlockIndex_ = blockIndex; }
156
157     /*! \brief Return sums for each coordinate dimension. */
158     const std::vector<CoordData>& coordData() const { return coordData_; }
159
160     /*! \brief Return the correlation integral tensor. */
161     const std::vector<double>& correlationIntegral() const { return correlationIntegral_; }
162
163 private:
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). */
171
172     /* Sums for each coordinate dimension. */
173     std::vector<CoordData> coordData_; /**< Array with sums for each coordinate dimension. */
174
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. */
177 };
178
179 /*! \internal
180  * \brief Correlation data for computing the correlation tensor of one grid point.
181  *
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.
192  */
193 class CorrelationTensor
194 {
195 public:
196     /*! \brief 64 blocks is a good trade-off between signal and noise */
197     static constexpr int c_numCorrelationBlocks = 64;
198
199     /*! \brief Constructor.
200      *
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.
204      */
205     CorrelationTensor(int numDim, int numBlockData, double blockLengthInit);
206
207     /*! \brief Get a const reference to the list of block data.
208      */
209     const std::vector<CorrelationBlockData>& blockDataList() const { return blockDataList_; }
210
211     /*! \brief
212      * Get the total weight of the data in the correlation matrix.
213      *
214      * \returns the weight of the added data.
215      */
216     double getWeight() const
217     {
218         /* The last blockdata has only 1 block containing all data */
219         return blockDataList().back().blockSumWeight();
220     }
221
222     /*! \brief Restore a correlation element from history.
223      *
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.
226      */
227     void restoreFromHistory(const std::vector<CorrelationBlockDataHistory>& blockDataBuffer,
228                             size_t*                                         bufferIndex);
229
230 private:
231     /*! \brief Updates the block length by doubling.
232      *
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).
237      */
238     void doubleBlockLengths();
239
240     /*! \brief Updates the block length such that data fits.
241      *
242      * \param[in] samplingLength  Sampling length of all data, in time or weight.
243      */
244     void updateBlockLengths(double samplingLength);
245
246 public:
247     /*! \brief Adds a weighted data vector to one point in the correlation grid.
248      *
249      * \note To avoid rounding noise, data with weight smaller than 1e-6
250      *       is ignored.
251      *
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
255      *                                 in time.
256      * \param[in] t                    The simulation time.
257      */
258     void addData(double weight, gmx::ArrayRef<const double> data, bool blockLengthInWeight, double t);
259
260     /*! \brief Returns an element of the time integrated correlation tensor at a given point in the grid.
261      *
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.
264      *
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.)
269      *
270      * \param[in] tensorIndex  Index in the tensor.
271      * \param[in] dtSample     The sampling interval length.
272      * \returns the integral.
273      */
274     double getTimeIntegral(int tensorIndex, double dtSample) const;
275
276     /*! \brief
277      * Returns the volume element of the correlation metric.
278      *
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.
282      *
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).
285      *
286      * \param[in] dtSample  The sampling interval length.
287      * \returns the volume element.
288      */
289     double getVolumeElement(double dtSample) const;
290
291 private:
292     std::vector<CorrelationBlockData> blockDataList_; /**< The block data for different, consecutively doubling block lengths. */
293 };
294
295 } // namespace gmx
296
297 #endif /* GMX_AWH_CORRELATIONTENSOR_H */