debbb846728fff2816f39523c2fe3e6b1b26326c
[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 "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/utility/basedefinitions.h"
56
57 namespace gmx
58 {
59
60 //! Target distribution enum.
61 enum class AwhTargetType : int
62 {
63     Constant,
64     Cutoff,
65     Boltzmann,
66     LocalBoltzmann,
67     Count,
68     Default = Constant
69 };
70 //! String for target distribution.
71 const char* enumValueToString(AwhTargetType enumValue);
72
73 //! Weight histogram growth enum.
74 enum class AwhHistogramGrowthType : int
75 {
76     ExponentialLinear,
77     Linear,
78     Count,
79     Default = ExponentialLinear
80 };
81 //! String for weight histogram growth
82 const char* enumValueToString(AwhHistogramGrowthType enumValue);
83
84 //! AWH potential type enum.
85 enum class AwhPotentialType : int
86 {
87     Convolved,
88     Umbrella,
89     Count,
90     Default = Convolved
91 };
92 //! String for AWH potential type
93 const char* enumValueToString(AwhPotentialType enumValue);
94
95 //! AWH bias reaction coordinate provider
96 enum class AwhCoordinateProviderType : int
97 {
98     Pull,
99     FreeEnergyLambda,
100     Count,
101     Default = Pull
102 };
103 //! String for AWH bias reaction coordinate provider.
104 const char* enumValueToString(AwhCoordinateProviderType enumValue);
105
106 /*! \cond INTERNAL */
107
108 //! Parameters for an AWH coordinate dimension.
109 struct AwhDimParams
110 {
111     AwhCoordinateProviderType eCoordProvider; /**< The module providing the reaction coordinate. */
112     int                       coordIndex;     /**< Index of reaction coordinate in the provider. */
113     double                    origin;         /**< Start value of the interval. */
114     double                    end;            /**< End value of the interval. */
115     double                    period; /**< The period of this dimension (= 0 if not periodic). */
116     double                    forceConstant; /**< The force constant in kJ/mol/nm^2, kJ/mol/rad^2 */
117     double diffusion; /**< Estimated diffusion constant in units of nm^2/ps, rad^2/ps or ps^-1. */
118     double coordValueInit; /**< The initial coordinate value. */
119     double coverDiameter; /**< The diameter that needs to be sampled around a point before it is considered covered. */
120 };
121
122 //! Parameters for an AWH bias.
123 struct AwhBiasParams
124 {
125     // TODO: Turn dimParams into a std::vector when moved into AWH module
126     int           ndim;       /**< Dimension of the coordinate space. */
127     AwhDimParams* dimParams;  /**< AWH parameters per dimension. */
128     AwhTargetType eTarget;    /**< Type of target distribution. */
129     double targetBetaScaling; /**< Beta scaling value for Boltzmann type target distributions. */
130     double targetCutoff; /**< Free energy cutoff value for cutoff type target distribution in kJ/mol.*/
131     AwhHistogramGrowthType eGrowth; /**< How the biasing histogram grows. */
132     bool     bUserData;    /**< Is there a user-defined initial PMF estimate and target estimate? */
133     double   errorInitial; /**< Estimated initial free energy error in kJ/mol. */
134     int      shareGroup; /**< When >0, the bias is shared with biases of the same group and across multiple simulations when shareBiasMultisim=true */
135     gmx_bool equilibrateHistogram; /**< True if the simulation starts out by equilibrating the histogram. */
136 };
137
138 //! Parameters for AWH.
139 struct AwhParams
140 {
141     // TODO: Turn awhBiasParams into a std::vector when moved into AWH module
142     int            numBias;       /**< The number of AWH biases.*/
143     AwhBiasParams* awhBiasParams; /**< AWH bias parameters.*/
144     int64_t        seed;          /**< Random seed.*/
145     int            nstOut;        /**< Output step interval.*/
146     int nstSampleCoord; /**< Number of samples per coordinate sample (also used for PMF) */
147     int numSamplesUpdateFreeEnergy; /**< Number of samples per free energy update. */
148     AwhPotentialType ePotential;    /**< Type of potential. */
149     gmx_bool         shareBiasMultisim; /**< When true, share biases with shareGroup>0 between multi-simulations */
150 };
151
152 /*! \endcond */
153
154 } // namespace gmx
155
156 #endif /* GMX_MDTYPES_AWH_PARAMS_H */