Skip to content

Commit db486b9

Browse files
authored
Fix diagonal() strides calculation for empty results (#2814)
Previously `dpnp.diagonal()` returned incorrect strides when the diagonal was empty (e.g., when offset >= array width). The stride was set to itemsize instead of the correct `st_n + st_m`. This PR: * fixes stride calculation to consistently use (st_n + st_m) * simplifies implementation using unified formula to calculate resulting shape, strides, and offset * extends tests coverage including scenarios with empty diagonals, views, and non-contiguous arrays The PR closes #2761.
1 parent 54d2109 commit db486b9

3 files changed

Lines changed: 77 additions & 21 deletions

File tree

CHANGELOG.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -77,6 +77,7 @@ Also, that release drops support for Python 3.9, making Python 3.10 the minimum
7777
* Resolved an issue where `dpnp` always returns the base allocation pointer, when the view start is expected [#2651](https://github.com/IntelPython/dpnp/pull/2651)
7878
* Fixed an issue causing an exception in `dpnp.geomspace` and `dpnp.logspace` when called with explicit `device` keyword but any input array is allocated on another device [#2723](https://github.com/IntelPython/dpnp/pull/2723)
7979
* Fixed `.data.ptr` property on array views to correctly return the pointer to the view's data location instead of the base allocation pointer [#2812](https://github.com/IntelPython/dpnp/pull/2812)
80+
* Resolved an issue with strides calculation in `dpnp.diagonal` to return correct values for empty diagonals [#2814](https://github.com/IntelPython/dpnp/pull/2814)
8081

8182
### Security
8283

dpnp/dpnp_iface_indexing.py

Lines changed: 11 additions & 17 deletions
Original file line numberDiff line numberDiff line change
@@ -719,24 +719,18 @@ def diagonal(a, offset=0, axis1=0, axis2=1):
719719
offset = -offset
720720

721721
a_shape = a.shape
722-
a_straides = a.strides
722+
a_strides = a.strides
723723
n, m = a_shape[-2:]
724-
st_n, st_m = a_straides[-2:]
725-
726-
# Compute shape, strides and offset of the resulting diagonal array
727-
# based on the input offset
728-
if offset == 0:
729-
out_shape = a_shape[:-2] + (min(n, m),)
730-
out_strides = a_straides[:-2] + (st_n + st_m,)
731-
out_offset = 0
732-
elif 0 < offset < m:
733-
out_shape = a_shape[:-2] + (min(n, m - offset),)
734-
out_strides = a_straides[:-2] + (st_n + st_m,)
735-
out_offset = st_m // a.itemsize * offset
736-
else:
737-
out_shape = a_shape[:-2] + (0,)
738-
out_strides = a_straides[:-2] + (a.itemsize,)
739-
out_offset = 0
724+
st_n, st_m = a_strides[-2:]
725+
726+
# Compute the diagonal array as a view:
727+
# - stride: sum of row and column strides (diag advances in both dimensions)
728+
# - shape: determined by diagonal size using max(0, min(n, m - offset))
729+
# - offset: starting position in buffer for non-zero offsets
730+
diag_size = max(0, min(n, m - offset))
731+
out_shape = a_shape[:-2] + (diag_size,)
732+
out_strides = a_strides[:-2] + (st_n + st_m,)
733+
out_offset = st_m // a.itemsize * offset
740734

741735
return dpnp_array(
742736
out_shape, buffer=a, strides=out_strides, offset=out_offset

dpnp/tests/test_indexing.py

Lines changed: 65 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -18,6 +18,7 @@
1818
from dpnp.exceptions import AxisError, ExecutionPlacementError
1919

2020
from .helper import (
21+
generate_random_numpy_array,
2122
get_abs_array,
2223
get_all_dtypes,
2324
get_array,
@@ -44,7 +45,9 @@ def wrapped(a, axis, **kwargs):
4445

4546

4647
class TestDiagonal:
47-
@pytest.mark.parametrize("dtype", get_all_dtypes(no_bool=True))
48+
@pytest.mark.parametrize(
49+
"dtype", get_all_dtypes(no_none=True, no_bool=True)
50+
)
4851
@pytest.mark.parametrize("offset", [-3, -1, 0, 1, 3])
4952
@pytest.mark.parametrize(
5053
"shape",
@@ -58,7 +61,7 @@ class TestDiagonal:
5861
"(2, 2, 2, 3)",
5962
],
6063
)
61-
def test_diagonal_offset(self, shape, dtype, offset):
64+
def test_offset(self, shape, dtype, offset):
6265
a = numpy.arange(numpy.prod(shape), dtype=dtype).reshape(shape)
6366
a_dp = dpnp.array(a)
6467
expected = numpy.diagonal(a, offset)
@@ -74,7 +77,7 @@ def test_diagonal_offset(self, shape, dtype, offset):
7477
((4, 3, 5, 2), [(0, 1), (1, 2), (2, 3), (0, 3)]),
7578
],
7679
)
77-
def test_diagonal_axes(self, shape, axis_pairs, dtype):
80+
def test_axes(self, shape, axis_pairs, dtype):
7881
a = numpy.arange(numpy.prod(shape), dtype=dtype).reshape(shape)
7982
a_dp = dpnp.array(a)
8083
for axis1, axis2 in axis_pairs:
@@ -91,7 +94,7 @@ def test_linalg_diagonal(self, offset):
9194
result = dpnp.linalg.diagonal(a_dp, offset=offset)
9295
assert_array_equal(expected, result)
9396

94-
def test_diagonal_errors(self):
97+
def test_errors(self):
9598
a = dpnp.arange(12).reshape(3, 4)
9699

97100
# unsupported type
@@ -115,6 +118,64 @@ def test_diagonal_errors(self):
115118
assert_raises(ValueError, a.diagonal, axis1=1, axis2=1)
116119
assert_raises(ValueError, a.diagonal, axis1=1, axis2=-1)
117120

121+
@pytest.mark.parametrize("dt", get_all_dtypes(no_none=True))
122+
@pytest.mark.parametrize(
123+
"shape, offset",
124+
[
125+
((2, 5), 5), # offset >= m
126+
((2, 5), 10), # offset >> m
127+
((4, 5), 6), # offset >= m
128+
((2, 5), -5), # negative offset >= n
129+
((3, 3, 4), 5), # 3D array, offset >= m
130+
],
131+
)
132+
def test_empty_strides(self, dt, shape, offset):
133+
a = generate_random_numpy_array(shape=shape, dtype=dt)
134+
ia = dpnp.array(a)
135+
136+
expected = numpy.diagonal(a, offset)
137+
result = dpnp.diagonal(ia, offset)
138+
139+
# Check both shape and strides match NumPy
140+
assert expected.shape == result.shape
141+
assert expected.strides == result.strides
142+
assert_array_equal(expected, result)
143+
144+
@pytest.mark.parametrize("dt", get_all_dtypes(no_none=True))
145+
def test_view(self, dt):
146+
a = generate_random_numpy_array(shape=(3, 4), dtype=dt)
147+
a = dpnp.array(a)
148+
ia = a.copy()
149+
150+
diag = dpnp.diagonal(a)
151+
diag[1] = 17 # modify a diagonal element
152+
ia[1, 1] = 17 # do the same in original copy of the array
153+
154+
assert (a == ia).all()
155+
156+
@pytest.mark.parametrize("dt", get_all_dtypes(no_none=True))
157+
@pytest.mark.parametrize(
158+
"slice_spec, offset",
159+
[
160+
((slice(None), slice(None, None, 2)), 0), # skip columns
161+
((slice(None, None, 2), slice(None)), 1), # skip rows
162+
((slice(None, None, 2), slice(None, None, 2)), 0), # skip both
163+
],
164+
)
165+
def test_noncontiguous(self, dt, slice_spec, offset):
166+
a = generate_random_numpy_array(shape=(4, 6), dtype=dt)
167+
a_sliced = a[slice_spec]
168+
ia = dpnp.array(a)
169+
ia_sliced = ia[slice_spec]
170+
171+
expected = numpy.diagonal(a_sliced, offset=offset)
172+
result = dpnp.diagonal(ia_sliced, offset=offset)
173+
174+
# Check strides match for non-contiguous arrays
175+
assert expected.shape == result.shape
176+
assert expected.strides == result.strides
177+
assert_array_equal(expected, result)
178+
118179

119180
class TestExtins:
120181
@pytest.mark.parametrize("dt", get_all_dtypes(no_none=True))

0 commit comments

Comments
 (0)