SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / mdtypes / awh_params.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,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 /*! \libinternal \file
38  *
39  * \brief
40  * Declares AWH parameter data types.
41  *
42  * Besides internal use by the AWH module, the AWH parameters are needed
43  * for reading the user input (mdp) file and for reading and writing the
44  * parameters to the mdrun input (tpr) file.
45  *
46  * \author Viveca Lindahl
47  * \inlibraryapi
48  * \ingroup module_mdtypes
49  */
50
51 #ifndef GMX_MDTYPES_AWH_PARAMS_H
52 #define GMX_MDTYPES_AWH_PARAMS_H
53
54 #include <vector>
55
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/basedefinitions.h"
59 #include "gromacs/utility/classhelpers.h"
60
61 struct t_inpfile;
62 struct t_inputrec;
63 struct pull_params_t;
64 using warninp_t = struct warninp*;
65
66 namespace gmx
67 {
68
69 class ISerializer;
70 //! Target distribution enum.
71 enum class AwhTargetType : int
72 {
73     Constant,
74     Cutoff,
75     Boltzmann,
76     LocalBoltzmann,
77     Count,
78     Default = Constant
79 };
80 //! String for target distribution.
81 const char* enumValueToString(AwhTargetType enumValue);
82
83 //! Weight histogram growth enum.
84 enum class AwhHistogramGrowthType : int
85 {
86     ExponentialLinear,
87     Linear,
88     Count,
89     Default = ExponentialLinear
90 };
91 //! String for weight histogram growth
92 const char* enumValueToString(AwhHistogramGrowthType enumValue);
93
94 //! AWH potential type enum.
95 enum class AwhPotentialType : int
96 {
97     Convolved,
98     Umbrella,
99     Count,
100     Default = Convolved
101 };
102 //! String for AWH potential type
103 const char* enumValueToString(AwhPotentialType enumValue);
104
105 //! AWH bias reaction coordinate provider
106 enum class AwhCoordinateProviderType : int
107 {
108     Pull,
109     FreeEnergyLambda,
110     Count,
111     Default = Pull
112 };
113 //! String for AWH bias reaction coordinate provider.
114 const char* enumValueToString(AwhCoordinateProviderType enumValue);
115
116 class AwhDimParams
117 {
118 public:
119     //! Constructor from input file.
120     AwhDimParams(std::vector<t_inpfile>* inp, const std::string& prefix, warninp_t wi, bool bComment);
121     //! Constructor to generate from file reading.
122     explicit AwhDimParams(ISerializer* serializer);
123
124     //! Move constructor.
125     AwhDimParams(AwhDimParams&&) = default;
126     //! Move assignment operator.
127     AwhDimParams& operator=(AwhDimParams&&) = default;
128     //! Delete copy constructor.
129     AwhDimParams(const AwhDimParams&) = delete;
130     //! Delete copy assignment.
131     AwhDimParams& operator=(const AwhDimParams&) = delete;
132
133     //! Which module is providing the reaction coordinate.
134     AwhCoordinateProviderType coordinateProvider() const { return eCoordProvider_; }
135     //! Index for reaction coordinate in provider.
136     int coordinateIndex() const { return coordIndex_; }
137     //! Start value for interval.
138     double origin() const { return origin_; }
139     //! End value for interval.
140     double end() const { return end_; }
141     //! Period for the dimension.
142     double period() const { return period_; }
143     //! Set period value dependent on state.
144     void setPeriod(double period) { period_ = period; }
145     //! Force constant for this dimension.
146     double forceConstant() const { return forceConstant_; }
147     //! Estimated diffusion constant.
148     double diffusion() const { return diffusion_; }
149     //! Initial value for coordinate.
150     double initialCoordinate() const { return coordValueInit_; }
151     //! Set initial coordinate value dependent on state.
152     void setInitialCoordinate(double initialCoordinate) { coordValueInit_ = initialCoordinate; }
153     //! Diameter needed to be sampled.
154     double coverDiameter() const { return coverDiameter_; }
155     //! Write datastructure.
156     void serialize(ISerializer* serializer);
157
158 private:
159     //! The module providing the reaction coordinate.
160     AwhCoordinateProviderType eCoordProvider_;
161     //! Index of reaction coordinate in the provider.
162     int coordIndex_ = 0;
163     //! Start value of the interval.
164     double origin_ = 0.0;
165     //! End value of the interval.
166     double end_ = 0.0;
167     //! The period of this dimension (= 0 if not periodic).
168     double period_ = 0.0;
169     //! The force constant in kJ/mol/nm^2, kJ/mol/rad^2
170     double forceConstant_ = 0.0;
171     //! Estimated diffusion constant in units of nm^2/ps or rad^2/ps or ps^-1.
172     double diffusion_ = 0.0;
173     //! The initial coordinate value.
174     double coordValueInit_ = 0.0;
175     //! The diameter that needs to be sampled around a point before it is considered covered.
176     double coverDiameter_ = 0.0;
177 };
178
179 class AwhBiasParams
180 {
181 public:
182     //! Constructor from input file.
183     AwhBiasParams(std::vector<t_inpfile>* inp, const std::string& prefix, warninp_t wi, bool bComment);
184     //! Constructor to generate from file reading.
185     explicit AwhBiasParams(ISerializer* serializer);
186
187     //! Move constructor.
188     AwhBiasParams(AwhBiasParams&&) = default;
189     //! Move assignment operator.
190     AwhBiasParams& operator=(AwhBiasParams&&) = default;
191     //! Delete copy constructor.
192     AwhBiasParams(const AwhBiasParams&) = delete;
193     //! Delete copy assignment.
194     AwhBiasParams& operator=(const AwhBiasParams&) = delete;
195
196     //! Which target distribution is searched.
197     AwhTargetType targetDistribution() const { return eTarget_; }
198     //! Beta scaling to reach target distribution.
199     double targetBetaScaling() const { return targetBetaScaling_; }
200     //! Cutoff for target.
201     double targetCutoff() const { return targetCutoff_; }
202     //! Which kind of growth to use.
203     AwhHistogramGrowthType growthType() const { return eGrowth_; }
204     //! User provided PMF estimate.
205     bool userPMFEstimate() const { return bUserData_; }
206     //! Estimated initial free energy error in kJ/mol.
207     double initialErrorEstimate() const { return errorInitial_; }
208     //! Dimensions of coordinate space.
209     int ndim() const { return dimParams_.size(); }
210     //! Number of groups to share this bias with.
211     int shareGroup() const { return shareGroup_; }
212     //! If the simulation starts with equilibrating histogram.
213     bool equilibrateHistogram() const { return equilibrateHistogram_; }
214     //! Access to dimension parameters.
215     ArrayRef<AwhDimParams> dimParams() { return dimParams_; }
216     //! Const access to dimension parameters.
217     ArrayRef<const AwhDimParams> dimParams() const { return dimParams_; }
218     //! Write datastructure.
219     void serialize(ISerializer* serializer);
220
221 private:
222     //! AWH parameters per dimension.
223     std::vector<AwhDimParams> dimParams_;
224     //! Type of target distribution.
225     AwhTargetType eTarget_;
226     //! Beta scaling value for Boltzmann type target distributions.
227     double targetBetaScaling_;
228     //! Free energy cutoff value for cutoff type target distribution in kJ/mol.
229     double targetCutoff_;
230     //! How the biasing histogram grows.
231     AwhHistogramGrowthType eGrowth_;
232     //! Is there a user-defined initial PMF estimate and target estimate?
233     bool bUserData_;
234     //! Estimated initial free energy error in kJ/mol.
235     double errorInitial_;
236     //! When >0, the bias is shared with biases of the same group and across multiple simulations when shareBiasMultisim=true
237     int shareGroup_;
238     //! True if the simulation starts out by equilibrating the histogram.
239     bool equilibrateHistogram_;
240 };
241 /*! \internal
242  * \brief Structure holding parameter information for AWH.
243  */
244 class AwhParams
245 {
246 public:
247     //! Constructor from input file.
248     AwhParams(std::vector<t_inpfile>* inp, warninp_t wi);
249     //! Constructor used to generate awh parameter from file reading.
250     explicit AwhParams(ISerializer* serializer);
251
252     //! Move constructor.
253     AwhParams(AwhParams&&) = default;
254     //! Move assignment operator.
255     AwhParams& operator=(AwhParams&&) = default;
256     //! Delete copy constructor.
257     AwhParams(const AwhParams&) = delete;
258     //! Delete copy assignment.
259     AwhParams& operator=(const AwhParams&) = delete;
260
261     //! Get number of biases.
262     int numBias() const { return awhBiasParams_.size(); }
263     //! Get access to bias parameters.
264     ArrayRef<AwhBiasParams> awhBiasParams() { return awhBiasParams_; }
265     //! Const access to bias parameters.
266     ArrayRef<const AwhBiasParams> awhBiasParams() const { return awhBiasParams_; }
267     //! What king of potential is being used. \todo should use actual enum class.
268     AwhPotentialType potential() const { return potentialEnum_; }
269     //! Seed used for starting AWH.
270     int64_t seed() const { return seed_; }
271     //! Output step interval.
272     int nstout() const { return nstOut_; }
273     //! Number of samples per coordinate sample.
274     int nstSampleCoord() const { return nstSampleCoord_; }
275     //! Number of samples per free energy update.
276     int numSamplesUpdateFreeEnergy() const { return numSamplesUpdateFreeEnergy_; }
277     //! If biases are shared in multisim.
278     bool shareBiasMultisim() const { return shareBiasMultisim_; }
279     //! Serialize awh parameters.
280     void serialize(ISerializer* serializer);
281
282 private:
283     //! AWH bias parameters.
284     std::vector<AwhBiasParams> awhBiasParams_;
285     //! Random seed.
286     int64_t seed_;
287     //! Output step interval.
288     int nstOut_;
289     //! Number of samples per coordinate sample (also used for PMF)
290     int nstSampleCoord_;
291     //! Number of samples per free energy update.
292     int numSamplesUpdateFreeEnergy_;
293     //! Type of potential.
294     AwhPotentialType potentialEnum_;
295     //! Whether to share biases with shareGroup>0 between multi-simulations.
296     bool shareBiasMultisim_;
297 };
298
299 } // namespace gmx
300
301 #endif /* GMX_MDTYPES_AWH_PARAMS_H */