Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
Show all changes
22 commits
Select commit Hold shift + click to select a range
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
17 changes: 11 additions & 6 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -802,23 +802,29 @@ if(opencoarrays_aware_compiler)
add_caf_test(strided_sendget 3 strided_sendget)
add_caf_test(get_with_vector_index 4 get_with_vector_index)


# Collective subroutine tests
add_caf_test(co_sum 4 co_sum_test)
add_caf_test(co_broadcast 4 co_broadcast_test)
add_caf_test(co_broadcast_derived_type 4 co_broadcast_derived_type_test)
if((gfortran_compiler AND (NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 10.0.0)) OR (CAF_RUN_DEVELOPER_TESTS OR $ENV{OPENCOARRAYS_DEVELOPER}))
add_caf_test(co_broadcast_allocatable_components 4 co_broadcast_allocatable_components_test)
endif()
add_caf_test(co_min 4 co_min_test)
add_caf_test(co_max 4 co_max_test)
add_caf_test(syncall 8 syncall)
add_caf_test(syncimages 8 syncimages)
add_caf_test(syncimages2 8 syncimages2)
add_caf_test(duplicate_syncimages 8 duplicate_syncimages)
add_caf_test(co_reduce 4 co_reduce_test)
add_caf_test(co_reduce_res_im 4 co_reduce_res_im)
add_caf_test(co_reduce_string 4 co_reduce_string)
add_caf_test(syncimages_status 8 syncimages_status)
add_caf_test(sync_ring_abort_np3 3 sync_image_ring_abort_on_stopped_image)
add_caf_test(sync_ring_abort_np7 7 sync_image_ring_abort_on_stopped_image)
add_caf_test(simpleatomics 8 atomics)

# Synchronization tests
add_caf_test(syncall 8 syncall)
add_caf_test(syncimages 8 syncimages)
add_caf_test(syncimages2 8 syncimages2)
add_caf_test(duplicate_syncimages 8 duplicate_syncimages)

# possible logic error in the following test
# add_caf_test(increment_my_neighbor 32 increment_my_neighbor)

Expand All @@ -828,7 +834,6 @@ if(opencoarrays_aware_compiler)
add_caf_test(co_heat 2 co_heat)
add_caf_test(asynchronous_hello_world 3 asynchronous_hello_world)


# Regression tests based on reported issues
if((gfortran_compiler AND (NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 7.0.0)) OR (CAF_RUN_DEVELOPER_TESTS OR $ENV{OPENCOARRAYS_DEVELOPER}))
if( CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 7.0.0 )
Expand Down
2 changes: 0 additions & 2 deletions src/mpi/mpi_caf.c
Original file line number Diff line number Diff line change
Expand Up @@ -7380,8 +7380,6 @@ PREFIX(co_broadcast) (gfc_descriptor_t *a, int source_image, int *stat,
size *= dimextent;
}

printf("DTYPE Size: %zd\n",GFC_DESCRIPTOR_SIZE(a));

