Skip to content
Open
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
4 changes: 2 additions & 2 deletions docs/implementation-status.md
Original file line number Diff line number Diff line change
Expand Up @@ -62,12 +62,12 @@ all of the PRIF-specific constants.

| Procedure | Status | Notes |
|-----------|--------|-------|
| `prif_allocate_coarray` | **YES** | |
| `prif_allocate_coarray` | **YES** | includes ucobound relaxation expected in PRIF 0.8 |
| `prif_allocate` | **YES** | |
| `prif_deallocate_coarray` | **YES** | `final_func` support requires flang 20+ |
| `prif_deallocate_coarrays` | **YES** | `final_func` support requires flang 20+ |
| `prif_deallocate` | **YES** | |
| `prif_alias_create` | **YES** | |
| `prif_alias_create` | **YES** | includes ucobound relaxation expected in PRIF 0.8 |
| `prif_alias_destroy` | **YES** | |

---
Expand Down
28 changes: 23 additions & 5 deletions src/caffeine/alias_s.F90
Original file line number Diff line number Diff line change
Expand Up @@ -12,10 +12,20 @@
contains

module procedure prif_alias_create
integer(c_int) :: corank

! validate inputs
call_assert(coarray_handle_check(source_handle))
corank = size(alias_lcobounds)
call_assert(corank > 0)
if (size(alias_ucobounds) == corank) then
call_assert(all(alias_lcobounds <= alias_ucobounds))
call_assert(product(alias_ucobounds - alias_lcobounds + 1) >= current_team%info%num_images)
else
call_assert(size(alias_ucobounds) == corank - 1)
call_assert(all(alias_lcobounds(1:corank-1) <= alias_ucobounds))
end if

call_assert(size(alias_lcobounds) == size(alias_ucobounds))
call_assert(product(alias_ucobounds - alias_lcobounds + 1) >= current_team%info%num_images)

allocate(alias_handle%info)
! start with a copy of the source descriptor
Expand All @@ -27,9 +37,17 @@
# endif

! apply provided cobounds
alias_handle%info%corank = size(alias_lcobounds)
alias_handle%info%lcobounds(1:size(alias_lcobounds)) = alias_lcobounds
alias_handle%info%ucobounds(1:size(alias_ucobounds)) = alias_ucobounds
alias_handle%info%corank = corank
alias_handle%info%lcobounds(1:corank) = alias_lcobounds
alias_handle%info%ucobounds(1:corank-1) = alias_ucobounds(1:corank-1)
call compute_coshape_epp(alias_lcobounds, alias_ucobounds, &
alias_handle%info%coshape_epp(1:corank))
# if ASSERTIONS
! The following entries are dead, but initialize them to help detect defects
alias_handle%info%lcobounds(corank+1:15) = huge(0_c_int64_t)
alias_handle%info%ucobounds(corank:14) = -huge(0_c_int64_t)
alias_handle%info%coshape_epp(corank+1:15) = 0
# endif

! reset some fields that are unused in aliases
alias_handle%info%reserved = c_null_ptr
Expand Down
25 changes: 20 additions & 5 deletions src/caffeine/allocation_s.F90
Original file line number Diff line number Diff line change
Expand Up @@ -19,11 +19,19 @@
type(c_ptr) :: whole_block
integer(c_ptrdiff_t) :: block_offset
integer(c_size_t) :: descriptor_size, total_size
integer(c_int) :: corank
type(prif_coarray_descriptor) :: unused
type(prif_coarray_descriptor), pointer :: unused2(:)

call_assert(size(lcobounds) == size(ucobounds))
call_assert(product(ucobounds - lcobounds + 1) >= current_team%info%num_images)
corank = size(lcobounds)
call_assert(corank > 0)
if (size(ucobounds) == corank) then
call_assert(all(lcobounds <= ucobounds))
call_assert(product(ucobounds - lcobounds + 1) >= current_team%info%num_images)
else
call_assert(size(ucobounds) == corank - 1)
call_assert(all(lcobounds(1:corank-1) <= ucobounds))
end if

