Add coolquotes
[gromacs.git] / src / gromacs / ewald / pme-gpu-timings.cpp
blob5c218c14104d13657f49d1e51e22c385365ffc85
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
36 /*! \internal \file
37 * \brief Implements PME GPU timing events wrappers.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
43 #include "gmxpre.h"
45 #include "pme-gpu-timings.h"
47 #include "gromacs/utility/gmxassert.h"
49 #include "pme-gpu-internal.h"
50 #include "pme-gpu-types-host.h"
51 #include "pme-gpu-types-host-impl.h"
53 /*! \brief
54 * Tells if CUDA-based performance tracking is enabled for PME.
56 * \param[in] pmeGpu The PME GPU data structure.
57 * \returns True if timings are enabled, false otherwise.
59 inline bool pme_gpu_timings_enabled(const PmeGpu *pmeGpu)
61 return pmeGpu->archSpecific->useTiming;
64 void pme_gpu_start_timing(const PmeGpu *pmeGpu, size_t PMEStageId)
66 if (pme_gpu_timings_enabled(pmeGpu))
68 GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
69 pmeGpu->archSpecific->timingEvents[PMEStageId].openTimingRegion(pmeGpu->archSpecific->pmeStream);
73 CommandEvent *pme_gpu_fetch_timing_event(const PmeGpu *pmeGpu, size_t PMEStageId)
75 CommandEvent *timingEvent = nullptr;
76 if (pme_gpu_timings_enabled(pmeGpu))
78 GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
79 timingEvent = pmeGpu->archSpecific->timingEvents[PMEStageId].fetchNextEvent();
81 return timingEvent;
84 void pme_gpu_stop_timing(const PmeGpu *pmeGpu, size_t PMEStageId)
86 if (pme_gpu_timings_enabled(pmeGpu))
88 GMX_ASSERT(PMEStageId < pmeGpu->archSpecific->timingEvents.size(), "Wrong PME GPU timing event index");
89 pmeGpu->archSpecific->timingEvents[PMEStageId].closeTimingRegion(pmeGpu->archSpecific->pmeStream);
93 void pme_gpu_get_timings(const PmeGpu *pmeGpu, gmx_wallclock_gpu_pme_t *timings)
95 if (pme_gpu_timings_enabled(pmeGpu))
97 GMX_RELEASE_ASSERT(timings, "Null GPU timing pointer");
98 for (size_t i = 0; i < pmeGpu->archSpecific->timingEvents.size(); i++)
100 timings->timing[i].t = pmeGpu->archSpecific->timingEvents[i].getTotalTime();
101 timings->timing[i].c = pmeGpu->archSpecific->timingEvents[i].getCallCount();
106 void pme_gpu_update_timings(const PmeGpu *pmeGpu)
108 if (pme_gpu_timings_enabled(pmeGpu))
110 for (const size_t &activeTimer : pmeGpu->archSpecific->activeTimers)
112 pmeGpu->archSpecific->timingEvents[activeTimer].getLastRangeTime();
117 void pme_gpu_reinit_timings(const PmeGpu *pmeGpu)
119 if (pme_gpu_timings_enabled(pmeGpu))
121 pmeGpu->archSpecific->activeTimers.clear();
122 pmeGpu->archSpecific->activeTimers.insert(gtPME_SPLINEANDSPREAD);
123 // TODO: no separate gtPME_SPLINE and gtPME_SPREAD as they are not used currently
124 if (pme_gpu_performs_FFT(pmeGpu))
126 pmeGpu->archSpecific->activeTimers.insert(gtPME_FFT_C2R);
127 pmeGpu->archSpecific->activeTimers.insert(gtPME_FFT_R2C);
129 if (pme_gpu_performs_solve(pmeGpu))
131 pmeGpu->archSpecific->activeTimers.insert(gtPME_SOLVE);
133 if (pme_gpu_performs_gather(pmeGpu))
135 pmeGpu->archSpecific->activeTimers.insert(gtPME_GATHER);
140 void pme_gpu_reset_timings(const PmeGpu *pmeGpu)
142 if (pme_gpu_timings_enabled(pmeGpu))
144 for (size_t i = 0; i < pmeGpu->archSpecific->timingEvents.size(); i++)
146 pmeGpu->archSpecific->timingEvents[i].reset();