if (rank == 0)
{
if( datatype == MPI_BYTE)
Expand Down
3 changes: 3 additions & 0 deletions src/tests/unit/collectives/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -1,6 +1,9 @@
caf_compile_executable(co_sum_test co_sum.F90)
caf_compile_executable(co_broadcast_test co_broadcast.F90)
caf_compile_executable(co_broadcast_derived_type_test co_broadcast_derived_type.f90)
if((gfortran_compiler AND (NOT CMAKE_Fortran_COMPILER_VERSION VERSION_LESS 10.0.0)) OR (CAF_RUN_DEVELOPER_TESTS OR $ENV{OPENCOARRAYS_DEVELOPER}))
caf_compile_executable(co_broadcast_allocatable_components_test co_broadcast_allocatable_components.f90)
endif()
caf_compile_executable(co_min_test co_min.F90)
caf_compile_executable(co_max_test co_max.F90)
caf_compile_executable(co_reduce_test co_reduce.F90)
Expand Down
62 changes: 62 additions & 0 deletions src/tests/unit/collectives/co_broadcast_allocatable_components.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,62 @@
program main
!! author: Damian Rouson
!!
!! Test co_broadcast with derived-type actual arguments
implicit none

integer, parameter :: sender=1 !! co_broadcast source_image
character(len=*), parameter :: text="text" !! character message data

type dynamic
character(len=:), allocatable :: string
character(len=len(text)), allocatable :: string_array(:)
complex, allocatable :: scalar
integer, allocatable :: vector(:)
logical, allocatable :: matrix(:,:)
real, allocatable :: superstring(:,:,:, :,:,:, :,:,:, :,:,:, :,:,: )
end type

type(dynamic) alloc_message, alloc_content

associate(me=>this_image())

alloc_content = dynamic( &
string=text, &
string_array=[text], &
scalar=(0.,1.), &
vector=reshape( [integer::], [0]), &
matrix=reshape( [.true.], [1,1]), &
superstring=reshape([1,2,3,4], [2,1,2, 1,1,1, 1,1,1, 1,1,1, 1,1,1 ]) &
)

alloc_message = alloc_content

if(me /= sender) then
alloc_message%vector = 0
alloc_message%matrix = .false.
alloc_message%superstring = 0
endif

sync all

call co_broadcast(alloc_message,source_image=sender)

associate( failures => [ &
alloc_message%string /= alloc_content%string, &
alloc_message%string_array /= alloc_content%string_array, &
alloc_message%scalar /= alloc_content%scalar, &
alloc_message%vector /= alloc_content%vector, &
alloc_message%matrix .neqv. alloc_content%matrix, &
alloc_message%superstring /= alloc_content%superstring &
] )

if ( any(failures) ) error stop "Test failed."

end associate

sync all ! Wait for each image to pass the test
if (me==sender) print *,"Test passed."

end associate

end program main
176 changes: 90 additions & 86 deletions src/tests/unit/collectives/co_broadcast_derived_type.f90
Original file line number Diff line number Diff line change
@@ -1,95 +1,99 @@
module object_interface
program main
!! author: Damian Rouson
!!
!! Test co_broadcast with derived-type actual arguments
implicit none
private
public :: object

type object
private
integer :: foo=0
logical :: bar=.false.
contains
procedure :: initialize
procedure :: co_broadcast_me
procedure :: not_equal
procedure :: copy
generic :: operator(/=)=>not_equal
generic :: assignment(=)=>copy
end type

integer, parameter :: sender=1 !! co_broadcast source_image
character(len=*), parameter :: text="text" !! character message data

interface
elemental impure module subroutine initialize(this,foo_,bar_)
implicit none
class(object), intent(out) :: this
integer, intent(in) :: foo_
logical, intent(in) :: bar_
end subroutine

elemental impure module subroutine co_broadcast_me(this,source_image)
implicit none
class(object), intent(inout) :: this
integer, intent(in) :: source_image
end subroutine

elemental module function not_equal(lhs,rhs) result(lhs_ne_rhs)
implicit none
class(object), intent(in) :: lhs,rhs
logical lhs_ne_rhs
end function

elemental impure module subroutine copy(lhs,rhs)
implicit none
class(object), intent(inout) :: lhs
class(object), intent(in) :: rhs
end subroutine
function f(x) result(y)
real x, y
end function
end interface

end module
type parent
integer :: heritable=0
end type

submodule(object_interface) object_implementation
implicit none
contains
module procedure co_broadcast_me
call co_broadcast(this%foo,source_image)
call co_broadcast(this%bar,source_image)
end procedure

module procedure initialize
this%foo = foo_
this%bar = bar_
end procedure

module procedure not_equal
lhs_ne_rhs = (lhs%foo /= rhs%foo) .or. (lhs%bar .neqv. rhs%bar)
end procedure

module procedure copy
lhs%foo = rhs%foo
lhs%bar = rhs%bar
end procedure
end submodule
type component
integer :: subcomponent=0
end type

program main
use object_interface, only : object
implicit none
type(object) message

call message%initialize(foo_=1,bar_=.true.)

emulate_co_broadcast: block
type(object) foobar
if (this_image()==1) foobar = message
call foobar%co_broadcast_me(source_image=1)
if ( foobar /= message ) error stop "Test failed."
end block emulate_co_broadcast

desired_co_broadcast: block
type(object) barfoo
if (this_image()==1) barfoo = message
call co_broadcast(barfoo,source_image=1) ! OpenCoarrays terminates here with the message "Unsupported data type"
if ( barfoo /= message ) error stop "Test failed."
end block desired_co_broadcast

sync all ! Wait for each image to pass the test
if (this_image()==1) print *,"Test passed."
type, extends(parent) :: child

! Scalar and array derived-type components
type(component) a, b(1,2,1, 1,1,1, 1)

! Scalar and array intrinsic-type components
character(len=len(text)) :: c="", z(0)
complex :: i=(0.,0.), j(1)=(0.,0.)
integer :: k=0, l(2,3)=0
logical :: r=.false., s(1,2,3, 1,2,3, 1)=.false.
real :: t=0., u(3,2,1)=0.

! Scalar and array pointer components
character(len=len(text)), pointer :: &
char_ptr=>null(), char_ptr_maxdim(:,:,:, :,:,:, :)=>null()
complex, pointer :: cplx_ptr=>null(), cplx_ptr_maxdim(:,:,:, :,:,:, :)=>null()
integer, pointer :: int_ptr =>null(), int_ptr_maxdim (:,:,:, :,:,:, :)=>null()
logical, pointer :: bool_ptr=>null(), bool_ptr_maxdim(:,:,:, :,:,:, :)=>null()
real, pointer :: real_ptr=>null(), real_ptr_maxdim(:,:,:, :,:,:, :)=>null()
procedure(f), pointer :: procedure_pointer=>null()
end type

type(child) message
type(child) :: content = child( & ! define content using the insrinsic structure constructor
parent=parent(heritable=-4), & ! parent
a=component(-3), b=reshape([component(-2),component(-1)], [1,2,1, 1,1,1, 1]), & ! derived types
c=text, z=[character(len=len(text))::], i=(0.,1.), j=(2.,3.), k=4, l=5, r=.true., s=.true., t=7., u=8. & ! intrinsic types
)

associate(me=>this_image())

if (me==sender) then
message = content
allocate(message%char_ptr, message%char_ptr_maxdim(1,1,2, 1,1,1, 1), source=text )
allocate(message%cplx_ptr, message%cplx_ptr_maxdim(1,1,1, 1,1,2, 1), source=(0.,1.))
allocate(message%int_ptr , message%int_ptr_maxdim (1,1,1, 1,1,1, 1), source=2 )
allocate(message%bool_ptr, message%bool_ptr_maxdim(1,1,1, 1,2,1, 1), source=.true. )
allocate(message%real_ptr, message%real_ptr_maxdim(1,1,1, 1,1,1, 1), source=3. )
end if

call co_broadcast(message,source_image=sender)

if (me==sender) then
deallocate(message%char_ptr, message%char_ptr_maxdim)
deallocate(message%cplx_ptr, message%cplx_ptr_maxdim)
deallocate(message%int_ptr , message%int_ptr_maxdim )
deallocate(message%bool_ptr, message%bool_ptr_maxdim)
deallocate(message%real_ptr, message%real_ptr_maxdim)
end if

!! Verify correct broadcast of all non-pointer components (pointers become undefined on the receiving image).
associate( failures => [ &
message%parent%heritable /= content%parent%heritable, &
message%a%subcomponent /= content%a%subcomponent, &
message%c /= content%c, &
message%z /= content%z, &
message%i /= content%i, &
message%j /= content%j, &
message%k /= content%k, &
message%l /= content%l, &
message%r .neqv. content%r, &
message%s .neqv. content%s, &
message%t /= content%t, &
any( message%u /= content%u ) &
] )

if ( any(failures) ) error stop "Test failed. "

end associate

sync all ! Wait for each image to pass the test
if (me==sender) print *,"Test passed."

end associate

end program main