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