Improve AWH bias sharing documentation
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / biassharing.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,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  *
38  * \brief
39  * Declares functions to check bias sharing properties.
40  *
41  * This actual sharing of biases is currently implemeted in BiasState.
42  *
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_awh
45  */
46
47 #ifndef GMX_AWH_BIASSHARING_H
48 #define GMX_AWH_BIASSHARING_H
49
50 #include <cstddef>
51
52 #include <memory>
53 #include <vector>
54
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/classhelpers.h"
57 #include "gromacs/utility/gmxmpi.h"
58
59 struct gmx_multisim_t;
60 struct t_commrec;
61
62 namespace gmx
63 {
64
65 template<typename>
66 class ArrayRef;
67 class AwhParams;
68
69 class BiasSharing
70 {
71 public:
72     /*! \brief Constructor
73      *
74      * \param[in] awhParams              Parameters for all biases in this simulation
75      * \param[in] commRecord             Intra-simulation communication record
76      * \param[in] simulationMastersComm  MPI communicator for all master ranks of all simulations that share this bias
77      */
78     BiasSharing(const AwhParams& awhParams, const t_commrec& commRecord, MPI_Comm simulationMastersComm);
79
80     ~BiasSharing();
81
82     //! Returns the number of simulations sharing bias \p biasIndex
83     int numSharingSimulations(int biasIndex) const { return numSharingSimulations_[biasIndex]; }
84
85     //! Returns the index of our simulation in the simulations sharing bias \p biasIndex
86     int sharingSimulationIndex(int biasIndex) const { return sharingSimulationIndices_[biasIndex]; }
87
88     //! Sums data of type int over the master ranks of all simulations sharing bias \p biasIndex
89     void sumOverSharingMasterRanks(ArrayRef<int> data, int biasIndex) const;
90
91     //! Sums data of type long over the master ranks of all simulations sharing bias \p biasIndex
92     void sumOverSharingMasterRanks(ArrayRef<long> data, int biasIndex) const;
93
94     /*! \brief Sums data of type int over all simulations sharing bias \p biasIndex
95      *
96      * The summing is performed over the master ranks of the simulations sharing bias \p biasIndex
97      * and the result is broadcasted to the other ranks within each simulation.
98      */
99     void sumOverSharingSimulations(ArrayRef<int> data, int biasIndex) const;
100
101     /*! \brief Sums data of type long over all simulations sharing bias \p biasIndex
102      *
103      * The summing is performed over the master ranks of the simulations sharing bias \p biasIndex
104      * and the result is broadcasted to the other ranks within each simulation.
105      */
106     void sumOverSharingSimulations(ArrayRef<double> data, int biasIndex) const;
107
108 private:
109     //! The number of simulations sharing for each bias
110     std::vector<int> numSharingSimulations_;
111     //! The index of our simulations in the simulations for each bias
112     std::vector<int> sharingSimulationIndices_;
113     //! Reference to the intra-simulation communication record
114     const t_commrec& commRecord_;
115
116     //! Communicator between master ranks sharing a bias, for each bias
117     std::vector<MPI_Comm> multiSimCommPerBias_;
118     //! List of MPI communicators created by this object so we can destroy them on distruction
119     std::vector<MPI_Comm> createdCommList_;
120
121     GMX_DISALLOW_COPY_AND_ASSIGN(BiasSharing);
122 };
123
124 /*! \brief Returns if any bias is sharing within a simulation.
125  *
126  * \param[in] awhParams  The AWH parameters.
127  */
128 bool haveBiasSharingWithinSimulation(const AwhParams& awhParams);
129
130 /*! \brief Checks whether biases are compatible for sharing between simulations, throws when not.
131  *
132  * Should be called simultaneously on the master rank of every simulation.
133  * Note that this only checks for technical compatibility. It is up to
134  * the user to check that the sharing physically makes sense.
135  * Throws an exception when shared biases are not compatible.
136  *
137  * \param[in] awhParams     The AWH parameters.
138  * \param[in] pointSize     Vector of grid-point sizes for each bias.
139  * \param[in] biasSharing   Object for communication for sharing bias data over simulations.
140  */
141 void biasesAreCompatibleForSharingBetweenSimulations(const AwhParams&       awhParams,
142                                                      ArrayRef<const size_t> pointSize,
143                                                      const BiasSharing&     biasSharing);
144
145 } // namespace gmx
146
147 #endif /* GMX_AWH_BIASSHARING_H */