Skip to content

Commit 0e2f624

Browse files
authored
Merge branch 'master' into stop-quiet
2 parents f703ab3 + 8777474 commit 0e2f624

File tree

6 files changed

+53
-31
lines changed

6 files changed

+53
-31
lines changed

CMakeLists.txt

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -657,8 +657,8 @@ if(opencoarrays_aware_compiler)
657657
add_caf_test(comp_allocated_2 2 comp_allocated_2)
658658
add_caf_test(alloc_comp_get_convert_nums 2 alloc_comp_get_convert_nums)
659659
if(NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 8)
660-
add_caf_test(team_number 2 ${tests_root}/unit/teams/team_number)
661-
add_caf_test(get_communicator 3 ${tests_root}/unit/teams/get_communicator)
660+
add_caf_test(team_number 2 team_number)
661+
add_caf_test(get_communicator 3 get_communicator)
662662
endif()
663663
endif()
664664

src/mpi/mpi_caf.c

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

4632-
/* Only increase the dim counter, when in an array ref, and
4633-
MODE != CAF_ARR_REF_SINGLE (delta == 1) see caf_array_ref_t. */
4634-
if (in_array_ref && dst_cur_dim < dst_rank && delta != 1)
4635-
++dst_cur_dim;
4632+
/* Only increase the dim counter, when in an array ref */
4633+
if (in_array_ref && dst_cur_dim < dst_rank)
4634+
{
4635+
// Mode != CAF_ARR_REF_SINGLE (delta == 1), and no rank reduction
4636+
if (!(delta == 1 && dst_rank != GFC_DESCRIPTOR_RANK(src)))
4637+
++dst_cur_dim;
4638+
}
46364639
}
46374640
size *= (ptrdiff_t)delta;
46384641
}
@@ -4778,10 +4781,13 @@ PREFIX (get_by_ref) (caf_token_t token, int image_index,
47784781
dst->dim[dst_cur_dim]._stride = size;
47794782
}
47804783
}
4781-
/* Only increase the dim counter, when in an array ref, and
4782-
MODE != CAF_ARR_REF_SINGLE (delta == 1) see caf_array_ref_t. */
4783-
if (in_array_ref && dst_cur_dim < dst_rank && delta != 1)
4784-
++dst_cur_dim;
4784+
/* Only increase the dim counter, when in an array ref */
4785+
if (in_array_ref && dst_cur_dim < dst_rank)
4786+
{
4787+
// Mode != CAF_ARR_REF_SINGLE (delta == 1), and no rank reduction
4788+
if (!(delta == 1 && dst_rank != GFC_DESCRIPTOR_RANK(src)))
4789+
++dst_cur_dim;
4790+
}
47854791
}
47864792
size *= (ptrdiff_t)delta;
47874793
}
@@ -5352,7 +5358,7 @@ PREFIX (send_by_ref) (caf_token_t token, int image_index,
53525358
"rank out of range.\n";
53535359
const char extentoutofrange[] = "libcaf_mpi::caf_send_by_ref(): "
53545360
"extent out of range.\n";
5355-
const char cannotallocdst[] = "libcaf_mpi::caf_get_by_ref(): "
5361+
const char cannotallocdst[] = "libcaf_mpi::caf_send_by_ref(): "
53565362
"can not allocate %d bytes of memory.\n";
53575363
const char unabletoallocdst[] = "libcaf_mpi::caf_send_by_ref(): "
53585364
"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

src/tests/unit/teams/CMakeLists.txt

Lines changed: 2 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -1,9 +1,2 @@
1-
add_executable(team_number team-number.f90)
2-
target_link_libraries(team_number OpenCoarrays)
3-
#TODO: Remove when gfortran supports team_number()
4-
target_include_directories(team_number PRIVATE ${CMAKE_BINARY_DIR}/include)
5-
6-
#add_executable(get_communicator get-communicator.f90)
7-
#target_link_libraries(get_communicator OpenCoarrays)
8-
#TODO: Remove when gfortran supports team_number()
9-
#target_include_directories(get_communicator PRIVATE ${CMAKE_BINARY_DIR}/include)
1+
caf_compile_executable(team_number team-number.f90)
2+
caf_compile_executable(get_communicator get-communicator.f90)

src/tests/unit/teams/get-communicator.f90

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
! BSD 3-Clause License
22
!
3-
! Copyright (c) 2016, Sourcery Institute
3+
! Copyright (c) 2018, Sourcery Institute
44
! All rights reserved.
55
!
66
! Redistribution and use in source and binary forms, with or without
@@ -83,6 +83,7 @@ pure function destination_team(image,numTeams) result(team)
8383
end function
8484

8585
subroutine mpi_matches_caf(comm)
86+
use mpi
8687
use iso_c_binding, only : c_int
8788
integer(c_int), intent(in) :: comm
8889
!! MPI communicator

src/tests/unit/teams/team-number.f90

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
! BSD 3-Clause License
22
!
3-
! Copyright (c) 2016, Sourcery Institute
3+
! Copyright (c) 2018, Sourcery Institute
44
! All rights reserved.
55
!
66
! Redistribution and use in source and binary forms, with or without

0 commit comments

Comments
 (0)