TABLE OF CONTENTS


tests/future_ca_omp2 [ Unit tests ]

[ Top ] [ Unit tests ]

NAME

future_ca_omp2

SYNOPSIS

!$Id: future_ca_omp2.f90 552 2018-08-31 11:30:39Z mexas $

program future_ca_omp2

PURPOSE

The future* tests are not part of casup. These are to test emerging capabilities. This test checks coarrays inside an OpenMP parallel region with 2D HX.

DESCRIPTION

Run on 4 images only! This is just for demo purposes. This is a simple relaxation problem. A 2D real array coarray is set real( this_image() ) on all images. The halos are fixed boundary conditions. The 4-element kernel computes the average of the 4 neighbours, up, down, left, right. Run multiple iterations until convergence. The HX is implemented using sync images inside OMP.

AUTHOR

Anton Shterenlikht

COPYRIGHT

See LICENSE

USES

USED BY

SOURCE

use :: omp_lib
implicit none
integer, parameter :: n=10000, niter=100
real :: a(0:n+1,0:n+1)[2,*], b(0:n+1,0:n+1), time1, time2
integer :: grid(2), i, j, img, iter, tmp, tid, nthr, nimgs

  img = this_image()
 grid = this_image( a )
nimgs = num_images()

if ( nimgs .ne. 4 ) then
  write (*,*) "ERROR: this demo program runs only on 4 images"
  error stop
end if

! Set to image number
a = real( img )
b = a

sync all
if (img .eq. 1 ) call dump( n, a, "start.raw" )
sync all

  ! Start counter
  if ( img .eq. 1 ) call cpu_time( time1 )

main: do iter = 1, niter
  a = b
  sync all
  !$omp parallel do default(none) private(i,tmp,tid) &
  !$omp shared(img,a,b,nthr,grid)
  loop_j: do j = 1, n
  loop_i: do i = 1, n
    nthr = omp_get_num_threads()
    if ( i .eq. n .and. grid(1) .ne. 2 ) a(n+1,j) = a(1,j) [ grid(1)+1, grid(2) ]
    if ( i .eq. 1 .and. grid(1) .ne. 1 ) a(0,  j) = a(n,j) [ grid(1)-1, grid(2) ]
    if ( j .eq. n .and. grid(2) .ne. 2 ) a(i,n+1) = a(i,1) [ grid(1), grid(2)+1 ]
    if ( j .eq. 1 .and. grid(2) .ne. 1 ) a(i,  0) = a(i,n) [ grid(1), grid(2)-1 ]
    b(i,j) = 0.25 * ( a(i-1,j) + a(i+1,j) + a(i,j-1) + a(i,j+1) )
  end do loop_i
  end do loop_j
  !$omp end parallel do
  if (img .eq. 1 ) write (*,*) "iter:", iter
!  write (*,"(a,i0,tr1,i0,tr1,999i1)") "iter, img, b: ", iter, img, b(1:n)
end do main

sync all

  ! Start counter
  if ( img .eq. 1 ) then
    call cpu_time( time2 )
    write (*,*) "Time, s:", time2-time1
  end if

a=b
if (img .eq. 1 ) call dump( n, a, "end.raw" )

contains

! Call this sub only from image 1!
subroutine dump( n, arr, fname )
integer, intent(in) :: n
real, intent(in) :: arr(0:n+1,0:n+1)[2,*]
character(len=*), intent(in) :: fname
integer :: fu, coi1, coi2, j
open( newunit=fu, file=fname, status="replace", form="unformatted", &
      access="stream")

! Assume 2x2 images grid
! nested loops for writing in correct order from all images
! Do not write halos
 do coi2 = 1,2
  do j = 1, n
   do coi1 = 1,2
!    write (*,'(es10.2)') arr(1:n, j) [ coi1, coi2 ]
    write (fu) arr(1:n, j) [ coi1, coi2 ]
   end do
  end do
 end do

close( fu )
end subroutine dump

end program future_ca_omp2