Add SIMD for AWH
[gromacs.git] / src / gromacs / awh / coordstate.h
blob3df240ac1d30545929231ce29193a9c0d81dd203
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017, 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 CoordState class.
41 * It sets and holds the current coordinate value and corresponding closest
42 * grid point index. These are (re)set at every step.
43 * With umbrella potential type, this class also holds and updates the umbrella
44 * potential reference location, which is a state variable that presists over
45 * the duration of an AWH sampling interval.
47 * \author Viveca Lindahl
48 * \author Berk Hess <hess@kth.se>
49 * \ingroup module_awh
52 #ifndef GMX_AWH_COORDSTATE_H
53 #define GMX_AWH_COORDSTATE_H
55 #include <vector>
57 #include "gromacs/utility/arrayref.h"
59 #include "dimparams.h"
61 namespace gmx
64 struct AwhBiasParams;
65 struct AwhBiasStateHistory;
66 class BiasParams;
67 class Grid;
69 /*! \internal \brief Keeps track of the current coordinate value, grid index and umbrella location.
71 class CoordState
73 public:
74 /*! \brief Constructor.
76 * \param[in] awhBiasParams The Bias parameters from inputrec.
77 * \param[in] dimParams The dimension Parameters.
78 * \param[in] grid The grid.
80 CoordState(const AwhBiasParams &awhBiasParams,
81 const std::vector<DimParams> &dimParams,
82 const Grid &grid);
84 /*! \brief
85 * Sample a new umbrella reference point given the current coordinate value.
87 * It is assumed that the probability distribution has been updated.
89 * \param[in] grid The grid.
90 * \param[in] gridpointIndex The grid point, sets the neighborhood.
91 * \param[in] probWeightNeighbor Probability weights of the neighbors.
92 * \param[in] step Step number, needed for the random number generator.
93 * \param[in] seed Random seed.
94 * \param[in] indexSeed Second random seed, should be the bias Index.
95 * \returns the index of the sampled point.
97 void sampleUmbrellaGridpoint(const Grid &grid,
98 int gridpointIndex,
99 gmx::ArrayRef<const double> probWeightNeighbor,
100 gmx_int64_t step,
101 gmx_int64_t seed,
102 int indexSeed);
104 /*! \brief Update the coordinate value with coordValue.
106 * \param[in] grid The grid.
107 * \param[in] coordValue The new coordinate value.
109 void setCoordValue(const Grid &grid,
110 const awh_dvec coordValue);
112 /*! \brief Restores the coordinate state from history.
114 * \param[in] stateHistory The AWH bias state history.
116 void restoreFromHistory(const AwhBiasStateHistory &stateHistory);
118 /*! \brief Returns the current coordinate value.
120 const awh_dvec &coordValue() const
122 return coordValue_;
125 /*! \brief Returns the grid point index for the current coordinate value.
127 int gridpointIndex() const
129 return gridpointIndex_;
132 /*! \brief Returns the index for the current reference grid point.
134 int umbrellaGridpoint() const
136 return umbrellaGridpoint_;
139 private:
140 awh_dvec coordValue_; /**< Current coordinate value in (nm or rad) */
141 int gridpointIndex_; /**< The grid point index for the current coordinate value */
142 int umbrellaGridpoint_; /**< Index for the current reference grid point for the umbrella, only used with umbrella potential type */
145 } // namespace gmx
147 #endif /* GMX_AWH_COORDSTATE_H */