2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, 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.
35 /*! \libinternal \file
37 * \brief Declares the multi-simulation support routines.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_mdrunutility
43 #ifndef GMX_MDRUNUTILITY_MULTISIM_H
44 #define GMX_MDRUNUTILITY_MULTISIM_H
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/gmxmpi.h"
51 #include "gromacs/utility/mpiinplacebuffers.h"
54 * \brief Coordinate multi-simulation resources for mdrun
56 * \todo Change this to class
60 //! Default constructor
62 /*! \brief Constructor useful for mdrun simulations
64 * Splits the communicator into multidirs.size() separate
65 * simulations, if >1, and creates a communication structure
66 * between the master these simulations.
68 * Valid to call regardless of build configuration, but \c
69 * multidirs must be empty unless a real MPI build is used. */
70 gmx_multisim_t(MPI_Comm comm
,
71 gmx::ArrayRef
<const std::string
> multidirs
);
75 //! The number of simulations in the set of multi-simulations
77 //! The index of the simulation that owns this object within the set
79 //! The MPI Group between master ranks of simulations, valid only on master ranks.
80 MPI_Group mpi_group_masters
= MPI_GROUP_NULL
;
81 //! The MPI communicator between master ranks of simulations, valid only on master ranks.
82 MPI_Comm mpi_comm_masters
= MPI_COMM_NULL
;
83 //! Communication buffers needed if MPI_IN_PLACE isn't supported
84 mpi_in_place_buf_t
* mpb
= nullptr;
87 //! Calculate the sum over the simulations of an array of ints
88 void gmx_sumi_sim(int nr
, int r
[], const gmx_multisim_t
*ms
);
90 //! Calculate the sum over the simulations of an array of large ints
91 void gmx_sumli_sim(int nr
, int64_t r
[], const gmx_multisim_t
*ms
);
93 //! Calculate the sum over the simulations of an array of floats
94 void gmx_sumf_sim(int nr
, float r
[], const gmx_multisim_t
*ms
);
96 //! Calculate the sum over the simulations of an array of doubles
97 void gmx_sumd_sim(int nr
, double r
[], const gmx_multisim_t
*ms
);
99 /*! \brief Return a vector containing the gathered values of \c
100 * localValue found on the master rank of each simulation. */
101 std::vector
<int> gatherIntFromMultiSimulation(const gmx_multisim_t
*ms
,
104 /*! \brief Check if val is the same on all simulations for a mdrun
107 * The string name is used to print to the log file and in a fatal error
108 * if the val's don't match. If bQuiet is true and the check passes,
109 * no output is written. */
110 void check_multi_int(FILE *log
, const gmx_multisim_t
*ms
,
111 int val
, const char *name
,
113 /*! \copydoc check_multi_int() */
114 void check_multi_int64(FILE *log
, const gmx_multisim_t
*ms
,
115 int64_t val
, const char *name
,
119 //! Convenience define for sum of reals
120 #define gmx_sum_sim gmx_sumd_sim
122 //! Convenience define for sum of reals
123 #define gmx_sum_sim gmx_sumf_sim
126 //! Are we doing multiple independent simulations?
127 static bool inline isMultiSim(const gmx_multisim_t
*ms
)
129 return ms
!= nullptr;
132 /*! \brief Return whether this rank is the master rank of a
133 * simulation, using \c ms (if it is valid) and otherwise \c
135 bool findIsSimulationMasterRank(const gmx_multisim_t
*ms
,
136 MPI_Comm communicator
);
138 //! Are we the master simulation of a possible multi-simulation?
139 bool isMasterSim(const gmx_multisim_t
*ms
);
141 /*! \brief Are we the master rank (of the master simulation, for a multi-sim).
143 * This rank prints the remaining run time etc. */
144 bool isMasterSimMasterRank(const gmx_multisim_t
*ms
,
147 //! Make a barrier across all multi-simulation master ranks
148 void multiSimBarrier(const gmx_multisim_t
*ms
);