Apply clang-format to source tree
[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,2018,2019, 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 {
62     eawhtargetCONSTANT,
63     eawhtargetCUTOFF,
64     eawhtargetBOLTZMANN,
65     eawhtargetLOCALBOLTZMANN,
66     eawhtargetNR
67 };
68 //! String for target distribution.
69 extern const char* eawhtarget_names[eawhtargetNR + 1];
70 //! Macro for target distribution string.
71 #define EAWHTARGET(e) enum_name(e, gmx::eawhtargetNR, gmx::eawhtarget_names)
72
73 //! Weight histogram growth enum.
74 enum
75 {
76     eawhgrowthEXP_LINEAR,
77     eawhgrowthLINEAR,
78     eawhgrowthNR
79 };
80 //! String for weight histogram growth
81 extern const char* eawhgrowth_names[eawhgrowthNR + 1];
82 //! Macro for weight histogram growth string.
83 #define EAWHGROWTH(e) enum_name(e, gmx::eawhgrowthNR, gmx::eawhgrowth_names)
84
85 //! AWH potential type enum.
86 enum
87 {
88     eawhpotentialCONVOLVED,
89     eawhpotentialUMBRELLA,
90     eawhpotentialNR
91 };
92 //! String for AWH potential type
93 extern const char* eawhpotential_names[eawhpotentialNR + 1];
94 //! Macro for AWH potential type string.
95 #define EAWHPOTENTIAL(e) enum_name(e, gmx::eawhpotentialNR, gmx::eawhpotential_names)
96
97 //! AWH bias reaction coordinate provider
98 enum
99 {
100     eawhcoordproviderPULL,
101     eawhcoordproviderNR
102 };
103 //! String for AWH bias reaction coordinate provider.
104 extern const char* eawhcoordprovider_names[eawhcoordproviderNR + 1];
105 //! Macro for AWH bias reaction coordinate provider.
106 #define EAWHCOORDPROVIDER(e) enum_name(e, gmx::eawhcoordproviderNR, gmx::eawhcoordprovider_names)
107
108 /*! \cond INTERNAL */
109
110 //! Parameters for an AWH coordinate dimension.
111 struct AwhDimParams
112 {
113     int    eCoordProvider; /**< The module providing the reaction coordinate. */
114     int    coordIndex;     /**< Index of reaction coordinate in the provider. */
115     double origin;         /**< Start value of the interval. */
116     double end;            /**< End value of the interval. */
117     double period;         /**< The period of this dimension (= 0 if not periodic). */
118     double forceConstant;  /**< The force constant in kJ/mol/nm^2, kJ/mol/rad^2 */
119     double diffusion;      /**< Estimated diffusion constant in units of nm^2/ps or rad^2/ps. */
120     double coordValueInit; /**< The initial coordinate value. */
121     double coverDiameter; /**< The diameter that needs to be sampled around a point before it is considered covered. */
122 };
123
124 //! Parameters for an AWH bias.
125 struct AwhBiasParams
126 {
127     // TODO: Turn dimParams into a std::vector when moved into AWH module
128     int           ndim;         /**< Dimension of the coordinate space. */
129     AwhDimParams* dimParams;    /**< AWH parameters per dimension. */
130     int           eTarget;      /**< Type of target distribution. */
131     double   targetBetaScaling; /**< Beta scaling value for Boltzmann type target distributions. */
132     double   targetCutoff; /**< Free energy cutoff value for cutoff type target distribution in kJ/mol.*/
133     int      eGrowth;      /**< How the biasing histogram grows. */
134     int      bUserData;    /**< Is there a user-defined initial PMF estimate and target estimate? */
135     double   errorInitial; /**< Estimated initial free energy error in kJ/mol. */
136     int      shareGroup; /**< When >0, the bias is shared with biases of the same group and across multiple simulations when shareBiasMultisim=true */
137     gmx_bool equilibrateHistogram; /**< True if the simulation starts out by equilibrating the histogram. */
138 };
139
140 //! Parameters for AWH.
141 struct AwhParams
142 {
143     // TODO: Turn awhBiasParams into a std::vector when moved into AWH module
144     int            numBias;       /**< The number of AWH biases.*/
145     AwhBiasParams* awhBiasParams; /**< AWH bias parameters.*/
146     int64_t        seed;          /**< Random seed.*/
147     int            nstOut;        /**< Output step interval.*/
148     int      nstSampleCoord; /**< Number of samples per coordinate sample (also used for PMF) */
149     int      numSamplesUpdateFreeEnergy; /**< Number of samples per free energy update. */
150     int      ePotential;                 /**< Type of potential. */
151     gmx_bool shareBiasMultisim; /**< When true, share biases with shareGroup>0 between multi-simulations */
152 };
153
154 /*! \endcond */
155
156 } // namespace gmx
157
158 #endif /* GMX_MDTYPES_AWH_PARAMS_H */