50d2603778fe3f21325d0d4869202d9023abd9da
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / biassharing.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020,2021, 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 bias sharing checking functionality.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_awh
42  */
43
44 #include "gmxpre.h"
45
46 #include "biassharing.h"
47
48 #include <vector>
49
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/mdrunutility/multisim.h"
52 #include "gromacs/mdtypes/awh_params.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/stringutil.h"
57
58 namespace gmx
59 {
60
61 bool haveBiasSharingWithinSimulation(const AwhParams& awhParams)
62 {
63     bool haveSharing = false;
64
65     const auto awhBiasParams = awhParams.awhBiasParams();
66     for (auto awhBiasParamIt = awhBiasParams.begin(); awhBiasParamIt != awhBiasParams.end(); ++awhBiasParamIt)
67     {
68         int shareGroup = awhBiasParamIt->shareGroup();
69         if (shareGroup > 0)
70         {
71             for (auto awhBiasParamIt2 = awhBiasParamIt + 1; awhBiasParamIt2 != awhBiasParams.end();
72                  ++awhBiasParamIt2)
73             {
74                 if (awhBiasParamIt2->shareGroup() == shareGroup)
75                 {
76                     haveSharing = true;
77                 }
78             }
79         }
80     }
81
82     return haveSharing;
83 }
84
85 void biasesAreCompatibleForSharingBetweenSimulations(const AwhParams&       awhParams,
86                                                      ArrayRef<const size_t> pointSize,
87                                                      const gmx_multisim_t*  multiSimComm)
88 {
89     const int numSim = multiSimComm->numSimulations_;
90
91     /* We currently enforce subsequent shared biases to have consecutive
92      * share-group values starting at 1. This means we can reduce shared
93      * biases in order over the ranks and it does not restrict possibilities.
94      */
95     int numShare = 0;
96     for (const auto& awhBiasParam : awhParams.awhBiasParams())
97     {
98         int group = awhBiasParam.shareGroup();
99         if (group > 0)
100         {
101             numShare++;
102             if (group != numShare)
103             {
104                 GMX_THROW(
105                         InvalidInputError("AWH biases that are shared should use consequetive "
106                                           "share-group values starting at 1"));
107             }
108         }
109     }
110     std::vector<int> numShareAll(numSim);
111     numShareAll[multiSimComm->simulationIndex_] = numShare;
112     gmx_sumi_sim(numShareAll.size(), numShareAll.data(), multiSimComm);
113     for (int sim = 1; sim < numSim; sim++)
114     {
115         if (numShareAll[sim] != numShareAll[0])
116         {
117             GMX_THROW(InvalidInputError(
118                     "Different simulations attempt to share different number of biases"));
119         }
120     }
121
122     std::vector<int> intervals(numSim * 2);
123     intervals[numSim * 0 + multiSimComm->simulationIndex_] = awhParams.nstSampleCoord();
124     intervals[numSim * 1 + multiSimComm->simulationIndex_] = awhParams.numSamplesUpdateFreeEnergy();
125     gmx_sumi_sim(intervals.size(), intervals.data(), multiSimComm);
126     for (int sim = 1; sim < numSim; sim++)
127     {
128         if (intervals[sim] != intervals[0])
129         {
130             GMX_THROW(
131                     InvalidInputError("All simulations should have the same AWH sample interval"));
132         }
133         if (intervals[numSim + sim] != intervals[numSim])
134         {
135             GMX_THROW(InvalidInputError(
136                     "All simulations should have the same AWH free-energy update interval"));
137         }
138     }
139
140     /* Check the point sizes. This is a sufficient condition for running
141      * as shared multi-sim run. No physics checks are performed here.
142      */
143     const auto& awhBiasParams = awhParams.awhBiasParams();
144     for (int b = 0; b < gmx::ssize(awhBiasParams); b++)
145     {
146         if (awhBiasParams[b].shareGroup() > 0)
147         {
148             std::vector<int64_t> pointSizes(numSim);
149             pointSizes[multiSimComm->simulationIndex_] = pointSize[b];
150             gmx_sumli_sim(pointSizes.size(), pointSizes.data(), multiSimComm);
151             for (int sim = 1; sim < numSim; sim++)
152             {
153                 if (pointSizes[sim] != pointSizes[0])
154                 {
155                     GMX_THROW(InvalidInputError(
156                             gmx::formatString("Shared AWH bias %d has different grid sizes in "
157                                               "different simulations\n",
158                                               b + 1)));
159                 }
160             }
161         }
162     }
163 }
164
165 } // namespace gmx