2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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.
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.
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.
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.
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.
38 * Implements bias sharing checking functionality.
40 * \author Berk Hess <hess@kth.se>
46 #include "biassharing.h"
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/exceptions.h"
55 #include "gromacs/utility/stringutil.h"
60 bool haveBiasSharingWithinSimulation(const AwhParams& awhParams)
62 bool haveSharing = false;
64 for (int k = 0; k < awhParams.numBias; k++)
66 int shareGroup = awhParams.awhBiasParams[k].shareGroup;
69 for (int i = k + 1; i < awhParams.numBias; i++)
71 if (awhParams.awhBiasParams[i].shareGroup == shareGroup)
82 void biasesAreCompatibleForSharingBetweenSimulations(const AwhParams& awhParams,
83 const std::vector<size_t>& pointSize,
84 const gmx_multisim_t* multiSimComm)
86 const int numSim = multiSimComm->numSimulations_;
88 /* We currently enforce subsequent shared biases to have consecutive
89 * share-group values starting at 1. This means we can reduce shared
90 * biases in order over the ranks and it does not restrict possibilities.
93 for (int b = 0; b < awhParams.numBias; b++)
95 int group = awhParams.awhBiasParams[b].shareGroup;
99 if (group != numShare)
102 InvalidInputError("AWH biases that are shared should use consequetive "
103 "share-group values starting at 1"));
107 std::vector<int> numShareAll(numSim);
108 numShareAll[multiSimComm->simulationIndex_] = numShare;
109 gmx_sumi_sim(numShareAll.size(), numShareAll.data(), multiSimComm);
110 for (int sim = 1; sim < numSim; sim++)
112 if (numShareAll[sim] != numShareAll[0])
114 GMX_THROW(InvalidInputError(
115 "Different simulations attempt to share different number of biases"));
119 std::vector<int> intervals(numSim * 2);
120 intervals[numSim * 0 + multiSimComm->simulationIndex_] = awhParams.nstSampleCoord;
121 intervals[numSim * 1 + multiSimComm->simulationIndex_] = awhParams.numSamplesUpdateFreeEnergy;
122 gmx_sumi_sim(intervals.size(), intervals.data(), multiSimComm);
123 for (int sim = 1; sim < numSim; sim++)
125 if (intervals[sim] != intervals[0])
128 InvalidInputError("All simulations should have the same AWH sample interval"));
130 if (intervals[numSim + sim] != intervals[numSim])
132 GMX_THROW(InvalidInputError(
133 "All simulations should have the same AWH free-energy update interval"));
137 /* Check the point sizes. This is a sufficient condition for running
138 * as shared multi-sim run. No physics checks are performed here.
140 for (int b = 0; b < awhParams.numBias; b++)
142 if (awhParams.awhBiasParams[b].shareGroup > 0)
144 std::vector<int64_t> pointSizes(numSim);
145 pointSizes[multiSimComm->simulationIndex_] = pointSize[b];
146 gmx_sumli_sim(pointSizes.size(), pointSizes.data(), multiSimComm);
147 for (int sim = 1; sim < numSim; sim++)
149 if (pointSizes[sim] != pointSizes[0])
151 GMX_THROW(InvalidInputError(
152 gmx::formatString("Shared AWH bias %d has different grid sizes in "
153 "different simulations\n",