2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017, 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/mdtypes/awh-params.h"
52 #include "gromacs/mdtypes/commrec.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/stringutil.h"
59 bool haveBiasSharingWithinSimulation(const AwhParams &awhParams)
61 bool haveSharing = false;
63 for (int k = 0; k < awhParams.numBias; k++)
65 int shareGroup = awhParams.awhBiasParams[k].shareGroup;
68 for (int i = k + 1; i < awhParams.numBias; i++)
70 if (awhParams.awhBiasParams[i].shareGroup == shareGroup)
81 void biasesAreCompatibleForSharingBetweenSimulations(const AwhParams &awhParams,
82 const std::vector<size_t> &pointSize,
83 const gmx_multisim_t *multiSimComm)
85 const int numSim = multiSimComm->nsim;
87 /* We currently enforce subsequent shared biases to have consecutive
88 * share-group values starting at 1. This means we can reduce shared
89 * biases in order over the ranks and it does not restrict possibilities.
92 for (int b = 0; b < awhParams.numBias; b++)
94 int group = awhParams.awhBiasParams[b].shareGroup;
98 if (group != numShare)
100 GMX_THROW(InvalidInputError("AWH biases that are shared should use consequetive share-group values starting at 1"));
104 std::vector<int> numShareAll(numSim);
105 numShareAll[multiSimComm->sim] = numShare;
106 gmx_sumi_sim(numShareAll.size(), numShareAll.data(), multiSimComm);
107 for (int sim = 1; sim < numSim; sim++)
109 if (numShareAll[sim] != numShareAll[0])
111 GMX_THROW(InvalidInputError("Different simulations attempt to share different number of biases"));
115 std::vector<int> intervals(numSim*2);
116 intervals[numSim*0 + multiSimComm->sim] = awhParams.nstSampleCoord;
117 intervals[numSim*1 + multiSimComm->sim] = awhParams.numSamplesUpdateFreeEnergy;
118 gmx_sumi_sim(intervals.size(), intervals.data(), multiSimComm);
119 for (int sim = 1; sim < numSim; sim++)
121 if (intervals[sim] != intervals[0])
123 GMX_THROW(InvalidInputError("All simulations should have the same AWH sample interval"));
125 if (intervals[numSim + sim] != intervals[numSim])
127 GMX_THROW(InvalidInputError("All simulations should have the same AWH free-energy update interval"));
131 /* Check the point sizes. This is a sufficient condition for running
132 * as shared multi-sim run. No physics checks are performed here.
134 for (int b = 0; b < awhParams.numBias; b++)
136 if (awhParams.awhBiasParams[b].shareGroup > 0)
138 std::vector<gmx_int64_t> pointSizes(numSim);
139 pointSizes[multiSimComm->sim] = pointSize[b];
140 gmx_sumli_sim(pointSizes.size(), pointSizes.data(), multiSimComm);
141 for (int sim = 1; sim < numSim; sim++)
143 if (pointSizes[sim] != pointSizes[0])
145 GMX_THROW(InvalidInputError(gmx::formatString("Shared AWH bias %d has different grid sizes in different simulations\n", b + 1)));