me = current_team%info%this_image
if (caf_have_child_teams()) then
Expand Down Expand Up @@ -62,11 +70,18 @@
call c_f_pointer(whole_block, unused2, [2])

coarray_handle%info%coarray_data = c_loc(unused2(2))
coarray_handle%info%corank = size(lcobounds)
coarray_handle%info%corank = corank
coarray_handle%info%coarray_size = size_in_bytes
coarray_handle%info%final_func = final_func
coarray_handle%info%lcobounds(1:size(lcobounds)) = lcobounds
coarray_handle%info%ucobounds(1:size(ucobounds)) = ucobounds
coarray_handle%info%lcobounds(1:corank) = lcobounds
coarray_handle%info%ucobounds(1:corank-1) = ucobounds(1:corank-1)
call compute_coshape_epp(lcobounds, ucobounds, coarray_handle%info%coshape_epp(1:corank))
# if ASSERTIONS
! The following entries are dead, but initialize them to help detect defects
coarray_handle%info%lcobounds(corank+1:15) = huge(0_c_int64_t)
coarray_handle%info%ucobounds(corank:14) = -huge(0_c_int64_t)
coarray_handle%info%coshape_epp(corank+1:15) = 0
# endif
coarray_handle%info%previous_handle = c_null_ptr
coarray_handle%info%next_handle = c_null_ptr
call add_to_team_list(coarray_handle)
Expand Down
83 changes: 53 additions & 30 deletions src/caffeine/coarray_queries_s.F90
Original file line number Diff line number Diff line change
Expand Up @@ -19,29 +19,63 @@

module procedure prif_lcobound_no_dim
call_assert(coarray_handle_check(coarray_handle))
call_assert(size(lcobounds) == coarray_handle%info%corank)

lcobounds = coarray_handle%info%lcobounds(1:coarray_handle%info%corank)
end procedure

module procedure prif_ucobound_with_dim
call_assert(coarray_handle_check(coarray_handle))
call_assert(dim >= 1 .and. dim <= coarray_handle%info%corank)

ucobound = coarray_handle%info%ucobounds(dim)
associate (info => coarray_handle%info, corank => coarray_handle%info%corank)
call_assert(dim >= 1 .and. dim <= corank)

if (corank == 1) then ! common-case optimization
ucobound = info%lcobounds(1) + current_team%info%num_images - 1
elseif (dim < corank) then
ucobound = info%ucobounds(dim)
else ! compute trailing ucobound, based on current team size
call_assert(dim == corank)
associate (epp => info%coshape_epp(corank), num_imgs => current_team%info%num_images)
if (epp >= num_imgs) then ! optimization to skip a divide
ucobound = info%lcobounds(corank)
else
ucobound = info%lcobounds(corank) + (num_imgs + epp - 1) / epp - 1
end if
end associate
end if
end associate
end procedure

module procedure prif_ucobound_no_dim
call_assert(coarray_handle_check(coarray_handle))
call_assert(size(ucobounds) == coarray_handle%info%corank)

ucobounds = coarray_handle%info%ucobounds(1:coarray_handle%info%corank)
associate (corank => coarray_handle%info%corank)
ucobounds(1:corank-1) = coarray_handle%info%ucobounds(1:corank-1)
call prif_ucobound_with_dim(coarray_handle, corank, ucobounds(corank))
end associate
end procedure

module procedure prif_coshape
integer(c_int64_t) :: trailing_ucobound

call_assert(coarray_handle_check(coarray_handle))

associate(info => coarray_handle%info)
sizes = info%ucobounds(1:info%corank) - info%lcobounds(1:info%corank) + 1
call_assert(size(sizes) == coarray_handle%info%corank)

associate(info => coarray_handle%info, corank => coarray_handle%info%corank)
if (corank == 1) then ! common-case optimization
sizes(1) = current_team%info%num_images
else
sizes(1:corank-1) = info%ucobounds(1:corank-1) - info%lcobounds(1:corank-1) + 1
associate (epp => info%coshape_epp(corank), num_imgs => current_team%info%num_images)
if (epp >= num_imgs) then ! optimization to skip a divide
sizes(corank) = 1
else
sizes(corank) = (num_imgs + epp - 1) / epp
end if
end associate
end if
end associate
end procedure

