TABLE OF CONTENTS


ca_hx/ca_hx_1D [ Submodules ]

[ Top ] [ ca_hx ] [ Submodules ]

NAME

ca_hx_1D

SYNOPSIS

!$Id: ca_hx_1D.f90 558 2018-10-14 16:28:29Z mexas $

submodule ( ca_hx ) ca_hx_1D

DESCRIPTION

Submodule with routines for a single codimension CA coarrays. Whole array coarrays, not just halos.

AUTHOR

Anton Shterenlikht

COPYRIGHT

See LICENSE

CONTAINS

USES

USED BY

SOURCE

implicit none

! These are halo coarrays. Parts of the CA array are designated
! "halo", but those are not coarrays.
!
! "minu" - minus, "plus" - plus. 
integer( kind=iarr ), allocatable ::          &
  h1Dminu(:,:,:)[:], h1Dplus(:,:,:)[:]

integer :: nei1D_L, nei1D_R

contains

ca_hx_1D/ca_1D_halloc [ Subroutines ]

[ Top ] [ ca_hx_1D ] [ Subroutines ]

NAME

ca_1D_halloc

SYNOPSIS

module subroutine ca_1D_halloc

SIDE EFFECTS

halo coarrays (submodule variables) are allocated

DESCRIPTION

This routine allocates 2 halo coarrays for von Neumann 6-neighbourhood. The coarrays are submodule variables. Coarray allocation is an implicit sync. Halos have depth hdepth. Halo coarrays have the same cobounds on all images.

NOTES

All images must call this routine!

USES

USED BY

SOURCE

integer :: i,j,k

main: associate( do => hdepth, c => sub )

if ( this_image() .eq. 1 ) write (*,*) "1D_halloc: d:", do, "c:", c

allocate( h1Dminu( c(1), c(2), do ) [*], stat=ierr, errmsg=errmsg )
if ( ierr .ne. 0 ) then
  write (*,*) "ERROR: ca_hx_1D/ca_1D_halloc: allocate( h1Dminu )." //  &
              " ierr: ", ierr, "errmsg:", trim(errmsg)
  error stop
end if

allocate( h1Dplus( c(1), c(2), do ) [*], stat=ierr, errmsg=errmsg )
if ( ierr .ne. 0 ) then
  write (*,*) "ERROR: ca_hx_1D/ca_1D_halloc: allocate( h1Dplus )." //  &
              " ierr: ", ierr, "errmsg:", trim(errmsg)
  error stop
end if

! Calculate once and keep forever
! In 1D version space has a single codimension.
! To use the old code, we make the first 2 codimensions, and the
! first 2 cobounds are equal to 1.
    ci    = 1
  ucob    = 1
    ci(3) = this_image()
  ucob(3) = num_images()

! Now can set mask_array. The mask array must reflect
! the global CA space, i.e. must not be affected by partitioning
! of the model into images. This means the mask array must
! depend on coindex set of this image.
do concurrent( i=1:c(1), j=1:c(2), k=1:c(3) )
  mask_array(i,j,k) = int( mod( (i+j+k + ( ci(1)-1)*c(1) +             &
    (ci(2)-1)*c(2) + (ci(3)-1)*c(3) ) , 2 ), kind=iarr ) 
end do

end associate main

! Calculate image numbers for the 2 neighbours
  nei1D_L = this_image()-1
  nei1D_R = this_image()+1

end subroutine ca_1D_halloc

ca_hx_1D/ca_1D_hdalloc [ Subroutines ]

[ Top ] [ ca_hx_1D ] [ Subroutines ]

NAME

ca_1D_hdalloc

SYNOPSIS

module subroutine ca_1D_hdalloc

INPUT

!    none

OUTPUT

!    none

SIDE EFFECTS

halo coarrays (module variables) are deallocated

DESCRIPTION

This routine deallocates 2 halo coarrays for von Neumann 6-neighbourhood. The coarrays are module variables. Coarray deallocation is an implicit sync. Halo coarrays have the same cobounds on all images.

