2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,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.
39 * This file contains datatypes and function declarations necessary
40 * for AWH to interface with the grid code.
42 * The grid organizes spatial properties of the AWH coordinate points.
43 * This includes traversing points in a specific order, locating
44 * neighboring points and calculating distances. Multiple dimensions
45 * as well as periodic dimensions are supported.
47 * \todo: Replace this by a more generic grid class once that is available.
49 * \author Viveca Lindahl
50 * \author Berk Hess <hess@kth.se>
54 #ifndef GMX_AWH_GRID_H
55 #define GMX_AWH_GRID_H
60 #include "dimparams.h" /* This is needed for awh_dvec */
68 * \brief An axis, i.e. dimension, of the grid.
73 /*! \brief Constructor.
75 * The spacing and number of points are set such that we have
76 * at least the requested point density.
77 * Requesting 0 point density results in the minimum number
80 * \param[in] origin Starting value.
81 * \param[in] end End value.
82 * \param[in] period Period, pass 0 if not periodic.
83 * \param[in] pointDensity Requested number of point per unit of axis length.
85 GridAxis(double origin, double end, double period, double pointDensity);
87 /*! \brief Constructor.
89 * \param[in] origin Starting value.
90 * \param[in] end End value.
91 * \param[in] period Period, pass 0 if not periodic.
92 * \param[in] numPoints The number of points.
94 GridAxis(double origin, double end, double period, int numPoints);
96 /*! \brief Returns if the axis has periodic boundaries.
98 bool isPeriodic() const { return period_ > 0; }
100 /*! \brief Returns the period of the grid along the axis.
102 double period() const { return period_; }
104 /*! \brief Returns the grid origin along the axis.
106 double origin() const { return origin_; }
108 /*! \brief Returns the grid point spacing along the axis.
110 double spacing() const { return spacing_; }
112 /*! \brief Returns the number of grid points along the axis.
114 int numPoints() const { return numPoints_; }
116 /*! \brief Returns the period of the grid in points along the axis.
118 * Returns 0 if the axis is not periodic.
120 int numPointsInPeriod() const { return numPointsInPeriod_; }
122 /*! \brief Returns the length of the interval.
124 * This returns the distance obtained by connecting the origin point to
125 * the end point in the positive direction. Note that this is generally
126 * not the shortest distance. For period > 0, both origin and
127 * end are expected to take values in the same periodic interval.
129 double length() const { return length_; }
131 /*! \brief Map a value to the nearest point index along an axis.
133 * \param[in] value Value along the axis.
134 * \returns the index nearest to the value.
136 int nearestIndex(double value) const;
139 double origin_; /**< Interval start value */
140 double length_; /**< Interval length */
141 double period_; /**< The period, 0 if not periodic */
142 double spacing_; /**< Point spacing */
143 int numPoints_; /**< Number of points in the interval */
144 int numPointsInPeriod_; /**< Number of points in a period (0 if no periodicity) */
148 * \brief A point in the grid.
150 * A grid point has a coordinate value and a coordinate index of the same dimensionality as the
151 * grid. It knows the linear indices of its neighboring point (which are useful only when handed up
156 awh_dvec coordValue; /**< Multidimensional coordinate value of this point */
157 awh_ivec index; /**< Multidimensional point indices */
158 std::vector<int> neighbor; /**< Linear point indices of the neighboring points */
162 * \brief The grid, generally multidimensional and periodic.
164 * The grid discretizes a multidimensional space with some given resolution.
165 * Each dimension is represented by an axis which sets the spatial extent,
166 * point spacing and periodicity of the grid in that direction.
171 /*! \brief Initializes the grid points.
177 * The point density per sigma of the Gaussian distribution in an umbrella.
179 * This value should be at least 1 to uniformly cover the reaction coordinate
180 * range with density and having it larger than 1 does not add information.
182 static constexpr double c_numPointsPerSigma = 1.0;
184 //! Cut-off in sigma for considering points, neglects 4e-8 of the density.
185 static constexpr double c_scopeCutoff = 5.5;
187 /*! \brief Construct a grid using AWH input parameters.
189 * \param[in] dimParams Dimension parameters including the expected inverse variance of the coordinate living on the grid (determines the grid spacing).
190 * \param[in] awhDimParams Dimension params from inputrec.
192 Grid(const std::vector<DimParams>& dimParams, const AwhDimParams* awhDimParams);
194 /*! \brief Returns the number of points in the grid.
196 * \returns the number of points in the grid.
198 size_t numPoints() const { return point_.size(); }
200 /*! \brief Returns a reference to a point on the grid.
202 * \returns a constant reference to a point on the grid.
204 const GridPoint& point(size_t pointIndex) const { return point_[pointIndex]; }
206 /*! \brief Returns the dimensionality of the grid.
208 * \returns the dimensionality of the grid.
210 int numDimensions() const { return axis_.size(); }
212 /*! \brief Returns the grid axes.
214 * \returns a constant reference to the grid axes.
216 const std::vector<GridAxis>& axis() const { return axis_; }
218 /*! \brief Returns a grid axis.
220 * param[in] dim Dimension to return the grid axis for.
221 * \returns a constant reference to the grid axis.
223 const GridAxis& axis(int dim) const { return axis_[dim]; }
225 /*! \brief Find the grid point with value nearest to the given value.
227 * \param[in] value Value vector.
228 * \returns the grid point index.
230 int nearestIndex(const awh_dvec value) const;
232 /*! \brief Query if the value is in the grid.
234 * \param[in] value Value vector.
235 * \returns true if the value is in the grid.
236 * \note It is assumed that any periodicity of value has already been taken care of.
238 bool covers(const awh_dvec value) const;
241 std::vector<GridPoint> point_; /**< Points on the grid */
242 std::vector<GridAxis> axis_; /**< Axes, one for each dimension. */
247 /*! \brief Convert a multidimensional grid point index to a linear one.
249 * \param[in] grid The grid.
250 * \param[in] indexMulti Multidimensional grid point index to convert to a linear one.
251 * \returns the linear index.
253 int multiDimGridIndexToLinear(const Grid& grid, const awh_ivec indexMulti);
255 /*! \brief Convert multidimensional array index to a linear one.
257 * \param[in] indexMulti Multidimensional index to convert to a linear one.
258 * \param[in] numDim Number of dimensions of the array.
259 * \param[in] numPointsDim Number of points of the array.
260 * \returns the linear index.
261 * \note This function can be used without having an initialized grid.
263 int multiDimArrayIndexToLinear(const awh_ivec indexMulti, int numDim, const awh_ivec numPointsDim);
265 /*! \brief Convert a linear grid point index to a multidimensional one.
267 * \param[in] grid The grid.
268 * \param[in] indexLinear Linear grid point index to convert to a multidimensional one.
269 * \param[out] indexMulti The multidimensional index.
271 void linearGridindexToMultiDim(const Grid& grid, int indexLinear, awh_ivec indexMulti);
273 /*! \brief Convert a linear array index to a multidimensional one.
275 * \param[in] indexLinear Linear array index
276 * \param[in] ndim Number of dimensions of the array.
277 * \param[in] numPointsDim Number of points for each dimension.
278 * \param[out] indexMulti The multidimensional index.
280 void linearArrayIndexToMultiDim(int indexLinear, int ndim, const awh_ivec numPointsDim, awh_ivec indexMulti);
283 * Find the next grid point in the sub-part of the grid given a starting point.
285 * The given grid point index is updated to the next valid grid point index
286 * by traversing the sub-part of the grid, here termed the subgrid.
287 * Since the subgrid range might extend beyond the actual size of the grid,
288 * the subgrid is traversed until a point both in the subgrid and grid is
289 * found. If no point is found, the function returns false and the index is
290 * not modified. The starting point needs to be inside of the subgrid. However,
291 * if this index is not given, meaning < 0, then the search is initialized at
292 * the subgrid origin, i.e. in this case the "next" grid point index is
293 * defined to be the first common grid/subgrid point.
295 * \param[in] grid The grid.
296 * \param[in] subgridOrigin Vector locating the subgrid origin relative to the grid origin.
297 * \param[in] subgridNpoints Number of points along each subgrid dimension.
298 * \param[in,out] gridPointIndex Pointer to the starting/next grid point index.
299 * \returns true if the grid point was updated.
301 bool advancePointInSubgrid(const Grid& grid,
302 const awh_ivec subgridOrigin,
303 const awh_ivec subgridNpoints,
304 int* gridPointIndex);
306 /*! \brief Maps each point in the grid to a point in the data grid.
308 * This functions maps an AWH bias grid to a user provided input data grid.
309 * The value of data grid point i along dimension d is given by data[d][i].
310 * The number of dimensions of the data should equal that of the grid.
311 * A fatal error is thrown if extracting the data fails or the data does not cover the whole grid.
313 * \param[out] gridpointToDatapoint Array mapping each grid point to a data point index.
314 * \param[in] data 2D array in format ndim x ndatapoints with data grid point values.
315 * \param[in] numDataPoints Number of data points.
316 * \param[in] dataFilename The data filename.
317 * \param[in] grid The grid.
318 * \param[in] correctFormatMessage String to include in error message if extracting the data fails.
320 void mapGridToDataGrid(std::vector<int>* gridpointToDatapoint,
321 const double* const* data,
323 const std::string& dataFilename,
325 const std::string& correctFormatMessage);
328 * Get the deviation along one dimension from the given value to a point in the grid.
330 * \param[in] grid The grid.
331 * \param[in] dimIndex Dimensional index in [0, ndim -1].
332 * \param[in] pointIndex Grid point index.
333 * \param[in] value Value along the given dimension.
334 * \returns the deviation of the given value to the given point.
336 double getDeviationFromPointAlongGridAxis(const Grid& grid, int dimIndex, int pointIndex, double value);