Add wallcycle counting to StatePropagatorDataGpu
authorSzilárd Páll <pall.szilard@gmail.com>
Fri, 6 Dec 2019 18:57:24 +0000 (19:57 +0100)
committerArtem Zhmurov <zhmurov@gmail.com>
Tue, 17 Dec 2019 15:57:00 +0000 (16:57 +0100)
Launch overheads are counted in the main GPU launch overhead counter and
a separate subcounter is used for the launch and a main counter for the
CPU blocking wait timing.

Note that this chnge introduces mdtypes->timing->mdtypes cyclic dependency,
the warning on which is suppressed.

Fixes #3207

Change-Id: I3b69df9e4888800b43712a42b863958db80f5caa

docs/doxygen/cycle-suppressions.txt
src/gromacs/ewald/pme_only.cpp
src/gromacs/ewald/tests/pmetestcommon.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/mdtypes/state_propagator_data_gpu.h
src/gromacs/mdtypes/state_propagator_data_gpu_impl.cpp
src/gromacs/mdtypes/state_propagator_data_gpu_impl.h
src/gromacs/mdtypes/state_propagator_data_gpu_impl_gpu.cpp
src/gromacs/timing/wallcycle.cpp
src/gromacs/timing/wallcycle.h

index 577d39c85673b43b26e1627586ae70cc94512957..d49152a172b9f04049128182eee6b24a4a71a7e2 100644 (file)
@@ -26,3 +26,7 @@ utility -> math
 
 # modular simulator uses shellfc from mdrun, but is later included in mdrun by simulator builder
 modularsimulator -> mdrun
+
+# Cycle counters in timing use comrec for the set up, which is in the mdtypes. This introduces
+# cyclic dependencies if the cycle counting is used anywhere in mdtypes.
+timing -> mdtypes
\ No newline at end of file
index be7a977daf1313fbf6f982727bd2de86dde75a8a..4e66664ae22663f3c2afa15cd2bfaae033dc1a7b 100644 (file)
@@ -643,7 +643,7 @@ int gmx_pmeonly(struct gmx_pme_t*         pme,
         // TODO: Special PME-only constructor is used here. There is no mechanism to prevent from using the other constructor here.
         //       This should be made safer.
         stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
-                commandStream, deviceContext, GpuApiCallBehavior::Async, paddingSize);
+                commandStream, deviceContext, GpuApiCallBehavior::Async, paddingSize, wcycle);
     }
 
 
index a57195d1f6f6fe0729a025052bc33bc04d0abaae..fc68e96f4131af69cbc38c607be52f92230f1946 100644 (file)
@@ -166,7 +166,7 @@ std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme
     //       restrict one from using other constructor here.
     return std::make_unique<StatePropagatorDataGpu>(
             pme_gpu_get_device_stream(&pme), pme_gpu_get_device_context(&pme),
-            GpuApiCallBehavior::Sync, pme_gpu_get_padding_size(&pme));
+            GpuApiCallBehavior::Sync, pme_gpu_get_padding_size(&pme), nullptr);
 }
 
 //! PME initialization with atom data
index 2a1a0be8b0a3a5934290851b4abb8495fdbd3f93..1b699c4a4d1d1470bc51f2ab25764f46ef0c0505 100644 (file)
@@ -1589,7 +1589,7 @@ int Mdrunner::mdrunner()
                                                       : GpuApiCallBehavior::Sync;
 
             stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
-                    pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize);
+                    pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
             fr->stateGpu = stateGpu.get();
         }
 
index c6ae19c589c0101759ddd7b135f2161f282872e8..697f9d2c861053326f139df573bd4d39c273570e 100644 (file)
@@ -60,6 +60,7 @@
 #include "locality.h"
 
 class GpuEventSynchronizer;
