Fix random typos
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_timings.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,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  *  \brief Implements PME GPU timing events wrappers.
38  *
39  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42
43 #include "gmxpre.h"
44
45 #include "pme_gpu_timings.h"
46
47 #include "gromacs/utility/gmxassert.h"
48
49 #include "pme_gpu_internal.h"
50 #include "pme_gpu_types_host.h"
51 #include "pme_gpu_types_host_impl.h"
52
53 bool pme_gpu_timings_enabled(const PmeGpu* pmeGpu)
54 {
55     return pmeGpu->archSpecific->useTiming;
56 }
57
58 void pme_gpu_start_timing(const PmeGpu* pmeGpu, PmeStage pmeStageId)
59 {
60     if (pme_gpu_timings_enabled(pmeGpu))
61     {
62         GMX_ASSERT(pmeStageId < PmeStage::Count, "Wrong PME GPU timing event index");
63         pmeGpu->archSpecific->timingEvents[pmeStageId].openTimingRegion(pmeGpu->archSpecific->pmeStream_);
64     }
65 }
66
67 void pme_gpu_stop_timing(const PmeGpu* pmeGpu, PmeStage pmeStageId)
68 {
69     if (pme_gpu_timings_enabled(pmeGpu))
70     {
71         GMX_ASSERT(pmeStageId < PmeStage::Count, "Wrong PME GPU timing event index");
72         pmeGpu->archSpecific->timingEvents[pmeStageId].closeTimingRegion(pmeGpu->archSpecific->pmeStream_);
73     }
74 }
75
76 void pme_gpu_get_timings(const PmeGpu* pmeGpu, gmx_wallclock_gpu_pme_t* timings)
77 {
78     if (pme_gpu_timings_enabled(pmeGpu))
79     {
80         GMX_RELEASE_ASSERT(timings, "Null GPU timing pointer");
81         for (auto key : keysOf(timings->timing))
82         {
83             timings->timing[key].t = pmeGpu->archSpecific->timingEvents[key].getTotalTime();
84             timings->timing[key].c = pmeGpu->archSpecific->timingEvents[key].getCallCount();
85         }
86     }
87 }
88
89 void pme_gpu_update_timings(const PmeGpu* pmeGpu)
90 {
91     if (pme_gpu_timings_enabled(pmeGpu))
92     {
93         for (const auto& activeTimer : pmeGpu->archSpecific->activeTimers)
94         {
95             pmeGpu->archSpecific->timingEvents[activeTimer].getLastRangeTime();
96         }
97     }
98 }
99
100 void pme_gpu_reinit_timings(const PmeGpu* pmeGpu)
101 {
102     if (pme_gpu_timings_enabled(pmeGpu))
103     {
104         pmeGpu->archSpecific->activeTimers.clear();
105         pmeGpu->archSpecific->activeTimers.insert(PmeStage::SplineAndSpread);
106         const auto& settings = pme_gpu_settings(pmeGpu);
107         // TODO: no separate gtPME_SPLINE and gtPME_SPREAD as they are not used currently
108         if (settings.performGPUFFT)
109         {
110             pmeGpu->archSpecific->activeTimers.insert(PmeStage::FftTransformC2R);
111             pmeGpu->archSpecific->activeTimers.insert(PmeStage::FftTransformR2C);
112         }
113         if (settings.performGPUSolve)
114         {
115             pmeGpu->archSpecific->activeTimers.insert(PmeStage::Solve);
116         }
117         if (settings.performGPUGather)
118         {
119             pmeGpu->archSpecific->activeTimers.insert(PmeStage::Gather);
120         }
121     }
122 }
123
124 void pme_gpu_reset_timings(const PmeGpu* pmeGpu)
125 {
126     if (pme_gpu_timings_enabled(pmeGpu))
127     {
128         for (auto key : keysOf(pmeGpu->archSpecific->timingEvents))
129         {
130             pmeGpu->archSpecific->timingEvents[key].reset();
131         }
132     }
133 }