Expand All @@ -53,29 +87,24 @@ subroutine image_index_helper(coarray_handle, sub, num_images, image_index)
integer(c_int), intent(out) :: image_index

integer :: dim
integer(c_int) :: prior_size

call_assert(coarray_handle_check(coarray_handle))

associate (info => coarray_handle%info)
call_assert(size(sub) == info%corank)
if (sub(1) .lt. info%lcobounds(1) .or. sub(1) .gt. info%ucobounds(1)) then
associate (info => coarray_handle%info, corank => coarray_handle%info%corank)
call_assert(size(sub) == corank)
if (sub(1) .lt. info%lcobounds(1) .or. &
(corank > 1 .and. sub(1) .gt. info%ucobounds(1))) then
image_index = 0
return
end if
image_index = 1 + INT(sub(1) - info%lcobounds(1), c_int)
prior_size = 1
! Future work: values of prior_size are invariant across calls w/ the same coarray_handle
! We could store them in the coarray metadata at allocation rather than redundantly
! computing them here, which would accelerate calls with corank > 1 by removing
! corank multiply/add operations and the loop-carried dependence
do dim = 2, size(sub)
prior_size = prior_size * INT(info%ucobounds(dim-1) - info%lcobounds(dim-1) + 1, c_int)
if (sub(dim) .lt. info%lcobounds(dim) .or. sub(dim) .gt. info%ucobounds(dim)) then
if (sub(dim) .lt. info%lcobounds(dim) .or. &
(dim < corank .and. sub(dim) .gt. info%ucobounds(dim))) then
image_index = 0
return
end if
image_index = image_index + INT(sub(dim) - info%lcobounds(dim), c_int) * prior_size
image_index = image_index + INT(sub(dim) - info%lcobounds(dim), c_int) * info%coshape_epp(dim)
end do
end associate

Expand Down Expand Up @@ -112,23 +141,17 @@ subroutine initial_index_helper(coarray_handle, sub, team, initial_team_index)
integer(c_int), intent(out) :: initial_team_index

integer :: dim
integer(c_int) :: prior_size, image_index
integer(c_int) :: image_index

call_assert(coarray_handle_check(coarray_handle))

associate (info => coarray_handle%info)
call_assert(size(sub) == info%corank)
call_assert(sub(1) .ge. info%lcobounds(1) .and. sub(1) .le. info%ucobounds(1))
associate (info => coarray_handle%info, corank => coarray_handle%info%corank)
call_assert(size(sub) == corank)
call_assert(sub(1) .ge. info%lcobounds(1) .and. (corank == 1 .or. sub(1) .le. info%ucobounds(1)))
image_index = 1 + INT(sub(1) - info%lcobounds(1), c_int)
prior_size = 1
! Future work: values of prior_size are invariant across calls w/ the same coarray_handle
! We could store them in the coarray metadata at allocation rather than redundantly
! computing them here, which would accelerate calls with corank > 1 by removing
! corank multiply/add operations and the loop-carried dependence
do dim = 2, size(sub)
prior_size = prior_size * INT(info%ucobounds(dim-1) - info%lcobounds(dim-1) + 1, c_int)
call_assert(sub(dim) .ge. info%lcobounds(dim) .and. sub(dim) .le. info%ucobounds(dim))
image_index = image_index + INT(sub(dim) - info%lcobounds(dim), c_int) * prior_size
call_assert(sub(dim) .ge. info%lcobounds(dim) .and. (dim == corank .or. sub(dim) .le. info%ucobounds(dim)))
image_index = image_index + INT(sub(dim) - info%lcobounds(dim), c_int) * info%coshape_epp(dim)
end do
end associate

Expand Down
2 changes: 1 addition & 1 deletion src/caffeine/image_queries_s.F90
Original file line number Diff line number Diff line change
Expand Up @@ -41,6 +41,7 @@
integer :: dim

