TABLE OF CONTENTS
tests/hxvn_timing_mpi [ Unit tests ]
[ Top ] [ Unit tests ]
NAME
hxvn_timing_mpi
SYNOPSIS
!$Id: hxvn_timing_mpi.f90 548 2018-04-27 14:49:39Z mexas $ program hxvn_timing_mpi
PURPOSE
Time HX for HCA, WCA and MPI. No computation or any other work is done. Several HX sizes are used.
DESCRIPTION
- ca_halloc - user allocates halo coarrays once, obviously at the start of the simulation.
- ca_hx_all - a high level routine to do all necessary hx operations, with necessary sync.
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
USED BY
Part of casup test suite
SOURCE
use casup implicit none integer( kind=iarr ), parameter :: huge_iarr = huge(0_iarr) real( kind=rdef ) :: & qual, & ! quality bsz0(3), & ! the given "box" size bsz(3), & ! updated "box" size dm, & ! mean grain size, linear dim, phys units lres, & ! linear resolution, cells per unit of length res ! resolutions, cells per grain integer( kind=iarr ), allocatable :: space(:,:,:) integer( kind=idef ) :: ir(3), nimgs, img, ng, c(3), hsize, i, iter integer( kind=ilrg ) :: icells, mcells real, allocatable :: rnd(:,:,:) real( kind=kind(1.0d0) ) :: time1, time2 logical :: flag integer :: ierr !*********************************************************************72 ! first executable statement img = this_image() nimgs = num_images() ! do a check on image 1 if ( img .eq. 1 ) then write (*,*) "running on", nimgs, "images in a 3D grid" write (*,*) "iarr kind:", iarr, "huge(0_iarr):", huge_iarr ! In this test space is assigned only 0 or 1, and no collectives, ! so even 1 byte integers will do. end if ! each image calculates the coarray grid dimensions call cgca_gdim( nimgs, ir, qual ) dm = 1.0 ! cell size res = 1.0 ! resolution ! Loop over several box sizes, i.e. halo sizes main: do i=1,35 ! 50, 100, 150 bsz0 = (/ 100*i, 100*i, 100*i /) ! numbers of cells in CA space ! calculate the resolution and the actual phys dimensions of the box ! subroutine cgca_cadim( bsz, res, dm, ir, c, lres, ng ) ! c - coarray sizes ! ir - coarray grid sizes bsz = bsz0 call cgca_cadim( bsz, res, dm, ir, c, lres, ng ) ! Check that the partition is sane if ( img .eq. 1 ) then if ( any(int(bsz) .ne. int(bsz0) ) ) then write (*,*) & "ERROR: bad decomposition - use a 'nicer' number of images" write (*,*) "ERROR: wanted :", int(bsz0) write (*,*) "ERROR: but got instead:", int(bsz) error stop 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 ) ! Halo size hsize = c(1)*c(2) write ( *, "(8(a,i0),tr1,g10.3,tr1,g10.3,3(a,g10.3),a)" ) & "nimgs: ", nimgs, " (", c(1), "," , c(2), "," , c(3), ")[", & ir(1), "," , ir(2), "," , ir(3), "] ", ng, qual, lres, & " (", bsz(1), ",", bsz(2), ",", bsz(3), ")" write (*,'(a,i0,a)') "Each image has ",icells, " cells" write (*,'(a,i0,a)') "The model has ", mcells, " cells" write (*,'(a,i0,a)') "Halo size is ", hsize, " cells" end if ! Initialise MPI if not done already call MPI_INITIALIZED( flag, ierr) if ( .not. flag ) then call MPI_INIT( ierr ) if ( img .eq. 1 ) write (*,*) "MPI not initialised, doing now!" end if ! 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, 1 ) ! Very crude RND allocate( rnd( size(space, dim=1), size(space, dim=2), & size(space, dim=3) ) ) call random_number( rnd ) space = nint( rnd, kind=iarr ) ! Just 0 or 1 ! allocate hx array coarrays, implicit sync all ! ir(3) - codimensions call ca_halloc( ir ) ! Create MPI subarray types call ca_mpi_halo_type_create( space ) ! do several HX, remote ops, sync inside ! Make sure the timer covers the slowest image sync all time1 = cgca_benchtime() do iter=1,10 call ca_mpi_hx_all( space ) end do ! Make sure the timer covers the slowest image sync all time2 = cgca_benchtime() if (img .eq. 1) write (*,*) "Halo, cells:", hsize, ". Time, s:", & time2-time1 ! deallocate arrays call ca_hdalloc deallocate( space ) deallocate( rnd ) ! free halo types call ca_mpi_halo_type_free end do main if (img .eq. 1 ) write (*,*) "PASS" end program hxvn_timing_mpi