77f23bbc956af3d386e33bf3c7ee2708c0736f49
[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, 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 /*! \libinternal \file
37  *
38  * \brief
39  * Declares AWH parameter data types.
40  *
41  * Besides internal use by the AWH module, the AWH parameters are needed
42  * for reading the user input (mdp) file and for reading and writing the
43  * parameters to the mdrun input (tpr) file.
44  *
45  * \author Viveca Lindahl
46  * \inlibraryapi
47  * \ingroup module_mdtypes
48  */
49
50 #ifndef GMX_MDTYPES_AWH_PARAMS_H
51 #define GMX_MDTYPES_AWH_PARAMS_H
52
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/utility/basedefinitions.h"
55
56 namespace gmx
57 {
58
59 //! Target distribution enum.
60 enum {
61     eawhtargetCONSTANT, eawhtargetCUTOFF, eawhtargetBOLTZMANN, eawhtargetLOCALBOLTZMANN, eawhtargetNR
62 };
63 //! String for target distribution.
64 extern const char *eawhtarget_names[eawhtargetNR+1];
65 //! Macro for target distribution string.
66 #define EAWHTARGET(e)  enum_name(e, gmx::eawhtargetNR, gmx::eawhtarget_names)
67
68 //! Weight histogram growth enum.
69 enum {
70     eawhgrowthEXP_LINEAR, eawhgrowthLINEAR, eawhgrowthNR
71 };
72 //! String for weight histogram growth
73 extern const char *eawhgrowth_names[eawhgrowthNR+1];
74 //! Macro for weight histogram growth string.
75 #define EAWHGROWTH(e)    enum_name(e, gmx::eawhgrowthNR, gmx::eawhgrowth_names)
76
77 //! AWH potential type enum.
78 enum {
79     eawhpotentialCONVOLVED, eawhpotentialUMBRELLA, eawhpotentialNR
80 };
81 //! String for AWH potential type
82 extern const char *eawhpotential_names[eawhpotentialNR+1];
83 //! Macro for AWH potential type string.
84 #define EAWHPOTENTIAL(e)    enum_name(e, gmx::eawhpotentialNR, gmx::eawhpotential_names)
85
86 //! AWH bias reaction coordinate provider
87 enum {
88     eawhcoordproviderPULL, eawhcoordproviderNR
89 };
90 //! String for AWH bias reaction coordinate provider.
91 extern const char *eawhcoordprovider_names[eawhcoordproviderNR+1];
92 //! Macro for AWH bias reaction coordinate provider.
93 #define EAWHCOORDPROVIDER(e)    enum_name(e, gmx::eawhcoordproviderNR, gmx::eawhcoordprovider_names)
94
95 /*! \cond INTERNAL */
96
97 //! Parameters for an AWH coordinate dimension.
98 struct AwhDimParams
99 {
100     int    eCoordProvider; /**< The module providing the reaction coordinate. */
101     int    coordIndex;     /**< Index of reaction coordinate in the provider. */
102     double origin;         /**< Start value of the interval. */
103     double end;            /**< End value of the interval. */
104     double period;         /**< The period of this dimension (= 0 if not periodic). */
105     double forceConstant;  /**< The force constant in kJ/mol/nm^2, kJ/mol/rad^2 */
106     double diffusion;      /**< Estimated diffusion constant in units of nm^2/ps or rad^2/ps. */
107     double coordValueInit; /**< The initial coordinate value. */
108     double coverDiameter;  /**< The diameter that needs to be sampled around a point before it is considered covered. */
109 };
110
111 //! Parameters for an AWH bias.
112 struct AwhBiasParams
113 {
114     // TODO: Turn dimParams into a std::vector when moved into AWH module
115     int           ndim;                  /**< Dimension of the coordinate space. */
116     AwhDimParams *dimParams;             /**< AWH parameters per dimension. */
117     int           eTarget;               /**< Type of target distribution. */
118     double        targetBetaScaling;     /**< Beta scaling value for Boltzmann type target distributions. */
119     double        targetCutoff;          /**< Free energy cutoff value for cutoff type target distribution in kJ/mol.*/
120     int           eGrowth;               /**< How the biasing histogram grows. */
121     int           bUserData;             /**< Is there a user-defined initial PMF estimate and target estimate? */
122     double        errorInitial;          /**< Estimated initial free energy error in kJ/mol. */
123     int           shareGroup;            /**< When >0, the bias is shared with biases of the same group and across multiple simulations when shareBiasMultisim=true */
124     gmx_bool      equilibrateHistogram;  /**< True if the simulation starts out by equilibrating the histogram.  */
125 };
126
127 //! Parameters for AWH.
128 struct AwhParams
129 {
130     // TODO: Turn awhBiasParams into a std::vector when moved into AWH module
131     int            numBias;                    /**< The number of AWH biases.*/
132     AwhBiasParams *awhBiasParams;              /**< AWH bias parameters.*/
133     gmx_int64_t    seed;                       /**< Random seed.*/
134     int            nstOut;                     /**< Output step interval.*/
135     int            nstSampleCoord;             /**< Number of samples per coordinate sample (also used for PMF) */
136     int            numSamplesUpdateFreeEnergy; /**< Number of samples per free energy update. */
137     int            ePotential;                 /**< Type of potential. */
138     gmx_bool       shareBiasMultisim;          /**< When true, share biases with shareGroup>0 between multi-simulations */
139 };
140
141 /*! \endcond */
142
143 }      // namespace gmx
144
145 #endif /* GMX_MDTYPES_AWH_PARAMS_H */