DPNP C++ backend kernel library 0.20.0dev4
Data Parallel Extension for NumPy*
Loading...
Searching...
No Matches
common.hpp
1//*****************************************************************************
2// Copyright (c) 2023, Intel Corporation
3// All rights reserved.
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7// - Redistributions of source code must retain the above copyright notice,
8// this list of conditions and the following disclaimer.
9// - Redistributions in binary form must reproduce the above copyright notice,
10// this list of conditions and the following disclaimer in the documentation
11// and/or other materials provided with the distribution.
12// - Neither the name of the copyright holder nor the names of its contributors
13// may be used to endorse or promote products derived from this software
14// without specific prior written permission.
15//
16// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
17// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
18// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
19// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
20// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
21// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
22// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
23// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
24// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
25// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF
26// THE POSSIBILITY OF SUCH DAMAGE.
27//*****************************************************************************
28
29#pragma once
30
31#include <cstddef>
32#include <type_traits>
33#include <utility>
34#include <vector>
35
36#include <oneapi/mkl.hpp>
37#include <sycl/sycl.hpp>
38
39#include <dpctl4pybind11.hpp>
40#include <pybind11/pybind11.h>
41
42// utils extension header
43#include "ext/common.hpp"
44
45// dpctl tensor headers
46#include "utils/memory_overlap.hpp"
47#include "utils/type_dispatch.hpp"
48
55#ifndef __INTEL_MKL_2023_2_0_VERSION_REQUIRED
56#define __INTEL_MKL_2023_2_0_VERSION_REQUIRED 20230002L
57#endif
58
59static_assert(INTEL_MKL_VERSION >= __INTEL_MKL_2023_2_0_VERSION_REQUIRED,
60 "OneMKL does not meet minimum version requirement");
61
62namespace ext_ns = ext::common;
63namespace py = pybind11;
64namespace td_ns = dpctl::tensor::type_dispatch;
65
66namespace dpnp::extensions::vm::py_internal
67{
68template <typename output_typesT, typename contig_dispatchT>
69bool need_to_call_unary_ufunc(sycl::queue &exec_q,
70 const dpctl::tensor::usm_ndarray &src,
71 const dpctl::tensor::usm_ndarray &dst,
72 const output_typesT &output_type_vec,
73 const contig_dispatchT &contig_dispatch_vector)
74{
75 // check type_nums
76 int src_typenum = src.get_typenum();
77 int dst_typenum = dst.get_typenum();
78
79 auto array_types = td_ns::usm_ndarray_types();
80 int src_typeid = array_types.typenum_to_lookup_id(src_typenum);
81 int dst_typeid = array_types.typenum_to_lookup_id(dst_typenum);
82
83 // check that types are supported
84 int func_output_typeid = output_type_vec[src_typeid];
85 if (dst_typeid != func_output_typeid) {
86 return false;
87 }
88
89 // OneMKL VM functions perform a copy on host if no double type support
90 if (!exec_q.get_device().has(sycl::aspect::fp64)) {
91 return false;
92 }
93
94 // check that queues are compatible
95 if (!dpctl::utils::queues_are_compatible(exec_q, {src, dst})) {
96 return false;
97 }
98
99 // dimensions must be the same
100 int dst_nd = dst.get_ndim();
101 if (dst_nd != src.get_ndim()) {
102 return false;
103 }
104 else if (dst_nd == 0) {
105 // don't call OneMKL for 0d arrays
106 return false;
107 }
108
109 // shapes must be the same
110 const py::ssize_t *src_shape = src.get_shape_raw();
111 const py::ssize_t *dst_shape = dst.get_shape_raw();
112 bool shapes_equal(true);
113 size_t src_nelems(1);
114
115 for (int i = 0; i < dst_nd; ++i) {
116 src_nelems *= static_cast<size_t>(src_shape[i]);
117 shapes_equal = shapes_equal && (src_shape[i] == dst_shape[i]);
118 }
119 if (!shapes_equal) {
120 return false;
121 }
122
123 // if nelems is zero, return false
124 if (src_nelems == 0) {
125 return false;
126 }
127
128 // ensure that output is ample enough to accommodate all elements
129 auto dst_offsets = dst.get_minmax_offsets();
130 // destination must be ample enough to accommodate all elements
131 {
132 size_t range =
133 static_cast<size_t>(dst_offsets.second - dst_offsets.first);
134 if (range + 1 < src_nelems) {
135 return false;
136 }
137 }
138
139 // check memory overlap
140 auto const &overlap = dpctl::tensor::overlap::MemoryOverlap();
141 if (overlap(src, dst)) {
142 return false;
143 }
144
145 // support only contiguous inputs
146 bool is_src_c_contig = src.is_c_contiguous();
147 bool is_dst_c_contig = dst.is_c_contiguous();
148
149 bool all_c_contig = (is_src_c_contig && is_dst_c_contig);
150 if (!all_c_contig) {
151 return false;
152 }
153
154 // MKL function is not defined for the type
155 if (contig_dispatch_vector[src_typeid] == nullptr) {
156 return false;
157 }
158 return true;
159}
160
161template <typename output_typesT, typename contig_dispatchT>
162bool need_to_call_unary_two_outputs_ufunc(
163 sycl::queue &exec_q,
164 const dpctl::tensor::usm_ndarray &src,
165 const dpctl::tensor::usm_ndarray &dst1,
166 const dpctl::tensor::usm_ndarray &dst2,
167 const output_typesT &output_type_vec,
168 const contig_dispatchT &contig_dispatch_vector)
169{
170 // check type_nums
171 int src_typenum = src.get_typenum();
172 int dst1_typenum = dst1.get_typenum();
173 int dst2_typenum = dst2.get_typenum();
174
175 const auto &array_types = td_ns::usm_ndarray_types();
176 int src_typeid = array_types.typenum_to_lookup_id(src_typenum);
177 int dst1_typeid = array_types.typenum_to_lookup_id(dst1_typenum);
178 int dst2_typeid = array_types.typenum_to_lookup_id(dst2_typenum);
179
180 std::pair<int, int> func_output_typeids = output_type_vec[src_typeid];
181
182 // check that types are supported
183 if (dst1_typeid != func_output_typeids.first ||
184 dst2_typeid != func_output_typeids.second) {
185 return false;
186 }
187
188 // OneMKL VM functions perform a copy on host if no double type support
189 if (!exec_q.get_device().has(sycl::aspect::fp64)) {
190 return false;
191 }
192
193 // check that queues are compatible
194 if (!dpctl::utils::queues_are_compatible(exec_q, {src, dst1, dst2})) {
195 return false;
196 }
197
198 // dimensions must be the same
199 int src_nd = src.get_ndim();
200 int dst1_nd = dst1.get_ndim();
201 int dst2_nd = dst2.get_ndim();
202 if (src_nd != dst1_nd || src_nd != dst2_nd) {
203 return false;
204 }
205 else if (dst1_nd == 0 || dst2_nd == 0) {
206 // don't call OneMKL for 0d arrays
207 return false;
208 }
209
210 // shapes must be the same
211 const py::ssize_t *src_shape = src.get_shape_raw();
212 const py::ssize_t *dst1_shape = dst1.get_shape_raw();
213 const py::ssize_t *dst2_shape = dst2.get_shape_raw();
214 bool shapes_equal(true);
215 size_t src_nelems(1);
216
217 for (int i = 0; i < src_nd; ++i) {
218 src_nelems *= static_cast<std::size_t>(src_shape[i]);
219 shapes_equal = shapes_equal && (src_shape[i] == dst1_shape[i]) &&
220 (src_shape[i] == dst2_shape[i]);
221 }
222 if (!shapes_equal) {
223 return false;
224 }
225
226 // if nelems is zero, return false
227 if (src_nelems == 0) {
228 return false;
229 }
230
231 // ensure that outputs are ample enough to accommodate all elements
232 auto dst1_offsets = dst1.get_minmax_offsets();
233 auto dst2_offsets = dst2.get_minmax_offsets();
234 // destinations must be ample enough to accommodate all elements
235 {
236 size_t range1 =
237 static_cast<size_t>(dst1_offsets.second - dst1_offsets.first);
238 size_t range2 =
239 static_cast<size_t>(dst2_offsets.second - dst2_offsets.first);
240 if ((range1 + 1 < src_nelems) || (range2 + 1 < src_nelems)) {
241 return false;
242 }
243 }
244
245 // check memory overlap
246 auto const &overlap = dpctl::tensor::overlap::MemoryOverlap();
247 if (overlap(src, dst1) || overlap(src, dst2) || overlap(dst1, dst2)) {
248 return false;
249 }
250
251 // support only contiguous inputs
252 bool is_src_c_contig = src.is_c_contiguous();
253 bool is_dst1_c_contig = dst1.is_c_contiguous();
254 bool is_dst2_c_contig = dst2.is_c_contiguous();
255
256 bool all_c_contig =
257 (is_src_c_contig && is_dst1_c_contig && is_dst2_c_contig);
258 if (!all_c_contig) {
259 return false;
260 }
261
262 // MKL function is not defined for the type
263 if (contig_dispatch_vector[src_typeid] == nullptr) {
264 return false;
265 }
266 return true;
267}
268
269template <typename output_typesT, typename contig_dispatchT>
270bool need_to_call_binary_ufunc(sycl::queue &exec_q,
271 const dpctl::tensor::usm_ndarray &src1,
272 const dpctl::tensor::usm_ndarray &src2,
273 const dpctl::tensor::usm_ndarray &dst,
274 const output_typesT &output_type_table,
275 const contig_dispatchT &contig_dispatch_table)
276{
277 // check type_nums
278 int src1_typenum = src1.get_typenum();
279 int src2_typenum = src2.get_typenum();
280 int dst_typenum = dst.get_typenum();
281
282 auto array_types = td_ns::usm_ndarray_types();
283 int src1_typeid = array_types.typenum_to_lookup_id(src1_typenum);
284 int src2_typeid = array_types.typenum_to_lookup_id(src2_typenum);
285 int dst_typeid = array_types.typenum_to_lookup_id(dst_typenum);
286
287 // check that types are supported
288 int output_typeid = output_type_table[src1_typeid][src2_typeid];
289 if (output_typeid != dst_typeid) {
290 return false;
291 }
292
293 // types must be the same
294 if (src1_typeid != src2_typeid) {
295 return false;
296 }
297
298 // OneMKL VM functions perform a copy on host if no double type support
299 if (!exec_q.get_device().has(sycl::aspect::fp64)) {
300 return false;
301 }
302
303 // check that queues are compatible
304 if (!dpctl::utils::queues_are_compatible(exec_q, {src1, src2, dst})) {
305 return false;
306 }
307
308 // dimensions must be the same
309 int dst_nd = dst.get_ndim();
310 if (dst_nd != src1.get_ndim() || dst_nd != src2.get_ndim()) {
311 return false;
312 }
313 else if (dst_nd == 0) {
314 // don't call OneMKL for 0d arrays
315 return false;
316 }
317
318 // shapes must be the same
319 const py::ssize_t *src1_shape = src1.get_shape_raw();
320 const py::ssize_t *src2_shape = src2.get_shape_raw();
321 const py::ssize_t *dst_shape = dst.get_shape_raw();
322 bool shapes_equal(true);
323 size_t src_nelems(1);
324
325 for (int i = 0; i < dst_nd; ++i) {
326 src_nelems *= static_cast<size_t>(src1_shape[i]);
327 shapes_equal = shapes_equal && (src1_shape[i] == dst_shape[i] &&
328 src2_shape[i] == dst_shape[i]);
329 }
330 if (!shapes_equal) {
331 return false;
332 }
333
334 // if nelems is zero, return false
335 if (src_nelems == 0) {
336 return false;
337 }
338
339 // ensure that output is ample enough to accommodate all elements
340 auto dst_offsets = dst.get_minmax_offsets();
341 // destination must be ample enough to accommodate all elements
342 {
343 size_t range =
344 static_cast<size_t>(dst_offsets.second - dst_offsets.first);
345 if (range + 1 < src_nelems) {
346 return false;
347 }
348 }
349
350 // check memory overlap
351 auto const &overlap = dpctl::tensor::overlap::MemoryOverlap();
352 if (overlap(src1, dst) || overlap(src2, dst)) {
353 return false;
354 }
355
356 // support only contiguous inputs
357 bool is_src1_c_contig = src1.is_c_contiguous();
358 bool is_src2_c_contig = src2.is_c_contiguous();
359 bool is_dst_c_contig = dst.is_c_contiguous();
360
361 bool all_c_contig =
362 (is_src1_c_contig && is_src2_c_contig && is_dst_c_contig);
363 if (!all_c_contig) {
364 return false;
365 }
366
367 // MKL function is not defined for the type
368 if (contig_dispatch_table[src1_typeid] == nullptr) {
369 return false;
370 }
371 return true;
372}
373
379#define MACRO_POPULATE_DISPATCH_VECTORS(__name__) \
380 template <typename fnT, typename T> \
381 struct ContigFactory \
382 { \
383 fnT get() \
384 { \
385 if constexpr (std::is_same_v<typename OutputType<T>::value_type, \
386 void>) { \
387 return nullptr; \
388 } \
389 else { \
390 return __name__##_contig_impl<T>; \
391 } \
392 } \
393 }; \
394 \
395 template <typename fnT, typename T> \
396 struct TypeMapFactory \
397 { \
398 std::enable_if_t<std::is_same<fnT, int>::value, int> get() \
399 { \
400 using rT = typename OutputType<T>::value_type; \
401 return td_ns::GetTypeid<rT>{}.get(); \
402 } \
403 }; \
404 \
405 static void populate_dispatch_vectors(void) \
406 { \
407 ext_ns::init_dispatch_vector<int, TypeMapFactory>( \
408 output_typeid_vector); \
409 ext_ns::init_dispatch_vector<unary_contig_impl_fn_ptr_t, \
410 ContigFactory>(contig_dispatch_vector); \
411 };
412
418#define MACRO_POPULATE_DISPATCH_2OUTS_VECTORS(__name__) \
419 template <typename fnT, typename T> \
420 struct ContigFactory \
421 { \
422 fnT get() \
423 { \
424 if constexpr (std::is_same_v<typename OutputType<T>::value_type1, \
425 void> || \
426 std::is_same_v<typename OutputType<T>::value_type2, \
427 void>) { \
428 fnT fn = nullptr; \
429 return fn; \
430 } \
431 else { \
432 fnT fn = __name__##_contig_impl<T>; \
433 return fn; \
434 } \
435 } \
436 }; \
437 \
438 template <typename fnT, typename T> \
439 struct TypeMapFactory \
440 { \
441 std::enable_if_t<std::is_same<fnT, std::pair<int, int>>::value, \
442 std::pair<int, int>> \
443 get() \
444 { \
445 using rT1 = typename OutputType<T>::value_type1; \
446 using rT2 = typename OutputType<T>::value_type2; \
447 return std::make_pair(td_ns::GetTypeid<rT1>{}.get(), \
448 td_ns::GetTypeid<rT2>{}.get()); \
449 } \
450 }; \
451 \
452 static void populate_dispatch_vectors(void) \
453 { \
454 ext_ns::init_dispatch_vector<std::pair<int, int>, TypeMapFactory>( \
455 output_typeid_vector); \
456 ext_ns::init_dispatch_vector<unary_two_outputs_contig_impl_fn_ptr_t, \
457 ContigFactory>(contig_dispatch_vector); \
458 };
459
465#define MACRO_POPULATE_DISPATCH_TABLES(__name__) \
466 template <typename fnT, typename T1, typename T2> \
467 struct ContigFactory \
468 { \
469 fnT get() \
470 { \
471 if constexpr (std::is_same_v< \
472 typename OutputType<T1, T2>::value_type, \
473 void>) { \
474 return nullptr; \
475 } \
476 else { \
477 return __name__##_contig_impl<T1, T2>; \
478 } \
479 } \
480 }; \
481 \
482 template <typename fnT, typename T1, typename T2> \
483 struct TypeMapFactory \
484 { \
485 std::enable_if_t<std::is_same<fnT, int>::value, int> get() \
486 { \
487 using rT = typename OutputType<T1, T2>::value_type; \
488 return td_ns::GetTypeid<rT>{}.get(); \
489 } \
490 }; \
491 \
492 static void populate_dispatch_tables(void) \
493 { \
494 ext_ns::init_dispatch_table<int, TypeMapFactory>( \
495 output_typeid_vector); \
496 ext_ns::init_dispatch_table<binary_contig_impl_fn_ptr_t, \
497 ContigFactory>(contig_dispatch_vector); \
498 };
499} // namespace dpnp::extensions::vm::py_internal