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 733 impl random init #740

Merged
merged 6 commits into from
Oct 20, 2021
Merged
Show file tree
Hide file tree
Changes from 3 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
3 changes: 3 additions & 0 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -782,6 +782,9 @@ if(opencoarrays_aware_compiler)
message( AUTHOR_WARNING "Skipping the following test to GFortran < 7.4.0 lack of compatibility:
send-strided-self.f90")
endif()
if((NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 12.0.0) OR (CAF_RUN_DEVELOPER_TESTS OR $ENV{OPENCOARRAYS_DEVELOPER}))
add_caf_test(random_init 4 random_init)
endif()
endif()

# Pure get tests
Expand Down
2 changes: 2 additions & 0 deletions src/libcaf.h
Original file line number Diff line number Diff line change
Expand Up @@ -331,6 +331,8 @@ void PREFIX (event_post) (caf_token_t, size_t, int, int *, char *, charlen_t);
void PREFIX (event_wait) (caf_token_t, size_t, int, int *, char *, charlen_t);
void PREFIX (event_query) (caf_token_t, size_t, int, int *, int *);

void PREFIX (random_init) (bool, bool);

/* Language extension */
#ifdef HAVE_MPI
MPI_Fint PREFIX (get_communicator) (caf_team_t *);
Expand Down
70 changes: 70 additions & 0 deletions src/mpi/mpi_caf.c
Original file line number Diff line number Diff line change
Expand Up @@ -8491,3 +8491,73 @@ void PREFIX(sync_team) (caf_team_t *team , int unused __attribute__((unused)))

int ierr = MPI_Barrier(*tmp_comm); chk_err(ierr);
}

extern void _gfortran_random_seed_i4 (int32_t *size, gfc_dim1_descriptor_t *put,
gfc_dim1_descriptor_t *get);

void PREFIX(random_init) (bool repeatable, bool image_distinct)
{
static gfc_dim1_descriptor_t rand_seed;
static bool rep_needs_init = true, arr_needs_init = true;
static int32_t seed_size;

if (arr_needs_init)
{
_gfortran_random_seed_i4(&seed_size, NULL, NULL);
memset(&rand_seed, 0, sizeof(gfc_dim1_descriptor_t));
rand_seed.base.base_addr = malloc(seed_size * sizeof(int32_t)); // because using seed_i4
rand_seed.base.offset = -1;
rand_seed.base.dtype.elem_len = sizeof(int32_t);
rand_seed.base.dtype.rank = 1;
rand_seed.base.dtype.type = BT_INTEGER;
rand_seed.base.span = 0;
rand_seed.dim[0].lower_bound = 1;
rand_seed.dim[0]._ubound = seed_size;
rand_seed.dim[0]._stride = 1;

arr_needs_init = false;
}

if (repeatable)
{
if (rep_needs_init)
{
int32_t lcg_seed = 57911963;
if (image_distinct)
{
lcg_seed *= caf_this_image;
}
int32_t *curr = rand_seed.base.base_addr;
for (int i = 0; i < seed_size; ++i)
{
const int32_t a = 16087;
const int32_t m = INT32_MAX;
const int32_t q = 127773;
const int32_t r = 2836;
lcg_seed = a * (lcg_seed % q) - r * (lcg_seed / q);
if (lcg_seed <= 0) lcg_seed += m;
*curr = lcg_seed;
++curr;
}
rep_needs_init = false;
}
_gfortran_random_seed_i4(NULL, &rand_seed, NULL);
}
else if (image_distinct)
{
_gfortran_random_seed_i4(NULL, NULL, NULL);
}
else
{
if (caf_this_image == 0)
{
_gfortran_random_seed_i4(NULL, NULL, &rand_seed);
MPI_Bcast(rand_seed.base.base_addr, seed_size, MPI_INT32_T, 0, CAF_COMM_WORLD);
}
else
{
MPI_Bcast(rand_seed.base.base_addr, seed_size, MPI_INT32_T, 0, CAF_COMM_WORLD);
_gfortran_random_seed_i4(NULL, &rand_seed, NULL);
}
}
}
1 change: 1 addition & 0 deletions src/tests/unit/simple/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@

caf_compile_executable(increment_my_neighbor increment_neighbor.f90)
caf_compile_executable(atomics testAtomics.f90)
caf_compile_executable(random_init random_init.f90)

# C tests
#include(CMakeForceCompiler)
Expand Down
124 changes: 124 additions & 0 deletions src/tests/unit/simple/random_init.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,124 @@
! random init test
!
! Copyright (c) 2021-2021, Sourcery, Inc.
! All rights reserved.
!
! Redistribution and use in source and binary forms, with or without
! modification, are permitted provided that the following conditions are met:
! * Redistributions of source code must retain the above copyright
! notice, this list of conditions and the following disclaimer.
! * Redistributions in binary form must reproduce the above copyright
! notice, this list of conditions and the following disclaimer in the
! documentation and/or other materials provided with the distribution.
! * Neither the name of the Sourcery, Inc., nor the
! names of its contributors may be used to endorse or promote products
! derived from this software without specific prior written permission.
!
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
! ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
! WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
! DISCLAIMED. IN NO EVENT SHALL SOURCERY, INC., BE LIABLE FOR ANY
! DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
! (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
! ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
!

program test_random_init
implicit none
integer :: me,np
integer(kind=4), dimension(:), allocatable :: random_num, from_master
integer(kind=4) :: seed_size
integer :: seed_eq

me = this_image()
np = num_images()

if (np .lt. 1) then
error stop "Need at least two images."
end if

call random_seed(size=seed_size)
allocate(random_num(1:seed_size))
allocate(from_master(1:seed_size))

call random_init(.true., .true.)

sync all
call random_seed(get=random_num)
if (me .eq. 1) then
from_master = random_num
end if
call co_broadcast(from_master, 1)
if (me .eq. 1) then
seed_eq = 0
else
seed_eq = any(random_num .eq. from_master)
end if
call co_max(seed_eq, 1)

if (me .eq. 1 .and. seed_eq .eq. 1) then
error stop "Test failed. (T,T)"
end if

call random_init(.false., .true.)

sync all
call random_seed(get=random_num)
if (me .eq. 1) then
from_master = random_num
end if
call co_broadcast(from_master, 1)
if (me .eq. 1) then
seed_eq = 0
else
seed_eq = any(random_num .eq. from_master)
end if
call co_max(seed_eq, 1)

if (me .eq. 1 .and. seed_eq .eq. 1) then
error stop "Test failed. (F,T)"
end if

sync all

call random_init(.false., .false.)

sync all
call random_seed(get=random_num)
if (me .eq. 1) then
from_master = random_num
end if
call co_broadcast(from_master, 1)
seed_eq = all(random_num .eq. from_master)
call co_min(seed_eq, 1)

print *,"me=", me, ", rand_num=", random_num, ", from_master=", from_master, ", seed_eq=", seed_eq
if (me .eq. 1 .and. seed_eq .eq. 0) then
error stop "Test failed. (F,F)"
end if

sync all

call random_init(.true., .false.)

sync all
call random_seed(get=random_num)
if (me .eq. 1) then
from_master = random_num
end if
call co_broadcast(from_master, 1)
seed_eq = all(random_num .eq. from_master)
call co_min(seed_eq, 1)

if (me .eq. 1 .and. seed_eq .eq. 0) then
error stop "Test failed. (T,F)"
end if

sync all

if (me .eq. 1) print *,"Test passed."

end program test_random_init