Fix output string for GPU support
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / biasparams.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,2019,2020, 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 /*! \internal \file
37  * \brief
38  * Implements the initialization of the BiasParams class.
39  *
40  * \author Viveca Lindahl
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_awh
43  */
44
45 #include "gmxpre.h"
46
47 #include "biasparams.h"
48
49 #include <cmath>
50
51 #include <algorithm>
52
53 #include "gromacs/math/functions.h"
54 #include "gromacs/mdtypes/awh_params.h"
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/gmxassert.h"
58
59 #include "biasgrid.h"
60
61 namespace gmx
62 {
63
64 namespace
65 {
66
67 /*! \brief Determines the interval for updating the target distribution.
68  *
69  * The interval value is based on the target distrbution type
70  * (this could be made a user-option but there is most likely
71  * no big need for tweaking this for most users).
72  *
73  * \param[in] awhParams      AWH parameters.
74  * \param[in] awhBiasParams  Bias parameters.
75  * \returns the target update interval in steps.
76  */
77 int64_t calcTargetUpdateInterval(const AwhParams& awhParams, const AwhBiasParams& awhBiasParams)
78 {
79     int64_t numStepsUpdateTarget = 0;
80     /* Set the target update frequency based on the target distrbution type
81      * (this could be made a user-option but there is most likely no big need
82      * for tweaking this for most users).
83      */
84     switch (awhBiasParams.eTarget)
85     {
86         case eawhtargetCONSTANT: numStepsUpdateTarget = 0; break;
87         case eawhtargetCUTOFF:
88         case eawhtargetBOLTZMANN:
89             /* Updating the target generally requires updating the whole grid so to keep the cost
90                down we generally update the target less often than the free energy (unless the free
91                energy update step is set to > 100 samples). */
92             numStepsUpdateTarget = std::max(100 % awhParams.numSamplesUpdateFreeEnergy,
93                                             awhParams.numSamplesUpdateFreeEnergy)
94                                    * awhParams.nstSampleCoord;
95             break;
96         case eawhtargetLOCALBOLTZMANN:
97             /* The target distribution is set equal to the reference histogram which is updated every free energy update.
98                So the target has to be updated at the same time. This leads to a global update each time because it is
99                assumed that a target distribution update can take any form. This is a bit unfortunate for a "local"
100                target distribution. One could avoid the global update by making a local target update function (and
101                postponing target updates for non-local points as for the free energy update). We avoid such additions
102                for now and accept that this target type always does global updates. */
103             numStepsUpdateTarget = awhParams.numSamplesUpdateFreeEnergy * awhParams.nstSampleCoord;
104             break;
105         default: GMX_RELEASE_ASSERT(false, "Unknown AWH target type"); break;
106     }
107
108     return numStepsUpdateTarget;
109 }
110
111 /*! \brief Determines the step interval for checking for covering.
112  *
113  * \param[in] awhParams      AWH parameters.
114  * \param[in] dimParams         Parameters for the dimensions of the coordinate.
115  * \param[in] gridAxis          The BiasGrid axes.
116  * \returns the check interval in steps.
117  */
118 int64_t calcCheckCoveringInterval(const AwhParams&              awhParams,
119                                   const std::vector<DimParams>& dimParams,
120                                   const std::vector<GridAxis>&  gridAxis)
121 {
122     /* Each sample will have a width of sigma. To cover the axis a
123        minimum number of samples of width sigma is required. */
124     int minNumSamplesCover = 0;
125     for (size_t d = 0; d < gridAxis.size(); d++)
126     {
127         int numSamplesCover;
128         if (!dimParams[d].isFepLambdaDimension())
129         {
130             GMX_RELEASE_ASSERT(
131                     dimParams[d].betak > 0,
132                     "Inverse temperature (beta) and force constant (k) should be positive.");
133             double sigma = 1.0 / std::sqrt(dimParams[d].betak);
134
135             /* The additional sample is here because to cover a discretized
136             axis of length sigma one needs two samples, one for each
137             end point. */
138             GMX_RELEASE_ASSERT(gridAxis[d].length() / sigma < std::numeric_limits<int>::max(),
139                                "The axis length in units of sigma should fit in an int");
140             numSamplesCover = static_cast<int>(std::ceil(gridAxis[d].length() / sigma)) + 1;
141         }
142         else
143         {
144             numSamplesCover = dimParams[d].numFepLambdaStates;
145         }
146         /* The minimum number of samples needed for simultaneously
147         covering all axes is limited by the axis requiring most
148         samples. */
149         minNumSamplesCover = std::max(minNumSamplesCover, numSamplesCover);
150     }
151
152     /* Convert to number of steps using the sampling frequency. The
153        check interval should be a multiple of the update step
154        interval. */
155     int numStepsUpdate = awhParams.numSamplesUpdateFreeEnergy * awhParams.nstSampleCoord;
156     GMX_RELEASE_ASSERT(awhParams.numSamplesUpdateFreeEnergy > 0,
157                        "When checking for AWH coverings, the number of samples per AWH update need "
158                        "to be > 0.");
159     int numUpdatesCheck = std::max(1, minNumSamplesCover / awhParams.numSamplesUpdateFreeEnergy);
160     int numStepsCheck   = numUpdatesCheck * numStepsUpdate;
161
162     GMX_RELEASE_ASSERT(numStepsCheck % numStepsUpdate == 0,
163                        "Only check covering at free energy update steps");
164
165     return numStepsCheck;
166 }
167
168 /*! \brief
169  * Returns an approximation of the geometry factor used for initializing the AWH update size.
170  *
171  * The geometry factor is defined as the following sum of Gaussians:
172  * sum_{k!=0} exp(-0.5*(k*pi*x)^2)/(pi*k)^2,
173  * where k is a xArray.size()-dimensional integer vector with k_i in {0,1,..}.
174  *
175  * \param[in] xArray  Array to evaluate.
176  * \returns the geometry factor.
177  */
178 double gaussianGeometryFactor(gmx::ArrayRef<const double> xArray)
179 {
180     /* For convenience we give the geometry factor function a name: zeta(x) */
181     constexpr size_t                    tableSize  = 5;
182     std::array<const double, tableSize> xTabulated = { { 1e-5, 1e-4, 1e-3, 1e-2, 1e-1 } };
183     std::array<const double, tableSize> zetaTable1d = { { 0.166536811948, 0.16653116886, 0.166250075882,
184                                                           0.162701098306, 0.129272430287 } };
185     std::array<const double, tableSize> zetaTable2d = { { 2.31985974274, 1.86307292523, 1.38159772648,
186                                                           0.897554759158, 0.405578211115 } };
187
188     gmx::ArrayRef<const double> zetaTable;
189
190     if (xArray.size() == 1)
191     {
192         zetaTable = zetaTable1d;
193     }
194     else if (xArray.size() == 2)
195     { // NOLINT bugprone-branch-clone
196         zetaTable = zetaTable2d;
197     }
198     else
199     {
200         /* TODO... but this is anyway a rough estimate and > 2 dimensions is not so popular.
201          * Remove the above NOLINT when addressing this */
202         zetaTable = zetaTable2d;
203     }
204
205     /* TODO. Really zeta is a function of an ndim-dimensional vector x and we shoudl have a ndim-dimensional lookup-table.
206        Here we take the geometric average of the components of x which is ok if the x-components are not very different. */
207     double xScalar = 1;
208     for (const double& x : xArray)
209     {
210         xScalar *= x;
211     }
212
213     GMX_ASSERT(!xArray.empty(), "We should have a non-empty input array");
214     xScalar = std::pow(xScalar, 1.0 / xArray.size());
215
216     /* Look up zeta(x) */
217     size_t xIndex = 0;
218     while ((xIndex < xTabulated.size()) && (xScalar > xTabulated[xIndex]))
219     {
220         xIndex++;
221     }
222
223     double zEstimate;
224     if (xIndex == xTabulated.size())
225     {
226         /* Take last value */
227         zEstimate = zetaTable[xTabulated.size() - 1];
228     }
229     else if (xIndex == 0)
230     {
231         zEstimate = zetaTable[xIndex];
232     }
233     else
234     {
235         /* Interpolate */
236         double x0 = xTabulated[xIndex - 1];
237         double x1 = xTabulated[xIndex];
238         double w  = (xScalar - x0) / (x1 - x0);
239         zEstimate = w * zetaTable[xIndex - 1] + (1 - w) * zetaTable[xIndex];
240     }
241
242     return zEstimate;
243 }
244
245 /*! \brief
246  * Estimate a reasonable initial reference weight histogram size.
247  *
248  * \param[in] dimParams         Parameters for the dimensions of the coordinate.
249  * \param[in] awhBiasParams     Bias parameters.
250  * \param[in] gridAxis          The BiasGrid axes.
251  * \param[in] beta              1/(k_B T).
252  * \param[in] samplingTimestep  Sampling frequency of probability weights.
253  * \returns estimate of initial histogram size.
254  */
255 double getInitialHistogramSizeEstimate(const std::vector<DimParams>& dimParams,
256                                        const AwhBiasParams&          awhBiasParams,
257                                        const std::vector<GridAxis>&  gridAxis,
258                                        double                        beta,
259                                        double                        samplingTimestep)
260 {
261     /* Get diffusion factor */
262     double              crossingTime = 0.;
263     std::vector<double> x;
264     for (size_t d = 0; d < gridAxis.size(); d++)
265     {
266         double axisLength = gridAxis[d].isFepLambdaAxis() ? 1.0 : gridAxis[d].length();
267         if (axisLength > 0)
268         {
269             crossingTime += awhBiasParams.dimParams[d].diffusion / (axisLength * axisLength);
270             /* The sigma of the Gaussian distribution in the umbrella */
271             double sigma = 1.;
272             if (dimParams[d].betak != 0)
273             {
274                 sigma /= std::sqrt(dimParams[d].betak);
275             }
276             x.push_back(sigma / axisLength);
277         }
278     }
279     GMX_RELEASE_ASSERT(crossingTime > 0, "We need at least one dimension with non-zero length");
280     double errorInitialInKT = beta * awhBiasParams.errorInitial;
281     double histogramSize    = gaussianGeometryFactor(x)
282                            / (crossingTime * gmx::square(errorInitialInKT) * samplingTimestep);
283
284     return histogramSize;
285 }
286
287 /*! \brief Returns the number of simulations sharing bias updates.
288  *
289  * \param[in] awhBiasParams          Bias parameters.
290  * \param[in] numSharingSimulations  The number of simulations to share the bias across.
291  * \returns the number of shared updates.
292  */
293 int getNumSharedUpdate(const AwhBiasParams& awhBiasParams, int numSharingSimulations)
294 {
295     GMX_RELEASE_ASSERT(numSharingSimulations >= 1, "We should ''share'' at least with ourselves");
296
297     int numShared = 1;
298
299     if (awhBiasParams.shareGroup > 0)
300     {
301         /* We do not yet support sharing within a simulation */
302         int numSharedWithinThisSimulation = 1;
303         numShared                         = numSharingSimulations * numSharedWithinThisSimulation;
304     }
305
306     return numShared;
307 }
308
309 } // namespace
310
311 BiasParams::BiasParams(const AwhParams&              awhParams,
312                        const AwhBiasParams&          awhBiasParams,
313                        const std::vector<DimParams>& dimParams,
314                        double                        beta,
315                        double                        mdTimeStep,
316                        DisableUpdateSkips            disableUpdateSkips,
317                        int                           numSharingSimulations,
318                        const std::vector<GridAxis>&  gridAxis,
319                        int                           biasIndex) :
320     invBeta(beta > 0 ? 1 / beta : 0),
321     numStepsSampleCoord_(awhParams.nstSampleCoord),
322     numSamplesUpdateFreeEnergy_(awhParams.numSamplesUpdateFreeEnergy),
323     numStepsUpdateTarget_(calcTargetUpdateInterval(awhParams, awhBiasParams)),
324     numStepsCheckCovering_(calcCheckCoveringInterval(awhParams, dimParams, gridAxis)),
325     eTarget(awhBiasParams.eTarget),
326     freeEnergyCutoffInKT(beta * awhBiasParams.targetCutoff),
327     temperatureScaleFactor(awhBiasParams.targetBetaScaling),
328     idealWeighthistUpdate(eTarget != eawhtargetLOCALBOLTZMANN),
329     numSharedUpdate(getNumSharedUpdate(awhBiasParams, numSharingSimulations)),
330     updateWeight(numSamplesUpdateFreeEnergy_ * numSharedUpdate),
331     localWeightScaling(eTarget == eawhtargetLOCALBOLTZMANN ? temperatureScaleFactor : 1),
332     initialErrorInKT(beta * awhBiasParams.errorInitial),
333     initialHistogramSize(getInitialHistogramSizeEstimate(dimParams,
334                                                          awhBiasParams,
335                                                          gridAxis,
336                                                          beta,
337                                                          numStepsSampleCoord_ * mdTimeStep)),
338     convolveForce(awhParams.ePotential == eawhpotentialCONVOLVED),
339     biasIndex(biasIndex),
340     disableUpdateSkips_(disableUpdateSkips == DisableUpdateSkips::yes)
341 {
342     if (beta <= 0)
343     {
344         GMX_THROW(InvalidInputError("To use AWH, the beta=1/(k_B T) should be > 0"));
345     }
346
347     for (int d = 0; d < awhBiasParams.ndim; d++)
348     {
349         /* The spacing in FEP dimensions is 1. The calculated coverRadius will be in lambda states
350          * (cf points in other dimensions). */
351         double coverRadiusInNm = 0.5 * awhBiasParams.dimParams[d].coverDiameter;
352         double spacing         = gridAxis[d].spacing();
353         coverRadius_[d] = spacing > 0 ? static_cast<int>(std::round(coverRadiusInNm / spacing)) : 0;
354     }
355 }
356
357 } // namespace gmx