DPNP C++ backend kernel library 0.18.0rc1
Data Parallel Extension for NumPy*
Loading...
Searching...
No Matches
interpolate.hpp
1//*****************************************************************************
2// Copyright (c) 2025, 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//
13// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
14// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
15// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
16// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
17// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
18// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
19// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
20// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
21// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
22// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF
23// THE POSSIBILITY OF SUCH DAMAGE.
24//*****************************************************************************
25
26#pragma once
27
28#include <sycl/sycl.hpp>
29#include <vector>
30
31#include "ext/common.hpp"
32
34
35namespace dpnp::kernels::interpolate
36{
37template <typename TCoord, typename TValue, typename TIdx = std::int64_t>
38sycl::event interpolate_impl(sycl::queue &q,
39 const TCoord *x,
40 const TIdx *idx,
41 const TCoord *xp,
42 const TValue *fp,
43 const TValue *left,
44 const TValue *right,
45 TValue *out,
46 const std::size_t n,
47 const std::size_t xp_size,
48 const std::vector<sycl::event> &depends)
49{
50 // Selected over the work-group version
51 // due to simpler execution and slightly better performance.
52 return q.submit([&](sycl::handler &h) {
53 h.depends_on(depends);
54 h.parallel_for(sycl::range<1>(n), [=](sycl::id<1> i) {
55 TValue left_val = left ? *left : fp[0];
56 TValue right_val = right ? *right : fp[xp_size - 1];
57
58 TCoord x_val = x[i];
59 TIdx x_idx = idx[i] - 1;
60
61 if (IsNan<TCoord>::isnan(x_val)) {
62 out[i] = x_val;
63 }
64 else if (x_idx < 0) {
65 out[i] = left_val;
66 }
67 else if (x_val == xp[xp_size - 1]) {
68 out[i] = fp[xp_size - 1];
69 }
70 else if (x_idx >= static_cast<TIdx>(xp_size - 1)) {
71 out[i] = right_val;
72 }
73 else {
74 TValue slope =
75 (fp[x_idx + 1] - fp[x_idx]) / (xp[x_idx + 1] - xp[x_idx]);
76 TValue res = slope * (x_val - xp[x_idx]) + fp[x_idx];
77
78 if (IsNan<TValue>::isnan(res)) {
79 res = slope * (x_val - xp[x_idx + 1]) + fp[x_idx + 1];
80 if (IsNan<TValue>::isnan(res) &&
81 (fp[x_idx] == fp[x_idx + 1])) {
82 res = fp[x_idx];
83 }
84 }
85 out[i] = res;
86 }
87 });
88 });
89}
90
91} // namespace dpnp::kernels::interpolate