Skip to content

Commit bcd48b4

Browse files
author
neok-m4700
committed
update patch for issue 493 - improved regression scenario
1 parent b3d02c3 commit bcd48b4

File tree

2 files changed

+46
-18
lines changed

2 files changed

+46
-18
lines changed

src/mpi/mpi_caf.c

Lines changed: 15 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -4627,10 +4627,13 @@ PREFIX (get_by_ref) (caf_token_t token, int image_index,
46274627
}
46284628
}
46294629

4630-
/* Only increase the dim counter, when in an array ref, and
4631-
MODE != CAF_ARR_REF_SINGLE (delta == 1) see caf_array_ref_t. */
4632-
if (in_array_ref && dst_cur_dim < dst_rank && delta != 1)
4633-
++dst_cur_dim;
4630+
/* Only increase the dim counter, when in an array ref */
4631+
if (in_array_ref && dst_cur_dim < dst_rank)
4632+
{
4633+
// Mode != CAF_ARR_REF_SINGLE (delta == 1), and no rank reduction
4634+
if (!(delta == 1 && dst_rank != GFC_DESCRIPTOR_RANK(src)))
4635+
++dst_cur_dim;
4636+
}
46344637
}
46354638
size *= (ptrdiff_t)delta;
46364639
}
@@ -4776,10 +4779,13 @@ PREFIX (get_by_ref) (caf_token_t token, int image_index,
47764779
dst->dim[dst_cur_dim]._stride = size;
47774780
}
47784781
}
4779-
/* Only increase the dim counter, when in an array ref, and
4780-
MODE != CAF_ARR_REF_SINGLE (delta == 1) see caf_array_ref_t. */
4781-
if (in_array_ref && dst_cur_dim < dst_rank && delta != 1)
4782-
++dst_cur_dim;
4782+
/* Only increase the dim counter, when in an array ref */
4783+
if (in_array_ref && dst_cur_dim < dst_rank)
4784+
{
4785+
// Mode != CAF_ARR_REF_SINGLE (delta == 1), and no rank reduction
4786+
if (!(delta == 1 && dst_rank != GFC_DESCRIPTOR_RANK(src)))
4787+
++dst_cur_dim;
4788+
}
47834789
}
47844790
size *= (ptrdiff_t)delta;
47854791
}
@@ -5350,7 +5356,7 @@ PREFIX (send_by_ref) (caf_token_t token, int image_index,
53505356
"rank out of range.\n";
53515357
const char extentoutofrange[] = "libcaf_mpi::caf_send_by_ref(): "
53525358
"extent out of range.\n";
5353-
const char cannotallocdst[] = "libcaf_mpi::caf_get_by_ref(): "
5359+
const char cannotallocdst[] = "libcaf_mpi::caf_send_by_ref(): "
53545360
"can not allocate %d bytes of memory.\n";
53555361
const char unabletoallocdst[] = "libcaf_mpi::caf_send_by_ref(): "
53565362
"unable to allocate memory on remote image.\n";

src/tests/regression/reported/issue-493-coindex-slice.f90

Lines changed: 31 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -5,7 +5,9 @@ program slice
55

66
type(coarr) :: co
77
integer :: nimg, me, z, nx, ny, nz, north, south, mex, mey, mez, coords(3)
8-
real, allocatable :: bufxz(:, :) ! a plane (2d) slice, normal in the y direction
8+
integer :: shape2d(2), shape3d(3)
9+
real, allocatable :: buf3d(:, :, :) ! a plane slice as a rank 3 array with a single transverse layer
10+
real, allocatable :: buf2d(:, :) ! a plane (2d) slice, normal in the y direction
911

1012
nx = 6
1113
ny = 4
@@ -17,7 +19,12 @@ program slice
1719
if (nimg /= 8) stop
1820

1921
allocate(co % a(nx, ny, nz)[1:2, 1:2, *])
20-
allocate(bufxz(nx, nz))
22+
allocate(buf2d(nx, nz), buf3d(nx, 1, nz))
23+
24+
! this example should NOT reallocate buf2d nor buf3d
25+
! compare shapes before and after syncing
26+
shape2d = shape(buf2d)
27+
shape3d = shape(buf3d)
2128

2229
co % a = reshape([(z, z=1, nx * ny * nz)], shape(co % a))
2330

@@ -33,19 +40,34 @@ program slice
3340
if (north <= 2) then
3441
z = image_index(co % a, [mex, north, mez])
3542
sync images(z)
36-
bufxz = co % a(1:nx, 1, 1:nz)[mex, north, mez]
37-
co % a(1:nx, ny, 1:nz) = bufxz
43+
! no reduction on rank
44+
buf3d = co % a(1:nx, 1:1, 1:nz)[mex, north, mez]
45+
co % a(1:nx, ny:ny, 1:nz) = buf3d
46+
47+
! reduction along dim 2
48+
buf2d = co % a(1:nx, 1, 1:nz)[mex, north, mez]
49+
co % a(1:nx, ny, 1:nz) = buf2d
3850
end if
3951
if (south >= 1) then
4052
z = image_index(co % a, [mex, south, mez])
4153
sync images(z)
42-
bufxz = co % a(1:nx, 1, 1:nz)[mex, south, mez]
43-
co % a(1:nx, ny, 1:nz) = bufxz
54+
buf3d = co % a(1:nx, ny:ny, 1:nz)[mex, south, mez]
55+
co % a(1:nx, 1:1, 1:nz) = buf3d
56+
57+
buf2d = co % a(1:nx, ny, 1:nz)[mex, south, mez]
58+
co % a(1:nx, 1, 1:nz) = buf2d
4459
end if
4560
sync all
4661

47-
deallocate(co % a, bufxz)
48-
if(this_image() == 1) write(*,*) "Test passed."
62+
deallocate(co % a, buf2d, buf3d)
63+
64+
if (any(shape2d /= shape(buf2d)) .or. any(shape3d /= shape(buf3d))) then
65+
write(*, *) 'Test failed!'
66+
error stop 5
67+
else
68+
write(*, *) 'Test passed.'
69+
end if
70+
4971
! Regression would cause error message:
5072
! Fortran runtime error on image <...>: libcaf_mpi::caf_get_by_ref(): rank out of range.
51-
end program
73+
end program

0 commit comments

Comments
 (0)