SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / coordstate.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 CoordState class.
40  *
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.
46  *
47  * \author Viveca Lindahl
48  * \author Berk Hess <hess@kth.se>
49  * \ingroup module_awh
50  */
51
52 #ifndef GMX_AWH_COORDSTATE_H
53 #define GMX_AWH_COORDSTATE_H
54
55 #include <vector>
56
57 #include "dimparams.h"
58
59 namespace gmx
60 {
61
62 template<typename>
63 class ArrayRef;
64 class AwhBiasParams;
65 struct AwhBiasStateHistory;
66 class BiasParams;
67 class BiasGrid;
68
69 /*! \internal \brief Keeps track of the current coordinate value, grid index and umbrella location.
70  */
71 class CoordState
72 {
73 public:
74     /*! \brief Constructor.
75      *
76      * \param[in] awhBiasParams  The Bias parameters from inputrec.
77      * \param[in] dimParams      The dimension Parameters.
78      * \param[in] grid           The grid.
79      */
80     CoordState(const AwhBiasParams& awhBiasParams, ArrayRef<const DimParams> dimParams, const BiasGrid& grid);
81
82     /*! \brief
83      * Sample a new umbrella reference point given the current coordinate value.
84      *
85      * It is assumed that the probability distribution has been updated.
86      *
87      * \param[in] grid                The grid.
88      * \param[in] gridpointIndex      The grid point, sets the neighborhood.
89      * \param[in] probWeightNeighbor  Probability weights of the neighbors.
90      * \param[in] step                Step number, needed for the random number generator.
91      * \param[in] seed                Random seed.
92      * \param[in] indexSeed           Second random seed, should be the bias Index.
93      * \returns the index of the sampled point.
94      */
95     void sampleUmbrellaGridpoint(const BiasGrid&             grid,
96                                  int                         gridpointIndex,
97                                  gmx::ArrayRef<const double> probWeightNeighbor,
98                                  int64_t                     step,
99                                  int64_t                     seed,
100                                  int                         indexSeed);
101
102     /*! \brief Update the coordinate value with coordValue.
103      *
104      * \param[in] grid        The grid.
105      * \param[in] coordValue  The new coordinate value.
106      */
107     void setCoordValue(const BiasGrid& grid, const awh_dvec coordValue);
108
109     /*! \brief Restores the coordinate state from history.
110      *
111      * \param[in] stateHistory  The AWH bias state history.
112      */
113     void restoreFromHistory(const AwhBiasStateHistory& stateHistory);
114
115     /*! \brief Returns the current coordinate value.
116      */
117     const awh_dvec& coordValue() const { return coordValue_; }
118
119     /*! \brief Returns the grid point index for the current coordinate value.
120      */
121     int gridpointIndex() const { return gridpointIndex_; }
122
123     /*! \brief Returns the index for the current reference grid point.
124      */
125     int umbrellaGridpoint() const { return umbrellaGridpoint_; }
126
127     /*! \brief Sets the umbrella grid point to the current grid point
128      */
129     void setUmbrellaGridpointToGridpoint();
130
131 private:
132     awh_dvec coordValue_;        /**< Current coordinate value in (nm or rad) */
133     int      gridpointIndex_;    /**< The grid point index for the current coordinate value */
134     int      umbrellaGridpoint_; /**< Index for the current reference grid point for the umbrella, only used with umbrella potential type */
135 };
136
137 } // namespace gmx
138
139 #endif /* GMX_AWH_COORDSTATE_H */