NOTES

All images must call this routine!

USES

USED BY

SOURCE

deallocate( h1Dminu, stat=ierr, errmsg=errmsg )
if ( ierr .ne. 0 ) then
  write (*,*) "ERROR: ca_hx_1D/ca_1D_hdalloc: " //                   &
    "deallocate( h1Dminu ). ierr: ", ierr, "errmsg:", trim(errmsg)
  error stop
end if

deallocate( h1Dplus, stat=ierr, errmsg=errmsg )
if ( ierr .ne. 0 ) then
  write (*,*) "ERROR: ca_hx_1D/ca_1D_hdalloc: " //                   &
    "deallocate( h1Dplus ). ierr: ", ierr, "errmsg:", trim(errmsg)
  error stop
end if

end subroutine ca_1D_hdalloc

ca_hx_1D/ca_1D_hx_check [ Subroutines ]

[ Top ] [ ca_hx_1D ] [ Subroutines ]

NAME

ca_1D_hx_check

SYNOPSIS

module subroutine ca_1D_hx_check( space, flag )

INPUT

integer( kind=iarr ), intent(in), allocatable :: space(:,:,:)

!    space - CA array

OUTPUT

integer, intent( out ) :: flag

!    flag .eq. 0 - check passed
!    flag .ne. 0 - check failed

SIDE EFFECTS

none

DESCRIPTION

This routine is of very limited use. It checks hx code for only a single case - all images set all their cell values to this_image() and then a single hx step is done. In this case I'm sure what's in my halos. So just check for this and either return flag=0 or flag > 0 indicating which halo failed. Fail values:

      0  - pass, no failures
      1  - test 1 failed
      2  - test 2 failed
      3  - tests 1,2 failed

USES

USED BY

SOURCE

flag = 0

! Test 1
if ( ci(3) .ne. 1       ) then
  if ( any( space( 1:sub(1), 1:sub(2), lhsta(3):0 )                    &
                    .ne. nei1D_L ) ) flag = 1
end if

! Test 2
if ( ci(3) .ne. ucob(3) ) then
  if ( any( space( 1:sub(1), 1:sub(2), rhsta(3):rhend(3) )             &
                    .ne. nei1D_R ) ) flag = 2
end if

end subroutine ca_1D_hx_check

ca_hx_1D/ca_1D_hx_sall [ Subroutines ]

[ Top ] [ ca_hx_1D ] [ Subroutines ]

NAME

ca_1D_hx_sall

SYNOPSIS

module subroutine ca_1D_hx_sall( space )

INPUT

integer( kind=iarr ), intent(inout), allocatable :: space(:,:,:)

!    space - non coarray array with CA model

OUTPUT

!    none

SIDE EFFECTS

halo coarrays are updated

DESCRIPTION

HX is a 2-step process. In step 1 I copy parts of my space array, not space halos!, into my halo coarrays. Step 2 is remote comms - I pull neighbour halo coarrays into my local space array halo sections. sync all is used for synchronisation (sall in the name). Note

USES

USED BY

SOURCE

if ( ci(3) .ne. 1       ) then
  h1Dminu(:,:,:) = space( 1:sub(1) , 1:sub(2) ,        1 : hdepth )
end if
if ( ci(3) .ne. ucob(3) ) then
  h1Dplus(:,:,:) = space( 1:sub(1) , 1:sub(2) , ihsta(3) : sub(3) )
end if

sync all

if ( ci(3) .ne. 1 ) then
  space( 1:sub(1) , 1:sub(2) , lhsta(3) : 0 ) =                        &
    h1Dplus(:,:,:) [ nei1D_L ]
end if

if ( ci(3) .ne. ucob(3) ) then
  space( 1:sub(1) , 1:sub(2) , rhsta(3) : rhend(3) ) =                 &
    h1Dminu(:,:,:) [ nei1D_R ]
end if

end subroutine ca_1D_hx_sall