222f08c20b7b770376d21db0e679e3e594a5e549
[alexxy/gromacs.git] / src / gromacs / gpu_utils / devicebuffer_sycl.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 #ifndef GMX_GPU_UTILS_DEVICEBUFFER_SYCL_H
36 #define GMX_GPU_UTILS_DEVICEBUFFER_SYCL_H
37
38 /*! \libinternal \file
39  *  \brief Implements the DeviceBuffer type and routines for SYCL.
40  *  Should only be included directly by the main DeviceBuffer file devicebuffer.h.
41  *  TODO: the intent is for DeviceBuffer to become a class.
42  *
43  *  \author Artem Zhmurov <zhmurov@gmail.com>
44  *  \author Erik Lindahl <erik.lindahl@gmail.com>
45  *  \author Andrey Alekseenko <al42and@gmail.com>
46  *
47  *  \inlibraryapi
48  */
49
50 #include <utility>
51
52 #include "gromacs/gpu_utils/device_context.h"
53 #include "gromacs/gpu_utils/device_stream.h"
54 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
55 #include "gromacs/gpu_utils/gmxsycl.h"
56 #include "gromacs/gpu_utils/gpu_utils.h" //only for GpuApiCallBehavior
57 #include "gromacs/gpu_utils/gputraits_sycl.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/stringutil.h"
60
61 #ifndef DOXYGEN
62 template<typename T>
63 class DeviceBuffer<T>::ClSyclBufferWrapper : public cl::sycl::buffer<T, 1>
64 {
65     using cl::sycl::buffer<T, 1>::buffer; // Get all the constructors
66 };
67
68 template<typename T>
69 using ClSyclBufferWrapper = typename DeviceBuffer<T>::ClSyclBufferWrapper;
70
71 //! Constructor.
72 template<typename T>
73 DeviceBuffer<T>::DeviceBuffer() : buffer_(nullptr)
74 {
75 }
76
77 //! Destructor.
78 template<typename T>
79 DeviceBuffer<T>::~DeviceBuffer() = default;
80
81 //! Copy constructor (references the same underlying SYCL buffer).
82 template<typename T>
83 DeviceBuffer<T>::DeviceBuffer(DeviceBuffer<T> const& src) :
84     buffer_(new ClSyclBufferWrapper(*src.buffer_))
85 {
86 }
87
88 //! Move constructor.
89 template<typename T>
90 DeviceBuffer<T>::DeviceBuffer(DeviceBuffer<T>&& src) noexcept = default;
91
92 //! Copy assignment (references the same underlying SYCL buffer).
93 template<typename T>
94 DeviceBuffer<T>& DeviceBuffer<T>::operator=(DeviceBuffer<T> const& src)
95 {
96     buffer_.reset(new ClSyclBufferWrapper(*src.buffer_));
97     return *this;
98 }
99
100 //! Move assignment.
101 template<typename T>
102 DeviceBuffer<T>& DeviceBuffer<T>::operator=(DeviceBuffer<T>&& src) noexcept = default;
103
104 /*! \brief Dummy assignment operator to allow compilation of some cross-platform code.
105  *
106  * A hacky way to make SYCL implementation of DeviceBuffer compatible with details of CUDA and
107  * OpenCL implementations.
108  *
109  * \todo Should be removed after DeviceBuffer refactoring.
110  *
111  * \tparam T Type of buffer content.
112  * \param nullPtr \c std::nullptr. Not possible to assign any other pointers.
113  */
114 template<typename T>
115 DeviceBuffer<T>& DeviceBuffer<T>::operator=(std::nullptr_t nullPtr)
116 {
117     buffer_.reset(nullPtr);
118     return *this;
119 }
120
121
122 namespace gmx::internal
123 {
124 //! Shorthand alias to create a placeholder SYCL accessor with chosen data type and access mode.
125 template<class T, enum cl::sycl::access::mode mode>
126 using PlaceholderAccessor =
127         cl::sycl::accessor<T, 1, mode, cl::sycl::access::target::global_buffer, cl::sycl::access::placeholder::true_t>;
128 } // namespace gmx::internal
129
130 /** \brief
131  * Thin wrapper around placeholder accessor that allows implicit construction from \c DeviceBuffer.
132  *
133  * "Placeholder accessor" is an indicator of the intent to create an accessor for certain buffer
134  * of a certain type, that is not yet bound to a specific command group handler (device). Such
135  * accessors can be created outside SYCL kernels, which is helpful if we want to pass them as
136  * function arguments.
137  *
138  * \tparam T Type of buffer content.
139  * \tparam mode Access mode.
140  */
141 template<class T, enum cl::sycl::access::mode mode>
142 class DeviceAccessor : public gmx::internal::PlaceholderAccessor<T, mode>
143 {
144 public:
145     // Inherit all the constructors
146     using gmx::internal::PlaceholderAccessor<T, mode>::PlaceholderAccessor;
147     //! Construct Accessor from DeviceBuffer (must be initialized)
148     DeviceAccessor(DeviceBuffer<T>& buffer) :
149         gmx::internal::PlaceholderAccessor<T, mode>(getSyclBuffer(buffer))
150     {
151     }
152
153 private:
154     //! Helper function to get sycl:buffer object from DeviceBuffer wrapper, with a sanity check.
155     static inline cl::sycl::buffer<T, 1>& getSyclBuffer(DeviceBuffer<T>& buffer)
156     {
157         GMX_ASSERT(bool(buffer), "Trying to construct accessor from an uninitialized buffer");
158         return *buffer.buffer_;
159     }
160 };
161
162 namespace gmx::internal
163 {
164 //! A "blackhole" class to be used when we want to ignore an argument to a function.
165 struct EmptyClassThatIgnoresConstructorArguments
166 {
167     template<class... Args>
168     [[maybe_unused]] EmptyClassThatIgnoresConstructorArguments(Args&&... /*args*/)
169     {
170     }
171 };
172 } // namespace gmx::internal
173
174 /** \brief
175  * Helper class to be used as function argument. Will either correspond to a device accessor, or an empty class.
176  *
177  * Example usage:
178  * \code
179     template <bool doFoo>
180     void getBarKernel(handler& cgh, OptionalAccessor<float, mode::read, doFoo> a_fooPrms)
181     {
182         if constexpr (doFoo)
183             cgh.require(a_fooPrms);
184         // Can only use a_fooPrms if doFoo == true
185     }
186
187     template <bool doFoo>
188     void callBar(DeviceBuffer<float> b_fooPrms)
189     {
190         // If doFoo is false, b_fooPrms will be ignored (can be not initialized).
191         // Otherwise, an accessor will be built (b_fooPrms must be a valid buffer).
192         auto kernel = getBarKernel<doFoo>(b_fooPrms);
193         // If the accessor in not enabled, anything can be passed as its ctor argument.
194         auto kernel2 = getBarKernel<false>(nullptr_t);
195     }
196  * \endcode
197  *
198  * \tparam T Data type of the underlying buffer
199  * \tparam mode Access mode of the accessor
200  * \tparam enabled Compile-time flag indicating whether we want to actually create an accessor.
201  */
202 template<class T, enum cl::sycl::access::mode mode, bool enabled>
203 using OptionalAccessor =
204         std::conditional_t<enabled, DeviceAccessor<T, mode>, gmx::internal::EmptyClassThatIgnoresConstructorArguments>;
205
206 #endif // #ifndef DOXYGEN
207
208 /*! \brief Check the validity of the device buffer.
209  *
210  * Checks if the buffer is valid and if its allocation is big enough.
211  *
212  * \param[in] buffer        Device buffer to be checked.
213  * \param[in] requiredSize  Number of elements that the buffer will have to accommodate.
214  *
215  * \returns Whether the device buffer exists and has enough capacity.
216  */
217 template<typename T>
218 static gmx_unused bool checkDeviceBuffer(const DeviceBuffer<T>& buffer, int requiredSize)
219 {
220     return buffer.buffer_ && (static_cast<int>(buffer.buffer_->get_count()) >= requiredSize);
221 }
222
223 /*! \libinternal \brief
224  * Allocates a device-side buffer.
225  * It is currently a caller's responsibility to call it only on not-yet allocated buffers.
226  *
227  * \tparam        ValueType            Raw value type of the \p buffer.
228  * \param[in,out] buffer               Pointer to the device-side buffer.
229  * \param[in]     numValues            Number of values to accommodate.
230  * \param[in]     deviceContext        The buffer's device context-to-be.
231  */
232 template<typename ValueType>
233 void allocateDeviceBuffer(DeviceBuffer<ValueType>* buffer, size_t numValues, const DeviceContext& deviceContext)
234 {
235     /* SYCL does not require binding buffer to a specific context or device. The ::context_bound
236      * property only enforces the use of only given context, and possibly offers some optimizations */
237     const cl::sycl::property_list bufferProperties{ cl::sycl::property::buffer::context_bound(
238             deviceContext.context()) };
239     buffer->buffer_.reset(
240             new ClSyclBufferWrapper<ValueType>(cl::sycl::range<1>(numValues), bufferProperties));
241 }
242
243 /*! \brief
244  * Frees a device-side buffer.
245  * This does not reset separately stored size/capacity integers,
246  * as this is planned to be a destructor of DeviceBuffer as a proper class,
247  * and no calls on \p buffer should be made afterwards.
248  *
249  * \param[in] buffer  Pointer to the buffer to free.
250  */
251 template<typename ValueType>
252 void freeDeviceBuffer(DeviceBuffer<ValueType>* buffer)
253 {
254     buffer->buffer_.reset(nullptr);
255 }
256
257 /*! \brief
258  * Performs the host-to-device data copy, synchronous or asynchronously on request.
259  *
260  * Unlike in CUDA and OpenCL, synchronous call does not guarantee that all previously
261  * submitted operations are complete, only the ones that are required for \p buffer consistency.
262  *
263  * \tparam        ValueType            Raw value type of the \p buffer.
264  * \param[in,out] buffer               Pointer to the device-side buffer.
265  * \param[in]     hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType.
266  * \param[in]     startingOffset       Offset (in values) at the device-side buffer to copy into.
267  * \param[in]     numValues            Number of values to copy.
268  * \param[in]     deviceStream         GPU stream to perform asynchronous copy in.
269  * \param[in]     transferKind         Copy type: synchronous or asynchronous.
270  * \param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
271  *                                     Ignored in SYCL.
272  */
273 template<typename ValueType>
274 void copyToDeviceBuffer(DeviceBuffer<ValueType>* buffer,
275                         const ValueType*         hostBuffer,
276                         size_t                   startingOffset,
277                         size_t                   numValues,
278                         const DeviceStream&      deviceStream,
279                         GpuApiCallBehavior       transferKind,
280                         CommandEvent* gmx_unused timingEvent)
281 {
282     if (numValues == 0)
283     {
284         return; // such calls are actually made with empty domains
285     }
286     GMX_ASSERT(buffer, "needs a buffer pointer");
287     GMX_ASSERT(hostBuffer, "needs a host buffer pointer");
288
289     GMX_ASSERT(checkDeviceBuffer(*buffer, startingOffset + numValues),
290                "buffer too small or not initialized");
291
292     cl::sycl::buffer<ValueType>& syclBuffer = *buffer->buffer_;
293
294     cl::sycl::event ev = deviceStream.stream().submit([&](cl::sycl::handler& cgh) {
295         /* Here and elsewhere in this file, accessor constructor is user instead of a more common
296          * buffer::get_access, since the compiler (icpx 2021.1-beta09) occasionally gets confused
297          * by all the overloads */
298         auto d_bufferAccessor = cl::sycl::accessor<ValueType, 1, cl::sycl::access::mode::discard_write>{
299             syclBuffer, cgh, cl::sycl::range(numValues), cl::sycl::id(startingOffset)
300         };
301         cgh.copy(hostBuffer, d_bufferAccessor);
302     });
303     if (transferKind == GpuApiCallBehavior::Sync)
304     {
305         ev.wait_and_throw();
306     }
307 }
308
309 /*! \brief
310  * Performs the device-to-host data copy, synchronous or asynchronously on request.
311  *
312  * Unlike in CUDA and OpenCL, synchronous call does not guarantee that all previously
313  * submitted operations are complete, only the ones that are required for \p buffer consistency.
314  *
315  * \tparam        ValueType            Raw value type of the \p buffer.
316  * \param[in,out] hostBuffer           Pointer to the raw host-side memory, also typed \p ValueType
317  * \param[in]     buffer               Pointer to the device-side buffer.
318  * \param[in]     startingOffset       Offset (in values) at the device-side buffer to copy from.
319  * \param[in]     numValues            Number of values to copy.
320  * \param[in]     deviceStream         GPU stream to perform asynchronous copy in.
321  * \param[in]     transferKind         Copy type: synchronous or asynchronous.
322  * \param[out]    timingEvent          A pointer to the H2D copy timing event to be filled in.
323  *                                     Ignored in SYCL.
324  */
325 template<typename ValueType>
326 void copyFromDeviceBuffer(ValueType*               hostBuffer,
327                           DeviceBuffer<ValueType>* buffer,
328                           size_t                   startingOffset,
329                           size_t                   numValues,
330                           const DeviceStream&      deviceStream,
331                           GpuApiCallBehavior       transferKind,
332                           CommandEvent* gmx_unused timingEvent)
333 {
334     if (numValues == 0)
335     {
336         return; // such calls are actually made with empty domains
337     }
338     GMX_ASSERT(buffer, "needs a buffer pointer");
339     GMX_ASSERT(hostBuffer, "needs a host buffer pointer");
340
341     GMX_ASSERT(checkDeviceBuffer(*buffer, startingOffset + numValues),
342                "buffer too small or not initialized");
343
344     cl::sycl::buffer<ValueType>& syclBuffer = *buffer->buffer_;
345
346     cl::sycl::event ev = deviceStream.stream().submit([&](cl::sycl::handler& cgh) {
347         const auto d_bufferAccessor = cl::sycl::accessor<ValueType, 1, cl::sycl::access::mode::read>{
348             syclBuffer, cgh, cl::sycl::range(numValues), cl::sycl::id(startingOffset)
349         };
350         cgh.copy(d_bufferAccessor, hostBuffer);
351     });
352     if (transferKind == GpuApiCallBehavior::Sync)
353     {
354         ev.wait_and_throw();
355     }
356 }
357
358
359 namespace gmx::internal
360 {
361 /*! \brief Helper function to clear device buffer.
362  *
363  * Not applicable to GROMACS's float3 (a.k.a. gmx::RVec) and other custom types.
364  * From SYCL specs: "T must be a scalar value or a SYCL vector type."
365  */
366 template<typename ValueType>
367 cl::sycl::event fillSyclBufferWithNull(cl::sycl::buffer<ValueType, 1>& buffer,
368                                        size_t                          startingOffset,
369                                        size_t                          numValues,
370                                        cl::sycl::queue                 queue)
371 {
372     using cl::sycl::access::mode;
373     const cl::sycl::range<1> range(numValues);
374     const cl::sycl::id<1>    offset(startingOffset);
375     const ValueType pattern = ValueType(0); // SYCL vectors support initialization by scalar
376
377     return queue.submit([&](cl::sycl::handler& cgh) {
378         auto d_bufferAccessor =
379                 cl::sycl::accessor<ValueType, 1, mode::discard_write>{ buffer, cgh, range, offset };
380         cgh.fill(d_bufferAccessor, pattern);
381     });
382 }
383
384 //! \brief Helper function to clear device buffer of type float3.
385 template<>
386 inline cl::sycl::event fillSyclBufferWithNull(cl::sycl::buffer<Float3, 1>& buffer,
387                                               size_t                       startingOffset,
388                                               size_t                       numValues,
389                                               cl::sycl::queue              queue)
390 {
391     cl::sycl::buffer<float, 1> bufferAsFloat = buffer.reinterpret<float, 1>(buffer.get_count() * DIM);
392     return fillSyclBufferWithNull<float>(
393             bufferAsFloat, startingOffset * DIM, numValues * DIM, std::move(queue));
394 }
395 } // namespace gmx::internal
396
397 /*! \brief
398  * Clears the device buffer asynchronously.
399  *
400  * \tparam        ValueType       Raw value type of the \p buffer.
401  * \param[in,out] buffer          Pointer to the device-side buffer.
402  * \param[in]     startingOffset  Offset (in values) at the device-side buffer to start clearing at.
403  * \param[in]     numValues       Number of values to clear.
404  * \param[in]     deviceStream    GPU stream.
405  */
406 template<typename ValueType>
407 void clearDeviceBufferAsync(DeviceBuffer<ValueType>* buffer,
408                             size_t                   startingOffset,
409                             size_t                   numValues,
410                             const DeviceStream&      deviceStream)
411 {
412     if (numValues == 0)
413     {
414         return;
415     }
416     GMX_ASSERT(buffer, "needs a buffer pointer");
417
418     GMX_ASSERT(checkDeviceBuffer(*buffer, startingOffset + numValues),
419                "buffer too small or not initialized");
420
421     cl::sycl::buffer<ValueType>& syclBuffer = *(buffer->buffer_);
422
423     gmx::internal::fillSyclBufferWithNull<ValueType>(
424             syclBuffer, startingOffset, numValues, deviceStream.stream());
425 }
426
427 /*! \brief Create a texture object for an array of type ValueType.
428  *
429  * Creates the device buffer and copies read-only data for an array of type ValueType.
430  * Like OpenCL, does not really do anything with textures, simply creates a buffer
431  * and initializes it.
432  *
433  * \tparam      ValueType      Raw data type.
434  *
435  * \param[out]  deviceBuffer   Device buffer to store data in.
436  * \param[in]   hostBuffer     Host buffer to get date from.
437  * \param[in]   numValues      Number of elements in the buffer.
438  * \param[in]   deviceContext  GPU device context.
439  */
440 template<typename ValueType>
441 void initParamLookupTable(DeviceBuffer<ValueType>* deviceBuffer,
442                           DeviceTexture* /* deviceTexture */,
443                           const ValueType*     hostBuffer,
444                           int                  numValues,
445                           const DeviceContext& deviceContext)
446 {
447     GMX_ASSERT(hostBuffer, "Host buffer should be specified.");
448     GMX_ASSERT(deviceBuffer, "Device buffer should be specified.");
449
450     /* Constructing buffer with cl::sycl::buffer(T* data, size_t size) will take ownership
451      * of this memory region making it unusable, which might lead to side-effects.
452      * On the other hand, cl::sycl::buffer(InputIterator<T> begin, InputIterator<T> end) will
453      * initialize the buffer without affecting ownership of the memory, although
454      * it will consume extra memory on host. */
455     const cl::sycl::property_list bufferProperties{ cl::sycl::property::buffer::context_bound(
456             deviceContext.context()) };
457     deviceBuffer->buffer_.reset(new ClSyclBufferWrapper<ValueType>(
458             hostBuffer, hostBuffer + numValues, bufferProperties));
459 }
460
461 /*! \brief Release the OpenCL device buffer.
462  *
463  * \tparam        ValueType     Raw data type.
464  *
465  * \param[in,out] deviceBuffer  Device buffer to store data in.
466  */
467 template<typename ValueType>
468 void destroyParamLookupTable(DeviceBuffer<ValueType>* deviceBuffer, DeviceTexture& /* deviceTexture */)
469 {
470     deviceBuffer->buffer_.reset(nullptr);
471 }
472
473 #endif // GMX_GPU_UTILS_DEVICEBUFFER_SYCL_H