Move the PmeGpuSpecific structure into pme-gpu-types-host-impl.h
[gromacs.git] / src / gromacs / ewald / pme-timings.cu
blobdd6bc3ea386b3a3e3806c03561a16672035c168e
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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  */
36 /*! \internal \file
37  *  \brief Implements PME GPU timing events in CUDA.
38  *
39  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
40  */
42 #include "gmxpre.h"
44 #include "pme-timings.cuh"
46 #include "gromacs/utility/gmxassert.h"
48 #include "pme-gpu-internal.h"
49 #include "pme-gpu-types-host.h"
50 #include "pme-gpu-types-host-impl.h"
52 /*! \brief \internal
53  * Tells if CUDA-based performance tracking is enabled for PME.
54  *
55  * \param[in] pme            The PME data structure.
56  * \returns                  True if timings are enabled, false otherwise.
57  */
58 inline bool pme_gpu_timings_enabled(const PmeGpu *pmeGpu)
60     return pmeGpu->archSpecific->useTiming;
63 void pme_gpu_start_timing(const PmeGpu *pmeGpu, size_t PMEStageId)
65     if (pme_gpu_timings_enabled(pmeGpu))
66     {
67         GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
68         pmeGpu->archSpecific->timingEvents[PMEStageId].openTimingRegion(pmeGpu->archSpecific->pmeStream);
69     }
72 CommandEvent *pme_gpu_fetch_timing_event(const PmeGpu *pmeGpu, size_t PMEStageId)
74     CommandEvent *timingEvent = nullptr;
75     if (pme_gpu_timings_enabled(pmeGpu))
76     {
77         GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
78         timingEvent = pmeGpu->archSpecific->timingEvents[PMEStageId].fetchNextEvent();
79     }
80     return timingEvent;
83 void pme_gpu_stop_timing(const PmeGpu *pmeGpu, size_t PMEStageId)
85     if (pme_gpu_timings_enabled(pmeGpu))
86     {
87         GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
88         pmeGpu->archSpecific->timingEvents[PMEStageId].closeTimingRegion(pmeGpu->archSpecific->pmeStream);
89     }
92 void pme_gpu_get_timings(const PmeGpu *pmeGpu, gmx_wallclock_gpu_pme_t *timings)
94     if (pme_gpu_timings_enabled(pmeGpu))
95     {
96         GMX_RELEASE_ASSERT(timings, "Null GPU timing pointer");
97         for (size_t i = 0; i < pmeGpu->archSpecific->timingEvents.size(); i++)
98         {
99             timings->timing[i].t = pmeGpu->archSpecific->timingEvents[i].getTotalTime();
100             timings->timing[i].c = pmeGpu->archSpecific->timingEvents[i].getCallCount();
101         }
102     }
105 void pme_gpu_update_timings(const PmeGpu *pmeGpu)
107     if (pme_gpu_timings_enabled(pmeGpu))
108     {
109         for (const size_t &activeTimer : pmeGpu->archSpecific->activeTimers)
110         {
111             pmeGpu->archSpecific->timingEvents[activeTimer].getLastRangeTime();
112         }
113     }
116 void pme_gpu_reinit_timings(const PmeGpu *pmeGpu)
118     if (pme_gpu_timings_enabled(pmeGpu))
119     {
120         pmeGpu->archSpecific->activeTimers.clear();
121         pmeGpu->archSpecific->activeTimers.insert(gtPME_SPLINEANDSPREAD);
122         // TODO: no separate gtPME_SPLINE and gtPME_SPREAD as they are not used currently
123         if (pme_gpu_performs_FFT(pmeGpu))
124         {
125             pmeGpu->archSpecific->activeTimers.insert(gtPME_FFT_C2R);
126             pmeGpu->archSpecific->activeTimers.insert(gtPME_FFT_R2C);
127         }
128         if (pme_gpu_performs_solve(pmeGpu))
129         {
130             pmeGpu->archSpecific->activeTimers.insert(gtPME_SOLVE);
131         }
132         if (pme_gpu_performs_gather(pmeGpu))
133         {
134             pmeGpu->archSpecific->activeTimers.insert(gtPME_GATHER);
135         }
136     }
139 void pme_gpu_reset_timings(const PmeGpu *pmeGpu)
141     if (pme_gpu_timings_enabled(pmeGpu))
142     {
143         for (size_t i = 0; i < pmeGpu->archSpecific->timingEvents.size(); i++)
144         {
145             pmeGpu->archSpecific->timingEvents[i].reset();
146         }
147     }