Skip to content

Commit 342421c

Browse files
author
Damian Rouson
authored
Merge pull request #741 from sourceryinstitute/issue-739-fix-seg-fault
Issue 739 fix seg fault
2 parents bad7f50 + a92ba10 commit 342421c

File tree

4 files changed

+33
-9
lines changed

4 files changed

+33
-9
lines changed

CMakeLists.txt

+1
Original file line numberDiff line numberDiff line change
@@ -795,6 +795,7 @@ if(opencoarrays_aware_compiler)
795795
add_caf_test(get_with_offset_1d 2 get_with_offset_1d)
796796
add_caf_test(whole_get_array 2 whole_get_array)
797797
add_caf_test(strided_get 2 strided_get)
798+
add_caf_test(get_static_array 2 get_static_array)
798799

799800
# Pure send tests
800801
add_caf_test(send_array 2 send_array)

src/mpi/mpi_caf.c

+8-9
Original file line numberDiff line numberDiff line change
@@ -4673,8 +4673,8 @@ PREFIX(get_by_ref) (caf_token_t token, int image_index,
46734673
size = 1;
46744674
while (riter)
46754675
{
4676-
dprint("caf_ref = %p, offset = %zd, remote_mem = %p, global_win(data, desc)) = (%d, %d)\n",
4677-
riter, data_offset, remote_memptr, access_data_through_global_win,
4676+
dprint("caf_ref = %p, type = %d, offset = %zd, remote_mem = %p, global_win(data, desc)) = (%d, %d)\n",
4677+
riter, riter->type, data_offset, remote_memptr, access_data_through_global_win,
46784678
access_desc_through_global_win);
46794679
switch (riter->type)
46804680
{
@@ -4996,7 +4996,7 @@ case kind: \
49964996
delta = riter->u.a.dim[i].v.nvec;
49974997
#define KINDCASE(kind, type) \
49984998
case kind: \
4999-
remote_memptr += \
4999+
data_offset += \
50005000
((type *)riter->u.a.dim[i].v.vector)[0] * riter->item_size; \
50015001
break
50025002

@@ -5026,15 +5026,14 @@ case kind: \
50265026
riter->u.a.dim[i].s.stride,
50275027
riter->u.a.dim[i].s.start,
50285028
riter->u.a.dim[i].s.end);
5029-
remote_memptr += riter->u.a.dim[i].s.start
5030-
* riter->u.a.dim[i].s.stride
5031-
* riter->item_size;
5029+
data_offset += riter->u.a.dim[i].s.start
5030+
* riter->u.a.dim[i].s.stride
5031+
* riter->item_size;
50325032
break;
50335033
case CAF_ARR_REF_SINGLE:
50345034
delta = 1;
5035-
remote_memptr += riter->u.a.dim[i].s.start
5036-
* riter->u.a.dim[i].s.stride
5037-
* riter->item_size;
5035+
data_offset += riter->u.a.dim[i].s.start
5036+
* riter->item_size;
50385037
break;
50395038
case CAF_ARR_REF_OPEN_END:
50405039
/* This and OPEN_START are mapped to a RANGE and therefore can

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

+1
Original file line numberDiff line numberDiff line change
@@ -13,6 +13,7 @@ caf_compile_executable(strided_get strided_get.f90)
1313
caf_compile_executable(get_with_vector_index get_with_vector_index.f90)
1414
## Inquiry functions (these are gets that could be optimized in the future to communicate only the descriptors)
1515
caf_compile_executable(alloc_comp_multidim_shape alloc_comp_multidim_shape.F90)
16+
caf_compile_executable(get_static_array get_static_array.f90)
1617

1718
## Pure send() tests
1819
caf_compile_executable(send_array send_array_test.f90)
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,23 @@
1+
program get_static_array
2+
type :: container
3+
integer, allocatable :: stuff(:)
4+
end type
5+
6+
type(container) :: co_containers(10)[*]
7+
8+
if (this_image() == 1) then
9+
allocate(co_containers(2)%stuff(4))
10+
co_containers(2)%stuff = [1,2,3,4]
11+
end if
12+
13+
sync all
14+
15+
if (this_image() == 2) then
16+
if (any(co_containers(2)[1]%stuff /= [1,2,3,4])) then
17+
error stop "Test failed."
18+
else
19+
print *, "Test passed."
20+
end if
21+
end if
22+
end program
23+

0 commit comments

Comments
 (0)