+struct gmx_wallcycle;
 
 namespace gmx
 {
@@ -108,13 +109,15 @@ public:
      *  \param[in] deviceContext   Device context, nullptr allowed.
      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
      *  \param[in] paddingSize     Padding size for coordinates buffer.
+     *  \param[in] wcycle          Wall cycle counter data.
      */
     StatePropagatorDataGpu(const void*        pmeStream,
                            const void*        localStream,
                            const void*        nonLocalStream,
                            const void*        deviceContext,
                            GpuApiCallBehavior transferKind,
-                           int                paddingSize);
+                           int                paddingSize,
+                           gmx_wallcycle*     wcycle);
 
     /*! \brief Constructor to use in PME-only rank and in tests.
      *
@@ -130,11 +133,13 @@ public:
      *  \param[in] deviceContext   Device context, nullptr allowed for non-OpenCL builds.
      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
      *  \param[in] paddingSize     Padding size for coordinates buffer.
+     *  \param[in] wcycle          Wall cycle counter data.
      */
     StatePropagatorDataGpu(const void*        pmeStream,
                            const void*        deviceContext,
                            GpuApiCallBehavior transferKind,
-                           int                paddingSize);
+                           int                paddingSize,
+                           gmx_wallcycle*     wcycle);
 
     //! Move constructor
     StatePropagatorDataGpu(StatePropagatorDataGpu&& other) noexcept;
index 9b201b9376f2d35aa5f85af8b0a7c5a1966c08f8..1e2eadaff463886f9b683890f6993c07d3b5dd71 100644 (file)
@@ -59,7 +59,8 @@ StatePropagatorDataGpu::StatePropagatorDataGpu(const void* /* pmeStream       */
                                                const void* /* nonLocalStream  */,
                                                const void* /* deviceContext   */,
                                                GpuApiCallBehavior /* transferKind    */,
-                                               int /* paddingSize     */) :
+                                               int /* paddingSize     */,
+                                               gmx_wallcycle* /*   wcycle */) :
     impl_(nullptr)
 {
 }
@@ -67,7 +68,8 @@ StatePropagatorDataGpu::StatePropagatorDataGpu(const void* /* pmeStream       */
 StatePropagatorDataGpu::StatePropagatorDataGpu(const void* /* pmeStream       */,
                                                const void* /* deviceContext   */,
                                                GpuApiCallBehavior /* transferKind    */,
-                                               int /* paddingSize     */) :
+                                               int /* paddingSize     */,
+                                               gmx_wallcycle* /*   wcycle */) :
     impl_(nullptr)
 {
 }
index af073b6284dde35d568465192004b375be8b2dd0..9ac9711d298862052e1cc3c138f93ba8ea3a292c 100644 (file)
@@ -58,6 +58,8 @@
 #include "gromacs/utility/classhelpers.h"
 #include "gromacs/utility/enumerationhelpers.h"
 
+struct gmx_wallcycle;
+
 namespace gmx
 {
 
@@ -108,13 +110,15 @@ public:
      *  \param[in] deviceContext   Device context, nullptr allowed.
      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
      *  \param[in] paddingSize     Padding size for coordinates buffer.
+     *  \param[in] wcycle          Wall cycle counter data.
      */
     Impl(const void*        pmeStream,
          const void*        localStream,
          const void*        nonLocalStream,
          const void*        deviceContext,
          GpuApiCallBehavior transferKind,
-         int                paddingSize);
+         int                paddingSize,
+         gmx_wallcycle*     wcycle);
 
     /*! \brief Constructor to use in PME-only rank and in tests.
      *
@@ -130,8 +134,13 @@ public:
      *  \param[in] deviceContext   Device context, nullptr allowed for non-OpenCL builds.
      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
      *  \param[in] paddingSize     Padding size for coordinates buffer.
+     *  \param[in] wcycle          Wall cycle counter data.
      */
-    Impl(const void* pmeStream, const void* deviceContext, GpuApiCallBehavior transferKind, int paddingSize);
+    Impl(const void*        pmeStream,
+         const void*        deviceContext,
+         GpuApiCallBehavior transferKind,
+         int                paddingSize,
+         gmx_wallcycle*     wcycle);
 
     ~Impl();
 
@@ -406,6 +415,9 @@ private:
     //! Allocation size for the force buffer
     int d_fCapacity_ = -1;
 
