Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Issue 762 #763

Merged
merged 3 commits into from
Jul 29, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -925,6 +925,7 @@ if(opencoarrays_aware_compiler)
endif()

add_caf_test(issue-700-allow-multiple-scalar-dim-array-gets 2 issue-700-allow-multiple-scalar-dim-array-gets)
add_caf_test(issue-762-mpi-crashing-on-exit 2 issue-762-mpi-crashing-on-exit)

# IMAGE FAIL tests
if(NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 7.0.0)
Expand Down
26 changes: 17 additions & 9 deletions src/runtime-libraries/mpi/mpi_caf.c
Original file line number Diff line number Diff line change
Expand Up @@ -1034,14 +1034,17 @@ finalize_internal(int status_code)
while (cur_stok)
{
prev_stok = cur_stok->prev;
ierr = MPI_Win_detach(global_dynamic_win, cur_stok); chk_err(ierr);
dprint("freeing slave token %p for memory %p", cur_stok->token, cur_stok->token->memptr);
ierr = MPI_Win_detach(global_dynamic_win, cur_stok->token); chk_err(ierr);
if (cur_stok->token->memptr)
{
ierr = MPI_Win_detach(global_dynamic_win, cur_stok->token->memptr);
chk_err(ierr);
free(cur_stok->token->memptr);
ierr = MPI_Free_mem(cur_stok->token->memptr);
chk_err(ierr);
}
free(cur_stok->token);
ierr = MPI_Free_mem(cur_stok->token);
chk_err(ierr);
free(cur_stok);
cur_stok = prev_stok;
}
Expand All @@ -1064,10 +1067,13 @@ finalize_internal(int status_code)
#ifdef GCC_GE_7
/* Unregister the window to the descriptors when freeing the token. */
dprint("MPI_Win_free(%p)\n", p);
ierr = MPI_Win_free(p); chk_err(ierr);
free(cur_tok->token);
ierr = MPI_Win_free(p);
chk_err(ierr);
ierr = MPI_Free_mem(cur_tok->token);
chk_err(ierr);
#else // GCC_GE_7
ierr = MPI_Win_free(p); chk_err(ierr);
ierr = MPI_Win_free(p);
chk_err(ierr);
#endif // GCC_GE_7
free(cur_tok);
cur_tok = prev;
Expand Down Expand Up @@ -1213,7 +1219,7 @@ PREFIX(register) (size_t size, caf_register_t type, caf_token_t *token,
struct caf_allocated_slave_tokens_t *tmp =
malloc(sizeof(struct caf_allocated_slave_tokens_t));
tmp->prev = caf_allocated_slave_tokens;
tmp->token = *token;
tmp->token = slave_token;
caf_allocated_slave_tokens = tmp;
}
else // (type == CAF_REGTYPE_COARRAY_ALLOC_ALLOCATE_ONLY)
Expand Down Expand Up @@ -1528,7 +1534,8 @@ PREFIX(deregister) (caf_token_t *token, int *stat, char *errmsg,
{
ierr = MPI_Win_detach(global_dynamic_win, slave_token->memptr);
chk_err(ierr);
free(slave_token->memptr);
ierr = MPI_Free_mem(slave_token->memptr);
chk_err(ierr);
slave_token->memptr = NULL;
if (type == CAF_DEREGTYPE_COARRAY_DEALLOCATE_ONLY)
{
Expand All @@ -1545,7 +1552,8 @@ PREFIX(deregister) (caf_token_t *token, int *stat, char *errmsg,
caf_allocated_slave_tokens = prev_stok;

free(cur_stok);
free(*token);
ierr = MPI_Free_mem(*token);
chk_err(ierr);
return;
}

Expand Down
1 change: 1 addition & 0 deletions src/tests/regression/reported/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -16,3 +16,4 @@ caf_compile_executable(issue-700-allow-multiple-scalar-dim-array-gets issue-700-
if (gfortran_compiler AND (NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 7.0.0))
caf_compile_executable(issue-515-mimic-mpi-gatherv issue-515-mimic-mpi-gatherv.f90)
endif()
caf_compile_executable(issue-762-mpi-crashing-on-exit issue-762-mpi-crashing-on-exit.f90)
13 changes: 13 additions & 0 deletions src/tests/regression/reported/issue-762-mpi-crashing-on-exit.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,13 @@
program hello_coarrays
implicit none
type :: array_type
integer, allocatable :: values(:)
end type
type(array_type) :: array[*]
allocate(array%values(2), source=0)
array%values = this_image()
if (all(array%values(:) .EQ. this_image())) then
print *,"Test passed."
end if
end program