TABLE OF CONTENTS


tests/hxvn_1D [ Unit tests ]

[ Top ] [ Unit tests ]

NAME

hxvn_1D

SYNOPSIS

!$Id: hxvn_1D.f90 561 2018-10-14 20:48:19Z mexas $

program hxvn_1D

PURPOSE

Test ca_1D_hx_sall. The kernel is a simply copy. only halos are coarrays, the rest of the model is not.

DESCRIPTION

Must work on any number of images, except when a good decomposition cannot be made. The user needs know nothing about sync.

NOTE

AUTHOR

Anton Shterenlikht

COPYRIGHT

See LICENSE

USES

casup

USED BY

Part of casup test suite

SOURCE

use casup

implicit none

integer( kind=iarr ), parameter :: huge_iarr = huge(0_iarr)

real( kind=rdef ) :: bsz(3)            ! "box" size







integer( kind=iarr ), allocatable :: space(:,:,:),                     &
  space1(:,:,:)

integer( kind=idef ) :: nimgs, img, c(3) ! coarray dimensions

integer( kind=ilrg ) :: icells, mcells

integer :: ierr, do, run

! logical :: flag

!*********************************************************************72
! first executable statement

! Read the box size from command line
call ca_cmd_real( n=3, data=bsz )

!bsz = (/ 50.0, 50.0, 24000.0 /) ! numbers of cells in CA space
!bsz0 = (/ 4.0e1, 8.0e1, 6.0e1 /) ! for testing on FreeBSD laptop



  img = this_image()
nimgs = num_images()

! do a check on image 1
if ( img .eq. 1 ) then
  write (*,*) "running on", nimgs, "images in a 1D grid"
  write (*,*) "iarr kind:", iarr, "huge(0_iarr):", huge_iarr

  ! In this test space is assigned image numbers - must be big enough
  ! integer kind to avoid integer overflow.
  if ( nimgs .gt. huge_iarr ) then
    write (*,*) "ERROR: num_images(): ", nimgs,              &
      " is greater than huge(0_iarr)"
    error stop
  end if
end if

! Only a single codimension, corank 1
c = int(bsz)
c(3) = int(bsz(3))/nimgs

! Check that the partition is sane
if ( img .eq. 1 ) then
  if ( c(3)*nimgs .ne. int(bsz(3)) ) then
    write (*,*)                                                        &
      "ERROR: bad decomposition: bsz(3) must be divisible by nimgs"
    write (*,*) "ERROR: bsz(3):", int(bsz(3))
    write (*,*) "ERROR: nimgs :", nimgs
    error stop
  end if
end if

! total number of cells in a coarray
icells = int( c(1), kind=ilrg ) * int( c(2), kind=ilrg ) *             &
         int( c(3), kind=ilrg )

! total number of cells in the model
mcells = icells * int( nimgs, kind=ilrg )

if ( img .eq. 1 ) then
  write ( *, "(5(a,i0),3(a,g10.3),a)" )                                &
    "nimgs: ", nimgs, " (", c(1), "," , c(2), "," , c(3),              &
    ")[", nimgs, "] (", bsz(1), ",", bsz(2), ",", bsz(3), ")"
     write (*,'(a,i0,a)') "Each image has ",icells, " cells"
     write (*,'(a,i0,a)') "The model has ", mcells, " cells"
end if

! MPI
! init
! lines
! here
! in
! test_mpi_hxvn

! run=1 => ca_iter_tl
! run=2 => ca_iter_dc
! run=3 => ca_iter_omp
outer: do run=1,3

  if ( img .eq. 1 ) then
    select case( run )
    case(1)
      write (*,*) "Checking ca_iter_tl - triple loop"
    case(2)
      write (*,*) "Checking ca_iter_dc - do concurrent"
    case(3)
      write (*,*) "Checking ca_iter_omp - OpenMP"
    end select
  end if

  ! Loop over several halo depths
  ! The max halo depth is 1/4 of the min dimension
  ! of the space CA array
  main: do do=1, int( 0.25 * min( c(1), c(2), c(3) ) )
  
    ! allocate space array
    !    space - CA array to allocate, with halos!
    !        c - array with space dimensions
    !        d - depth of the halo layer

    call ca_spalloc( space,  c, do )
    call ca_spalloc( space1, c, do )
  
    ! Set space to my image number
    space  = int( img, kind=iarr )
    space1 = space

    ! allocate hx arrays, implicit sync all

    call ca_1D_halloc

    ! MPI subarray types in
    ! test_mpi_hxvn

    ! do hx, remote ops
    call ca_1D_hx_sall(     space )
  
    ! halo check, local ops
    ! space - space array, with halos
    ! flag - default integer
    call ca_1D_hx_check(    space=space, flag=ierr )
    if ( ierr .ne. 0 ) then
      write (*,*) "ERROR: ca_hx_check    failed: img:", img,           &
                  "flag:", ierr
      error stop
    end if

    ! CA iterations
    ! subroutine ca_run(    space, hx_sub, iter_sub, kernel, niter )
    select case( run )
    case(1)
      call ca_run(    space = space, hx_sub = ca_1D_hx_sall,           &
        iter_sub = ca_iter_tl,  kernel = ca_kernel_copy, niter = 13 )
    case(2)
      call ca_run(    space = space, hx_sub = ca_1D_hx_sall,           &
        iter_sub = ca_iter_dc,  kernel = ca_kernel_copy, niter = 13 )
    case(3)
      call ca_run(    space = space, hx_sub = ca_1D_hx_sall,           &
        iter_sub = ca_iter_omp, kernel = ca_kernel_copy, niter = 13 )
    end select
  
    ! Must be the same
    if ( any( space(  1:c(1), 1:c(2), 1:c(3) ) .ne.                    &
              space1( 1:c(1), 1:c(2), 1:c(3) ) ) ) then
      write (*,*) "img:", img, "FAIL: space .ne. space1"
      error stop
    end if

    ! deallocate halos, implicit sync all
    call ca_1D_hdalloc

    ! Free MPI types in
    ! test_mpi_hxvn

    ! deallocate space
    deallocate( space )
    deallocate( space1 )
  
    if (img .eq. 1 ) write (*,*) "PASS, halo depth:", do

  end do main

end do outer

end program hxvn_1D