call_assert(coarray_handle_check(coarray_handle))
call_assert(size(cosubscripts) == coarray_handle%info%corank)

if (present(team)) then
offset = team%info%this_image - 1
Expand All @@ -58,7 +59,6 @@
offset = offset / dsz
end do
cosubscripts(info%corank) = offset + info%lcobounds(info%corank)
call_assert(cosubscripts(info%corank) <= info%ucobounds(info%corank))
end associate

# if ASSERTIONS
Expand Down
35 changes: 30 additions & 5 deletions src/caffeine/prif_private_s.F90
Original file line number Diff line number Diff line change
Expand Up @@ -431,6 +431,28 @@ pure function optional_value(var) result(c_val)
end if
end function

pure subroutine compute_coshape_epp(lcobounds, ucobounds, coshape_epp)
!! Compute the exclusive prefix product of the coshape for the given cobounds
integer(c_int64_t), intent(in) :: lcobounds(:), ucobounds(:)
integer(c_int), intent(out) :: coshape_epp(:)
integer(c_int64_t) :: product
integer :: d

associate (corank => size(lcobounds))
call_assert(corank > 0)
call_assert(size(coshape_epp) == corank)
call_assert(size(ucobounds) == corank .or. size(ucobounds) == corank-1)

coshape_epp(1) = 1
product = 1
do d = 2, corank
product = product * (ucobounds(d-1) - lcobounds(d-1) + 1)
call_assert_describe(product < huge(0_c_int), "Overflow in cobounds. product(coshape(a)) must be < 2 billion")
coshape_epp(d) = int(product, c_int)
end do
end associate
end subroutine

! Report the provided error stat/msg using the provided optional stat/errmsg args
subroutine report_error(report_stat, report_msg, stat, errmsg, errmsg_alloc)
integer(c_int), intent(in) :: report_stat
Expand Down Expand Up @@ -460,16 +482,19 @@ elemental impure function coarray_handle_check(coarray_handle) result(result_)
implicit none
type(prif_coarray_handle), intent(in) :: coarray_handle
logical :: result_
integer(c_int) :: i
integer(c_int) :: i, epp(15)

call assert_always(associated(coarray_handle%info), "unassociated info pointer in prif_coarray_handle")
associate(info => coarray_handle%info)
call assert_always(info%corank >= 1, "invalid corank in prif_coarray_handle")
call assert_always(info%corank <= size(info%ucobounds), "invalid corank in prif_coarray_handle")
call assert_always(all([(info%lcobounds(i) <= info%ucobounds(i), i = 1, info%corank)]), &
associate(info => coarray_handle%info, corank => coarray_handle%info%corank)
call assert_always(corank >= 1, "invalid corank in prif_coarray_handle")
call assert_always(corank <= size(info%lcobounds), "invalid corank in prif_coarray_handle")
call assert_always(all([(info%lcobounds(i) <= info%ucobounds(i), i = 1, corank-1)]), &
"invalid cobounds in prif_coarray_handle")
call assert_always(info%coarray_size > 0, "invalid data size in prif_coarray_handle")
call assert_always(c_associated(info%coarray_data), "invalid data pointer in prif_coarray_handle")
call compute_coshape_epp(info%lcobounds(1:corank),info%ucobounds(1:corank-1),epp(1:corank))
call assert_always(all(info%coshape_epp(1:corank) == epp(1:corank)), &
"invalid coshape_epp in prif_coarray_handle")
end associate

result_ = .true.
Expand Down
3 changes: 2 additions & 1 deletion src/prif.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1179,7 +1179,8 @@ module subroutine prif_atomic_ref_logical_indirect(image_num, atom_remote_ptr, v
integer(c_size_t) :: coarray_size
type(c_funptr) :: final_func
type(c_ptr) :: previous_handle = c_null_ptr, next_handle = c_null_ptr
integer(c_int64_t) :: lcobounds(15), ucobounds(15)
integer(c_int64_t) :: lcobounds(15), ucobounds(14)
integer(c_int) :: coshape_epp(15)
type(c_ptr) :: p_context_data
type(c_ptr) :: reserved
end type
Expand Down
Loading