Skip to content

Commit 781a084

Browse files
authored
Merge pull request #255 from sourceryinstitute/bugfix-strided_sendget
Bugfix strided sendget Fixes #254
2 parents 4d6a0d7 + 0c792d3 commit 781a084

File tree

4 files changed

+121
-6
lines changed

4 files changed

+121
-6
lines changed

CMakeLists.txt

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -370,6 +370,7 @@ if(opencoarrays_aware_compiler)
370370
add_mpi_test(get_with_offset_1d 2 ${tests_root}/unit/send-get/get_with_offset_1d)
371371
add_mpi_test(whole_get_array 2 ${tests_root}/unit/send-get/whole_get_array)
372372
add_mpi_test(strided_get 2 ${tests_root}/unit/send-get/strided_get)
373+
add_mpi_test(strided_sendget 3 ${tests_root}/unit/send-get/strided_sendget)
373374
add_mpi_test(co_sum 4 ${tests_root}/unit/collectives/co_sum_test)
374375
add_mpi_test(co_broadcast 4 ${tests_root}/unit/collectives/co_broadcast_test)
375376
add_mpi_test(co_min 4 ${tests_root}/unit/collectives/co_min_test)

src/mpi/mpi_caf.c

Lines changed: 8 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -912,37 +912,39 @@ PREFIX (sendget) (caf_token_t token_s, size_t offset_s, int image_index_s,
912912
ptrdiff_t array_offset_dst = 0;
913913
ptrdiff_t stride = 1;
914914
ptrdiff_t extent = 1;
915+
ptrdiff_t tot_ext = 1;
915916
for (j = 0; j < rank-1; j++)
916917
{
917-
array_offset_dst += ((i / (extent*stride))
918+
array_offset_dst += ((i / tot_ext)
918919
% (dest->dim[j]._ubound
919920
- dest->dim[j].lower_bound + 1))
920921
* dest->dim[j]._stride;
921922
extent = (dest->dim[j]._ubound - dest->dim[j].lower_bound + 1);
922923
stride = dest->dim[j]._stride;
924+
tot_ext *= extent;
923925
}
924926

925-
extent = (dest->dim[rank-1]._ubound - dest->dim[rank-1].lower_bound + 1);
926-
array_offset_dst += (i / extent) * dest->dim[rank-1]._stride;
927+
array_offset_dst += (i / tot_ext) * dest->dim[rank-1]._stride;
927928
dst_offset = offset_s + array_offset_dst*GFC_DESCRIPTOR_SIZE (dest);
928929

929930
ptrdiff_t array_offset_sr = 0;
930931
if (GFC_DESCRIPTOR_RANK (src) != 0)
931932
{
932933
stride = 1;
933934
extent = 1;
935+
tot_ext = 1;
934936
for (j = 0; j < GFC_DESCRIPTOR_RANK (src)-1; j++)
935937
{
936-
array_offset_sr += ((i / (extent*stride))
938+
array_offset_sr += ((i / tot_ext)
937939
% (src->dim[j]._ubound
938940
- src->dim[j].lower_bound + 1))
939941
* src->dim[j]._stride;
940942
extent = (src->dim[j]._ubound - src->dim[j].lower_bound + 1);
941943
stride = src->dim[j]._stride;
944+
tot_ext *= extent;
942945
}
943946

944-
extent = (src->dim[rank-1]._ubound - src->dim[rank-1].lower_bound + 1);
945-
array_offset_sr += (i / extent) * src->dim[rank-1]._stride;
947+
array_offset_sr += (i / tot_ext) * src->dim[rank-1]._stride;
946948
array_offset_sr *= GFC_DESCRIPTOR_SIZE (src);
947949
}
948950
src_offset = offset_g + array_offset_sr;

src/tests/unit/send-get/CMakeLists.txt

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -18,3 +18,6 @@ target_link_libraries(whole_get_array OpenCoarrays)
1818

1919
add_executable(strided_get strided_get.f90)
2020
target_link_libraries(strided_get OpenCoarrays)
21+
22+
add_executable(strided_sendget strided_sendget.f90)
23+
target_link_libraries(strided_sendget OpenCoarrays)
Lines changed: 109 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,109 @@
1+
! Test that sendget with strides on either side (of the assignment) works
2+
! as expected.
3+
!
4+
! This test needs at least three images, because sendget has the potential
5+
! to check whether on image used in the communication is the current one.
6+
! More than three images do not pay, because there is no general code in
7+
! this test.
8+
!
9+
! Written by Andre Vehreschild
10+
11+
program stridedsendgettest
12+
13+
implicit none
14+
15+
integer, parameter :: src_image = 2, dst_image = 3, master_image = 1
16+
integer, save, dimension(4,6) :: srcmat[*], dstmat[*]
17+
integer, save, dimension(6) :: srcvec[*], dstvec[*]
18+
integer :: i
19+
logical :: test_passed = .true.
20+
21+
! Make sure that enough images are available for this test.
22+
! Everything less than dst_image == 3 may make sendget use an
23+
! optimized version saving a part of the communication, which is
24+
! not what the test should test.
25+
if (num_images() < dst_image) then
26+
print*, "Pretend that the test was run and passed, even though there are too few images to perform test:"
27+
print*, "Test passed"
28+
error stop "Need at least three images."
29+
end if
30+
31+
! On the src_image, set some defined values, to be able to distinguish
32+
! strides going wrong.
33+
if (this_image() == src_image) then
34+
srcvec = [(2 * i, i = 1, 6)]
35+
srcmat = reshape([(i * 2, i = 1, 4*6)], [4,6])
36+
! On the dst_image set values that enable to recognize unset values.
37+
elseif (this_image() == dst_image) then
38+
dstmat = -1
39+
dstvec = -2
40+
end if
41+
42+
! Make sure data is valid on all images.
43+
sync all
44+
45+
! master_image is the controller in this communication and therefore needs
46+
! to initiate the communication.
47+
if (this_image() == master_image) then
48+
! Transfer data from the src-vector to the dst-vector on image
49+
! dst_image. This is a transfer of a contingous block of data and here for
50+
! completeness only.
51+
dstvec(:)[dst_image] = srcvec(:)[src_image]
52+
! This statement uses a stride in the send phase of the communication.
53+
dstmat(3,:)[dst_image] = srcvec(:)[src_image]
54+
end if
55+
56+
! Make sure the communication has completed.
57+
sync all
58+
59+
! Check the result of communication on the dst_image.
60+
if (this_image() == dst_image) then
61+
! Check that transfering to the vector has succeeded.
62+
if (any(dstvec /= [2, 4, 6, 8, 10, 12])) error stop "SendGet vec/vec does not match."
63+
64+
! Check that transfering a vector into a matrix changes only the
65+
! values desired.
66+
if (any(dstmat /= reshape([-1, -1, 2, -1, &
67+
-1, -1, 4, -1, &
68+
-1, -1, 6, -1, &
69+
-1, -1, 8, -1, &
70+
-1, -1, 10, -1, &
71+
-1, -1, 12, -1], [4, 6]))) then
72+
error stop "SendGet matrow/vec does not match."
73+
end if
74+
! Reset the dst-buffers to enable new test.
75+
dstvec = -2
76+
dstmat = -1
77+
end if
78+
79+
! Wait for dst having done its tests.
80+
sync all
81+
if (this_image() == master_image) then
82+
! Execute strided get in sendget and store in a vector to just
83+
! test the get.
84+
dstvec(:)[dst_image] = srcmat(2,:)[src_image]
85+
! Test both strided get and strided send at once.
86+
dstmat(3,:)[dst_image] = srcmat(2,:)[src_image]
87+
end if
88+
89+
! Ensure that the communication is all done.
90+
sync all
91+
if (this_image() == dst_image) then
92+
! Check, that the strided get has the expected result.
93+
if (any(dstvec /= [4, 12, 20, 28, 36, 44])) error stop "SendGet vec/matrow does not match."
94+
95+
! And that both communications with stride work as expected.
96+
if (any(dstmat /= reshape([-1, -1, 4, -1, &
97+
-1, -1, 12, -1, &
98+
-1, -1, 20, -1, &
99+
-1, -1, 28, -1, &
100+
-1, -1, 36, -1, &
101+
-1, -1, 44, -1], [4, 6]))) then
102+
error stop "SendGet matrow/matrow does not match."
103+
end if
104+
105+
! Above checks would stop with an error on failure, so its save
106+
! to unguardedly print here, when all tests pass.
107+
print *, "Test passed"
108+
end if
109+
end program

0 commit comments

Comments
 (0)