Skip to content

Commit e33ff6c

Browse files
author
neok-m4700
committed
FIX issue 503
1 parent 8fb3bab commit e33ff6c

File tree

4 files changed

+46
-1
lines changed

4 files changed

+46
-1
lines changed

CMakeLists.txt

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -775,6 +775,7 @@ if(opencoarrays_aware_compiler)
775775
# issues reported by @neok-m4700
776776
add_caf_test(issue-488-multi-dim-cobounds-true 8 issue-488-multi-dim-cobounds true)
777777
add_caf_test(issue-488-multi-dim-cobounds-false 8 issue-488-multi-dim-cobounds false)
778+
add_caf_test(issue-503-multidim-array-broadcast 2 issue-503-multidim-array-broadcast)
778779

779780
# IMAGE FAIL tests
780781
if(NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 7.0.0)

src/mpi/mpi_caf.c

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -6921,7 +6921,7 @@ PREFIX (co_broadcast) (gfc_descriptor_t *a, int source_image, int *stat, char *e
69216921
extent = (a->dim[j]._ubound - a->dim[j].lower_bound + 1);
69226922
stride = a->dim[j]._stride;
69236923
}
6924-
array_offset_sr += (i / extent) * a->dim[rank-1]._stride;
6924+
array_offset_sr += (i / (extent * stride)) * a->dim[rank-1]._stride;
69256925
void *sr = (void *)((char *) a->base_addr
69266926
+ array_offset_sr*GFC_DESCRIPTOR_SIZE (a));
69276927

src/tests/regression/reported/CMakeLists.txt

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -7,3 +7,4 @@ caf_compile_executable(issue-493-coindex-slice issue-493-coindex-slice.f90)
77
caf_compile_executable(issue-422-send issue-422-send.F90)
88
caf_compile_executable(issue-422-send-get issue-422-send-get.F90)
99
caf_compile_executable(issue-488-multi-dim-cobounds issue-488-multi-dim-cobounds.f90)
10+
caf_compile_executable(issue-503-multidim-array-broadcast issue-503-multidim-array-broadcast.f90)
Lines changed: 43 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,43 @@
1+
! caf -o bcast bcast.f90
2+
! cafrun -np 2 ./bcast
3+
! mwe for issue github.com/sourceryinstitute/OpenCoarrays/issues/503
4+
5+
program main
6+
real, allocatable, dimension(:, :, :) :: arr1, arr2
7+
integer :: i, me, nimg
8+
real :: red1, sum1, red2, sum2
9+
10+
allocate(arr1(10, 20, 8))
11+
allocate(arr2(10, 20, 30))
12+
13+
me = this_image()
14+
nimg = num_images()
15+
16+
if (me == 1) then
17+
arr1 = reshape([(i, i=1, size(arr1))], shape(arr1))
18+
arr2 = reshape([(i, i=1, size(arr2))], shape(arr2))
19+
end if
20+
21+
call co_broadcast(arr1, source_image=1)
22+
call co_broadcast(arr2, source_image=1)
23+
24+
sum1 = sum(arr1)
25+
sum2 = sum(arr2)
26+
27+
red1 = sum1
28+
red2 = sum2
29+
30+
call co_sum(red1)
31+
call co_sum(red2)
32+
33+
sync all
34+
print *, me, ' sum1=', sum1, ' red1=', red1
35+
print *, me, ' sum2=', sum2, ' red2=', red2
36+
37+
if (abs(red1 - nimg * sum1) > epsilon(0.) .or. abs(red2 - nimg * sum2) > epsilon(0.)) then
38+
write(*,*) 'Test failed!'
39+
else
40+
write(*,*) 'Test passed.'
41+
end if
42+
deallocate(arr1, arr2)
43+
end program

0 commit comments

Comments
 (0)