Clang-11 fixes for applied forces
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / correlationgrid.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief
39  * Declares the CorrelationGrid class to collect correlation statistics on a grid, using several block lengths.
40  *
41  * \author Viveca Lindahl
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_awh
44  */
45
46 #ifndef GMX_AWH_CORRELATIONGRID_H
47 #define GMX_AWH_CORRELATIONGRID_H
48
49 #include <cstddef>
50
51 #include <vector>
52
53 #include "gromacs/utility/basedefinitions.h"
54 #include "gromacs/utility/gmxassert.h"
55
56 #include "correlationtensor.h"
57
58 namespace gmx
59 {
60
61 template<typename>
62 class ArrayRef;
63 struct CorrelationGridHistory;
64
65 /*! \internal
66  * \brief BiasGrid of local correlation tensors.
67  *
68  * This class provides the means for a bias to interaction with the grid
69  * of correlation tensors. The grid should have the same number of points
70  * and the same dimensionality as the bias grid.
71  */
72 class CorrelationGrid
73 {
74 public:
75     //! Enum that sets how we measure block length.
76     enum class BlockLengthMeasure
77     {
78         Time,  //!< Measure block length in time.
79         Weight //!< Measure block length in sampled weight.
80     };
81
82     /*! \brief Constructor.
83      *
84      * \param[in] numPoints           Number of points in the grid.
85      * \param[in] numDims             Number of dimensions of the grid.
86      * \param[in] blockLengthInit     Initial length of the blocks used for block averaging.
87      * \param[in] blockLengthMeasure  Sets how we measure block length.
88      * \param[in] dtSample            Time step for sampling correlations.
89      */
90     CorrelationGrid(int                numPoints,
91                     int                numDims,
92                     double             blockLengthInit,
93                     BlockLengthMeasure blockLengthMeasure,
94                     double             dtSample);
95
96     /*! \brief Adds a weighted data vector to one point in the correlation grid.
97      *
98      * \param[in] pointIndex  Index of the point to add data to.
99      * \param[in] weight      Weight to assign to the data.
100      * \param[in] data        One data point for each grid dimension.
101      * \param[in] t           The time when the data was sampled.
102      */
103     void addData(int pointIndex, double weight, gmx::ArrayRef<const double> data, double t)
104     {
105         tensors_[pointIndex].addData(weight, data, blockLengthMeasure == BlockLengthMeasure::Weight, t);
106     }
107
108     /*! \brief Restores the correlation grid state from the correlation grid history.
109      *
110      * The setup in the history should match that of this simulation.
111      * If this is not the case, an exception is thrown.
112      *
113      * \param[in] correlationGridHistory  The correlation grid state history.
114      */
115     void restoreStateFromHistory(const CorrelationGridHistory& correlationGridHistory);
116
117     /*! \brief Returns the number of elements in the tensor: dim*(dim+1)/2.
118      */
119     int tensorSize() const
120     {
121         GMX_RELEASE_ASSERT(!tensors_.empty(), "Should only call tensorSize on a valid grid");
122
123         return tensors_[0].blockDataList()[0].correlationIntegral().size();
124     }
125
126     /*! \brief Returns the size of the block data list.
127      */
128     int blockDataListSize() const
129     {
130         GMX_RELEASE_ASSERT(!tensors_.empty(), "Should only call tensorSize on a valid grid");
131
132         return tensors_[0].blockDataList().size();
133     }
134
135     /*! \brief Get a const reference to the correlation grid data.
136      */
137     const std::vector<CorrelationTensor>& tensors() const { return tensors_; }
138
139     /* Right now the below functions are only used for an initial log printing. */
140
141     /*! \brief Get the current blocklength.
142      */
143     double getBlockLength() const;
144
145     /*! \brief Get the current number of blocks.
146      *
147      * If we have a finite block span we have a constant number of blocks,
148      * otherwise we are always adding more blocks (and we don't keep
149      * track of the number), so we return -1.
150      */
151     int getNumBlocks() const;
152
153     const double             dtSample;           /**< Time in between samples. */
154     const BlockLengthMeasure blockLengthMeasure; /**< The measure for the block length. */
155 private:
156     std::vector<CorrelationTensor> tensors_; /**< Correlation tensor grid */
157 };
158
159 } // namespace gmx
160
161 #endif /* GMX_AWH_CORRELATIONGRID_H */