+    //! \brief Pointer to wallcycle structure.
+    gmx_wallcycle* wcycle_;
+
     /*! \brief Performs the copy of data from host to device buffer.
      *
      * \todo Template on locality.
index f42ad7230e9b220eaca5be3cae2ebb401fa8390f..29821a43133ddaacd8078a7c2ae2596c95aa2037 100644 (file)
 #    endif
 #    include "gromacs/math/vectypes.h"
 #    include "gromacs/mdtypes/state_propagator_data_gpu.h"
+#    include "gromacs/timing/wallcycle.h"
 #    include "gromacs/utility/classhelpers.h"
 
 #    include "state_propagator_data_gpu_impl.h"
 
+
 namespace gmx
 {
 
@@ -67,9 +69,11 @@ StatePropagatorDataGpu::Impl::Impl(const void*        pmeStream,
                                    const void*        nonLocalStream,
                                    const void*        deviceContext,
                                    GpuApiCallBehavior transferKind,
-                                   int                paddingSize) :
+                                   int                paddingSize,
+                                   gmx_wallcycle*     wcycle) :
     transferKind_(transferKind),
-    paddingSize_(paddingSize)
+    paddingSize_(paddingSize),
+    wcycle_(wcycle)
 {
     static_assert(GMX_GPU != GMX_GPU_NONE,
                   "This object should only be constructed on the GPU code-paths.");
@@ -131,9 +135,11 @@ StatePropagatorDataGpu::Impl::Impl(const void*        pmeStream,
 StatePropagatorDataGpu::Impl::Impl(const void*        pmeStream,
                                    const void*        deviceContext,
                                    GpuApiCallBehavior transferKind,
-                                   int                paddingSize) :
+                                   int                paddingSize,
+                                   gmx_wallcycle*     wcycle) :
     transferKind_(transferKind),
-    paddingSize_(paddingSize)
+    paddingSize_(paddingSize),
+    wcycle_(wcycle)
 {
     static_assert(GMX_GPU != GMX_GPU_NONE,
                   "This object should only be constructed on the GPU code-paths.");
@@ -171,6 +177,9 @@ StatePropagatorDataGpu::Impl::~Impl() {}
 
 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
 {
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     numAtomsLocal_ = numAtomsLocal;
     numAtomsAll_   = numAtomsAll;
 
@@ -203,6 +212,9 @@ void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
     {
         clearDeviceBufferAsync(&d_f_, 0, d_fCapacity_, localStream_);
     }
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
@@ -243,11 +255,13 @@ void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<float>
                                                 AtomLocality                         atomLocality,
                                                 CommandStream                        commandStream)
 {
-
     GMX_UNUSED_VALUE(dataSize);
 
     GMX_ASSERT(dataSize >= 0, "Trying to copy to device buffer before it was allocated.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     int atomsStartAt, numAtomsToCopy;
     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
 
@@ -264,6 +278,9 @@ void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<float>
         copyToDeviceBuffer(&d_data, reinterpret_cast<const float*>(&h_data.data()[atomsStartAt]),
                            elementsStartAt, numElementsToCopy, commandStream, transferKind_, nullptr);
     }
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
@@ -272,11 +289,13 @@ void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_dat
                                                   AtomLocality             atomLocality,
                                                   CommandStream            commandStream)
 {
-
     GMX_UNUSED_VALUE(dataSize);
 
     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     int atomsStartAt, numAtomsToCopy;
     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
 
@@ -293,6 +312,9 @@ void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_dat
         copyFromDeviceBuffer(reinterpret_cast<float*>(&h_data.data()[atomsStartAt]), &d_data,
                              elementsStartAt, numElementsToCopy, commandStream, transferKind_, nullptr);
     }
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getCoordinates()
@@ -308,6 +330,9 @@ void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<cons
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying positions with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, commandStream);
 
     // markEvent is skipped in OpenCL as:
@@ -318,6 +343,9 @@ void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<cons
     {
         xReadyOnDevice_[atomLocality].markEvent(commandStream);
     }
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 GpuEventSynchronizer*
@@ -349,8 +377,10 @@ StatePropagatorDataGpu::Impl::getCoordinatesReadyOnDeviceEvent(AtomLocality atom
 
 void StatePropagatorDataGpu::Impl::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
 {
+    wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
     xReadyOnDevice_[atomLocality].waitForEvent();
+    wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
 }
 
 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::xUpdatedOnDevice()
@@ -365,14 +395,22 @@ void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVe
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying positions with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, commandStream);
     // Note: unlike copyCoordinatesToGpu this is not used in OpenCL, and the conditional is not needed.
     xReadyOnHost_[atomLocality].markEvent(commandStream);
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 void StatePropagatorDataGpu::Impl::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
 {
+    wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
     xReadyOnHost_[atomLocality].waitForEvent();
+    wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
 }
 
 
@@ -389,8 +427,14 @@ void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying velocities with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, commandStream);
     vReadyOnDevice_[atomLocality].markEvent(commandStream);
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
@@ -406,13 +450,21 @@ void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying velocities with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, commandStream);
     vReadyOnHost_[atomLocality].markEvent(commandStream);
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 void StatePropagatorDataGpu::Impl::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
 {
+    wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
     vReadyOnHost_[atomLocality].waitForEvent();
+    wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
 }
 
 
@@ -429,8 +481,14 @@ void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying forces with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, commandStream);
     fReadyOnDevice_[atomLocality].markEvent(commandStream);
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
@@ -458,13 +516,21 @@ void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_
     GMX_ASSERT(commandStream != nullptr,
                "No stream is valid for copying forces with given atom locality.");
 
+    wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
+    wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+
     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, commandStream);
     fReadyOnHost_[atomLocality].markEvent(commandStream);
+
+    wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
+    wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
 }
 
 void StatePropagatorDataGpu::Impl::waitForcesReadyOnHost(AtomLocality atomLocality)
 {
+    wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
     fReadyOnHost_[atomLocality].waitForEvent();
+    wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
 }
 
 void* StatePropagatorDataGpu::Impl::getUpdateStream()
@@ -488,16 +554,18 @@ StatePropagatorDataGpu::StatePropagatorDataGpu(const void*        pmeStream,
                                                const void*        nonLocalStream,
                                                const void*        deviceContext,
                                                GpuApiCallBehavior transferKind,
-                                               int                paddingSize) :
-    impl_(new Impl(pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize))
+                                               int                paddingSize,
+                                               gmx_wallcycle*     wcycle) :
+    impl_(new Impl(pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle))
 {
 }
 
 StatePropagatorDataGpu::StatePropagatorDataGpu(const void*        pmeStream,
                                                const void*        deviceContext,
                                                GpuApiCallBehavior transferKind,
-                                               int                paddingSize) :
-    impl_(new Impl(pmeStream, deviceContext, transferKind, paddingSize))
+                                               int                paddingSize,
+                                               gmx_wallcycle*     wcycle) :
+    impl_(new Impl(pmeStream, deviceContext, transferKind, paddingSize, wcycle))
 {
 }
 
index 59a55b8bd1ae5ca52e4b40a47f1317d0c729ddfb..638b8dc475e7da59a742f77479e202089109fdb0 100644 (file)
@@ -132,6 +132,7 @@ static const char* wcn[ewcNR] = { "Run",
                                   "Reduce GPU PME F",
                                   "Wait GPU NB nonloc.",
                                   "Wait GPU NB local",
+                                  "Wait GPU state copy",
                                   "NB X/F buffer ops.",
                                   "Vsite spread",
                                   "COM pull force",
@@ -168,6 +169,7 @@ static const char* wcsn[ewcsNR] = {
     "Launch NB GPU tasks",
     "Launch Bonded GPU tasks",
     "Launch PME GPU tasks",
+    "Launch state copy",
     "Ewald F correction",
     "NB X buffer ops.",
     "NB F buffer ops.",
index fa65919c609a794fa2aded880ae370ec06ef4691..c2d5bc8bf51798a5dc8514c85377b5d12ac8cbd9 100644 (file)
@@ -81,6 +81,7 @@ enum
     ewcPME_GPU_F_REDUCTION,
     ewcWAIT_GPU_NB_NL,
     ewcWAIT_GPU_NB_L,
+    ewcWAIT_GPU_STATE_PROPAGATOR_DATA,
     ewcNB_XF_BUF_OPS,
     ewcVSITESPREAD,
     ewcPULLPOT,
@@ -120,6 +121,7 @@ enum
     ewcsLAUNCH_GPU_NONBONDED,
     ewcsLAUNCH_GPU_BONDED,
     ewcsLAUNCH_GPU_PME,
+    ewcsLAUNCH_STATE_PROPAGATOR_DATA,
     ewcsEWALD_CORRECTION,
     ewcsNB_X_BUF_OPS,
     ewcsNB_F_BUF_OPS,