Skip to content

Defect: same-image copies incorrect with teams #632

Closed
@nathanweeks

Description

@nathanweeks

Defect/Bug Report

  • OpenCoarrays Version: 2.4.0-6-gebc6f30
  • Fortran Compiler: gfortran 8.2.0
  • C compiler used for building lib: gcc 8.2.0
  • Installation method: CC=gcc FC=gfortran cmake .. -DCMAKE_BUILD_TYPE=Debug -DCMAKE_INSTALL_PREFIX=...
  • Output of uname -a: Linux ... 4.4.156-94.61.1.16335.0.PTF.1107299-default tests dis_transpose: test passed  #1 SMP ... x86_64 x86_64 x86_64 GNU/Linux
  • MPI library being used: MPICH 3.3
  • Machine architecture and number of physical cores:
  • Version of CMake: 3.11.4

Observed Behavior

*caf_get(), *caf_send(), and *caf_sendget() detect when the cosubscript(s) of the coarray(s) refer to the same image as the executing image; in that case, an optimized copy (avoiding MPI) is performed. However, this detection can be incorrect when the image indexes of the current team don't correspond to the cobounds of the coarray(s), leading to incorrect results.

Expected Behavior

Steps to Reproduce

The following examples form teams that have one image per team.

for *caf_get():

program teams_coarray_get
  use, intrinsic :: iso_fortran_env, only: team_type
  implicit none
  type(team_type) :: team
  integer, allocatable :: L(:)
  integer :: i, R[*]

  allocate(L(num_images()))
  R = this_image()

  form team (this_image(), team)
  change team (team)
    do i = 1, size(L)
      L(i) = R[i]
    end do
  end team

  write(*,*) this_image(), ':', L
end program teams_coarray_get

With 3 images, all images are expected to have L = [1,2,3]; however, only image 1 does:

$ caf teams_coarray_get.f90
$ cafrun -np 3 ./a.out
           1 :           1           2           3
           2 :           2           2           3
           3 :           3           2           3

*caf_send():

program teams_coarray_send
  use, intrinsic :: iso_fortran_env, only: team_type
  implicit none
  type(team_type) :: team
  integer, allocatable :: R(:)[:]
  integer :: i, initial_team_this_image

  allocate(R(num_images())[*], source=0)
  initial_team_this_image = this_image()

  form team (this_image(), team)
  change team (team)
    do i = 1, size(R)
      R(initial_team_this_image)[i] = initial_team_this_image
    end do
  end team

  write(*,*) this_image(), ':', R
end program teams_coarray_send

With 3 images, all are expected to have R = [1,2,3]; image 1's result is incorrect:

$ caf teams_coarray_send.f90
$ cafrun -n 3 ./a.out
           1 :           1           0           0
           2 :           1           2           3
           3 :           1           2           3

*caf_sendget():

program teams_coarray_sendget
  use, intrinsic :: iso_fortran_env, only: team_type
  implicit none
  type(team_type) :: team
  integer, allocatable :: R1(:,:)[:]
  integer :: i, j, R2[*]

  allocate(R1(num_images(), num_images())[*], source=0)
  R2 = this_image()

  form team (this_image(), team)
  change team (team)
    do concurrent (i = 1:ubound(R1,1), j = 1:ubound(R1,1))
      R1(R2,j)[i] = R2[j]
    end do
  end team

  write(*,*) this_image(), ':', R1
end program teams_coarray_sendget

On each image, R1 is expected to be a 3x3 matrix of the form:

1 2 3
1 2 3
1 2 3

("1 1 1 2 2 2 3 3 3" when output by the gfortran write(,) statement). However:

$ caf teams_coarray_sendget.f90
$ cafrun -np 3 ./a.out
           1 :           1           0           0           2           0           0           3           0           0
           2 :           1           2           3           2           2           2           3           3           3
           3 :           1           2           3           2           2           2           3           3           3

Metadata

Metadata

Assignees

No one assigned

    Labels

    No labels
    No labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions