SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / domdec / dlbtiming.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,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.
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 /*! \libinternal \file
37  *
38  * \brief This file declares functions for timing the load imbalance due to domain decomposition.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \inlibraryapi
42  * \ingroup module_domdec
43  */
44
45 #ifndef GMX_DOMDEC_DLBTIMING_H
46 #define GMX_DOMDEC_DLBTIMING_H
47
48 #include "gromacs/mdtypes/commrec.h"
49
50 struct BalanceRegion;
51 struct gmx_domdec_t;
52 struct t_nrnb;
53
54 /*! \brief Tells if we should open the balancing region */
55 enum class DdAllowBalanceRegionReopen
56 {
57     no, //!< Do not allow opening an already open region
58     yes //!< Allow opening an already open region
59 };
60
61 /*! \brief Tells if we had to wait for a GPU to finish computation */
62 enum class DdBalanceRegionWaitedForGpu
63 {
64     no, //!< The GPU finished computation before the CPU needed the result
65     yes //!< We had to wait for the GPU to finish computation
66 };
67
68 /*! \brief Re-open the, already opened, load balance timing region
69  *
70  * This function should be called after every MPI communication that occurs
71  * in the main MD loop.
72  * Note that the current setup assumes that all MPI communication acts like
73  * a global barrier. But if some ranks don't participate in communication
74  * or if some ranks communicate faster with neighbors than others,
75  * the obtained timings might not accurately reflect the computation time.
76  *
77  * \param[in,out] dd  The domain decomposition struct
78  */
79 void ddReopenBalanceRegionCpu(const gmx_domdec_t* dd);
80
81 /*! \libinternal
82  * \brief Manager for starting and stopping the dynamic load balancing region
83  */
84 class DDBalanceRegionHandler
85 {
86 public:
87     //! Constructor, pass a pointer to t_commrec or nullptr when not using domain decomposition
88     DDBalanceRegionHandler(const t_commrec* cr) :
89         useBalancingRegion_(cr != nullptr ? havePPDomainDecomposition(cr) : false),
90         dd_(cr != nullptr ? cr->dd : nullptr)
91     {
92     }
93
94     /*! \brief Returns whether were are actually using the balancing region
95      */
96     bool useBalancingRegion() const { return useBalancingRegion_; }
97
98     /*! \brief Open the load balance timing region on the CPU
99      *
100      * Opens the balancing region for timing how much time it takes to perform
101      * the (balancable part of) the MD step. This should be called right after
102      * the last communication during the previous step to maximize the region.
103      * In practice this means right after the force communication finished
104      * or just before neighbor search at search steps.
105      * It is assumed that computation done in the region either scales along
106      * with the domain size or takes constant time.
107      *
108      * \param[in] allowReopen  Allows calling with a potentially already opened region
109      */
110     void openBeforeForceComputationCpu(DdAllowBalanceRegionReopen allowReopen) const
111     {
112         if (useBalancingRegion_)
113         {
114             openRegionCpuImpl(allowReopen);
115         }
116     }
117
118     /*! \brief Open the load balance timing region for the CPU
119      *
120      * This can only be called within a region that is open on the CPU side.
121      */
122     void openBeforeForceComputationGpu() const
123     {
124         if (useBalancingRegion_)
125         {
126             openRegionGpuImpl();
127         }
128     }
129
130     /*! \brief Re-open the, already opened, load balance timing region
131      *
132      * This function should be called after every MPI communication that occurs
133      * in the main MD loop.
134      * Note that the current setup assumes that all MPI communication acts like
135      * a global barrier. But if some ranks don't participate in communication
136      * or if some ranks communicate faster with neighbors than others,
137      * the obtained timings might not accurately reflect the computation time.
138      */
139     void reopenRegionCpu() const
140     {
141         if (useBalancingRegion_)
142         {
143             ddReopenBalanceRegionCpu(dd_);
144         }
145     }
146
147     /*! \brief Close the load balance timing region on the CPU side
148      */
149     void closeAfterForceComputationCpu() const
150     {
151         if (useBalancingRegion_)
152         {
153             closeRegionCpuImpl();
154         }
155     }
156
157     /*! \brief Close the load balance timing region on the GPU side
158      *
159      * This should be called after the CPU receives the last (local) results
160      * from the GPU. The wait time for these results is estimated, depending
161      * on the \p waitedForGpu parameter.
162      * If called on an already closed region, this call does nothing.
163      *
164      * \param[in] waitCyclesGpuInCpuRegion  The time we waited for the GPU earlier, overlapping completely with the open CPU region
165      * \param[in] waitedForGpu              Tells if we waited for the GPU to finish now
166      */
167     void closeAfterForceComputationGpu(float                       waitCyclesGpuInCpuRegion,
168                                        DdBalanceRegionWaitedForGpu waitedForGpu) const
169     {
170         if (useBalancingRegion_)
171         {
172             closeRegionGpuImpl(waitCyclesGpuInCpuRegion, waitedForGpu);
173         }
174     }
175
176 private:
177     /*! \brief Open the load balance timing region on the CPU
178      *
179      * \param[in] allowReopen  Allows calling with a potentially already opened region
180      */
181     void openRegionCpuImpl(DdAllowBalanceRegionReopen allowReopen) const;
182
183     /*! \brief Open the load balance timing region for the GPU
184      *
185      * This can only be called within a region that is open on the CPU side.
186      */
187     void openRegionGpuImpl() const;
188
189     /*! \brief Close the load balance timing region on the CPU side
190      */
191     void closeRegionCpuImpl() const;
192
193     /*! \brief Close the load balance timing region on the GPU side
194      *
195      * \param[in] waitCyclesGpuInCpuRegion  The time we waited for the GPU earlier, overlapping completely with the open CPU region
196      * \param[in] waitedForGpu              Tells if we waited for the GPU to finish now
197      */
198     void closeRegionGpuImpl(float waitCyclesGpuInCpuRegion, DdBalanceRegionWaitedForGpu waitedForGpu) const;
199
200     //! Tells whether the balancing region should be active
201     bool useBalancingRegion_;
202     //! A pointer to the DD struct, only valid with useBalancingRegion_=true
203     gmx_domdec_t* dd_;
204 };
205
206 /*! \brief Returns a pointer to a constructed \p BalanceRegion struct
207  *
208  * Should be replaced by a proper constructor once BalanceRegion is a proper
209  * class (requires restructering in domdec.cpp).
210  */
211 BalanceRegion* ddBalanceRegionAllocate();
212
213 /*! \brief Start the force flop count */
214 void dd_force_flop_start(struct gmx_domdec_t* dd, t_nrnb* nrnb);
215
216 /*! \brief Stop the force flop count */
217 void dd_force_flop_stop(struct gmx_domdec_t* dd, t_nrnb* nrnb);
218
219 //! Clear the cycle counts used for tuning.
220 void clear_dd_cycle_counts(gmx_domdec_t* dd);
221
222 #endif