2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /*! \libinternal \file
40 * Declares AWH parameter data types.
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.
46 * \author Viveca Lindahl
48 * \ingroup module_mdtypes
51 #ifndef GMX_MDTYPES_AWH_PARAMS_H
52 #define GMX_MDTYPES_AWH_PARAMS_H
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"
64 using warninp_t = struct warninp*;
70 //! Target distribution enum.
71 enum class AwhTargetType : int
80 //! String for target distribution.
81 const char* enumValueToString(AwhTargetType enumValue);
83 //! Weight histogram growth enum.
84 enum class AwhHistogramGrowthType : int
89 Default = ExponentialLinear
91 //! String for weight histogram growth
92 const char* enumValueToString(AwhHistogramGrowthType enumValue);
94 //! AWH potential type enum.
95 enum class AwhPotentialType : int
102 //! String for AWH potential type
103 const char* enumValueToString(AwhPotentialType enumValue);
105 //! AWH bias reaction coordinate provider
106 enum class AwhCoordinateProviderType : int
113 //! String for AWH bias reaction coordinate provider.
114 const char* enumValueToString(AwhCoordinateProviderType enumValue);
119 //! Constructor from input file.
120 AwhDimParams(std::vector<t_inpfile>* inp,
121 const std::string& prefix,
122 const t_inputrec& ir,
125 //! Constructor to generate from file reading.
126 explicit AwhDimParams(ISerializer* serializer);
128 //! Move constructor.
129 AwhDimParams(AwhDimParams&&) = default;
130 //! Move assignment operator.
131 AwhDimParams& operator=(AwhDimParams&&) = default;
132 //! Delete copy constructor.
133 AwhDimParams(const AwhDimParams&) = delete;
134 //! Delete copy assignment.
135 AwhDimParams& operator=(const AwhDimParams&) = delete;
137 //! Which module is providing the reaction coordinate.
138 AwhCoordinateProviderType coordinateProvider() const { return eCoordProvider_; }
139 //! Index for reaction coordinate in provider.
140 int coordinateIndex() const { return coordIndex_; }
141 //! Start value for interval.
142 double origin() const { return origin_; }
143 //! End value for interval.
144 double end() const { return end_; }
145 //! Period for the dimension.
146 double period() const { return period_; }
147 //! Set period value dependent on state.
148 void setPeriod(double period) { period_ = period; }
149 //! Force constant for this dimension.
150 double forceConstant() const { return forceConstant_; }
151 //! Estimated diffusion constant.
152 double diffusion() const { return diffusion_; }
153 //! Initial value for coordinate.
154 double initialCoordinate() const { return coordValueInit_; }
155 //! Set initial coordinate value dependent on state.
156 void setInitialCoordinate(double initialCoordinate) { coordValueInit_ = initialCoordinate; }
157 //! Diameter needed to be sampled.
158 double coverDiameter() const { return coverDiameter_; }
159 //! Write datastructure.
160 void serialize(ISerializer* serializer);
163 //! The module providing the reaction coordinate.
164 AwhCoordinateProviderType eCoordProvider_;
165 //! Index of reaction coordinate in the provider.
167 //! Start value of the interval.
169 //! End value of the interval.
171 //! The period of this dimension (= 0 if not periodic).
173 //! The force constant in kJ/mol/nm^2, kJ/mol/rad^2
174 double forceConstant_;
175 //! Estimated diffusion constant in units of nm^2/ps or rad^2/ps or ps^-1.
177 //! The initial coordinate value.
178 double coordValueInit_;
179 //! The diameter that needs to be sampled around a point before it is considered covered.
180 double coverDiameter_;
186 //! Constructor from input file.
187 AwhBiasParams(std::vector<t_inpfile>* inp,
188 const std::string& prefix,
189 const t_inputrec& ir,
192 //! Constructor to generate from file reading.
193 explicit AwhBiasParams(ISerializer* serializer);
195 //! Move constructor.
196 AwhBiasParams(AwhBiasParams&&) = default;
197 //! Move assignment operator.
198 AwhBiasParams& operator=(AwhBiasParams&&) = default;
199 //! Delete copy constructor.
200 AwhBiasParams(const AwhBiasParams&) = delete;
201 //! Delete copy assignment.
202 AwhBiasParams& operator=(const AwhBiasParams&) = delete;
204 //! Which target distribution is searched.
205 AwhTargetType targetDistribution() const { return eTarget_; }
206 //! Beta scaling to reach target distribution.
207 double targetBetaScaling() const { return targetBetaScaling_; }
208 //! Cutoff for target.
209 double targetCutoff() const { return targetCutoff_; }
210 //! Which kind of growth to use.
211 AwhHistogramGrowthType growthType() const { return eGrowth_; }
212 //! User provided PMF estimate.
213 bool userPMFEstimate() const { return bUserData_; }
214 //! Estimated initial free energy error in kJ/mol.
215 double initialErrorEstimate() const { return errorInitial_; }
216 //! Dimensions of coordinate space.
217 int ndim() const { return dimParams_.size(); }
218 //! Number of groups to share this bias with.
219 int shareGroup() const { return shareGroup_; }
220 //! If the simulation starts with equilibrating histogram.
221 bool equilibrateHistogram() const { return equilibrateHistogram_; }
222 //! Access to dimension parameters.
223 ArrayRef<AwhDimParams> dimParams() { return dimParams_; }
224 //! Const access to dimension parameters.
225 ArrayRef<const AwhDimParams> dimParams() const { return dimParams_; }
226 //! Write datastructure.
227 void serialize(ISerializer* serializer);
230 //! AWH parameters per dimension.
231 std::vector<AwhDimParams> dimParams_;
232 //! Type of target distribution.
233 AwhTargetType eTarget_;
234 //! Beta scaling value for Boltzmann type target distributions.
235 double targetBetaScaling_;
236 //! Free energy cutoff value for cutoff type target distribution in kJ/mol.
237 double targetCutoff_;
238 //! How the biasing histogram grows.
239 AwhHistogramGrowthType eGrowth_;
240 //! Is there a user-defined initial PMF estimate and target estimate?
242 //! Estimated initial free energy error in kJ/mol.
243 double errorInitial_;
244 //! When >0, the bias is shared with biases of the same group and across multiple simulations when shareBiasMultisim=true
246 //! True if the simulation starts out by equilibrating the histogram.
247 bool equilibrateHistogram_;
250 * Structure holding parameter information for AWH.
255 //! Constructor from input file.
256 AwhParams(std::vector<t_inpfile>* inp, const t_inputrec& ir, warninp_t wi);
257 //! Constructor used to generate awh parameter from file reading.
258 explicit AwhParams(ISerializer* serializer);
260 //! Move constructor.
261 AwhParams(AwhParams&&) = default;
262 //! Move assignment operator.
263 AwhParams& operator=(AwhParams&&) = default;
264 //! Delete copy constructor.
265 AwhParams(const AwhParams&) = delete;
266 //! Delete copy assignment.
267 AwhParams& operator=(const AwhParams&) = delete;
269 //! Get number of biases.
270 int numBias() const { return awhBiasParams_.size(); }
271 //! Get access to bias parameters.
272 ArrayRef<AwhBiasParams> awhBiasParams() { return awhBiasParams_; }
273 //! Const access to bias parameters.
274 ArrayRef<const AwhBiasParams> awhBiasParams() const { return awhBiasParams_; }
275 //! What king of potential is being used. \todo should use actual enum class.
276 AwhPotentialType potential() const { return potentialEnum_; }
277 //! Seed used for starting AWH.
278 int64_t seed() const { return seed_; }
279 //! Output step interval.
280 int nstout() const { return nstOut_; }
281 //! Number of samples per coordinate sample.
282 int nstSampleCoord() const { return nstSampleCoord_; }
283 //! Number of samples per free energy update.
284 int numSamplesUpdateFreeEnergy() const { return numSamplesUpdateFreeEnergy_; }
285 //! If biases are shared in multisim.
286 bool shareBiasMultisim() const { return shareBiasMultisim_; }
287 //! Serialize awh parameters.
288 void serialize(ISerializer* serializer);
291 //! AWH bias parameters.
292 std::vector<AwhBiasParams> awhBiasParams_;
295 //! Output step interval.
297 //! Number of samples per coordinate sample (also used for PMF)
299 //! Number of samples per free energy update.
300 int numSamplesUpdateFreeEnergy_;
301 //! Type of potential.
302 AwhPotentialType potentialEnum_;
303 //! Whether to share biases with shareGroup>0 between multi-simulations.
304 bool shareBiasMultisim_;
309 #endif /* GMX_MDTYPES_AWH_PARAMS_H */