TABLE OF CONTENTS
tests/ising_perf_1D [ Unit tests ]
[ Top ] [ Unit tests ]
NAME
ising_perf_1D
SYNOPSIS
!$Id: ising_perf_1D.f90 561 2018-10-14 20:48:19Z mexas $ program ising_perf_1D
PURPOSE
Test performance of ising magnetisation with 1D (corank 1). This test uses coarray collectives. Reproducibility on any number of images requires a serial RND for the whole model. This is too slow. Here each image uses its own RND, meaning the results are *not* reproducible when the number of images change. This test is purely for performance analysis.
DESCRIPTION
See ca_kernel_ising and related routines for details.
AUTHOR
Anton Shterenlikht
COPYRIGHT
See LICENSE
USES
cgca
USED BY
Part of CGPACK test suite
SOURCE
use casup implicit none integer( kind=iarr ), parameter :: huge_iarr = huge(0_iarr) real( kind=rdef ) :: bsz(3) ! updated "box" size integer( kind=iarr ), allocatable :: space(:,:,:), & space0(:,:,:) integer( kind=idef ) :: ir(3), nimgs, img, ng, c(3) ! coarray dimensions integer( kind=ilrg ) :: icells, mcells !real, allocatable :: space_ini(:,:,:), rnd_array(:) integer :: i, iter, seed_size, run integer( kind=ilrg ) :: energy0, energy1, energy2, magnet0, magnet1, & magnet2 integer, allocatable :: seed_array(:) ! logical :: flag real :: time1, time2 real, allocatable :: rnd_arr(:,:,:) !*********************************************************************72 ! first executable statement ! Read the box size from command line call ca_cmd_real( n=3, data=bsz ) ! When dm=res=1, then bsz is simply CA dimensions in cells! !bsz = (/ 2.0e3, 2.0e3, 2.0e3 /) ! dimensions of the CA model !bsz = (/ 3.0e3, 3.0e3, 3.0e3 /) ! dimensions of the CA model !bsz = (/ 1.2e2, 1.2e2, 1.2e2 /) ! 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 ! total number of cells in a coarray icells = product( int( c, kind=ilrg ) ) ! total number of cells in the model mcells = icells * int( nimgs, kind=ilrg ) 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" ! In this test sum over all cells on an image is done, so the kind ! must be big enough to contain the total number of cells ! on an image. if ( icells .gt. huge_iarr ) then write (*,*) "ERROR: number of cells on an image:", icells, & "is greater than huge(0_iarr)" error stop end if end if ! allocate space arrays and set all values to zero ! space - CA array to allocate, with halos! ! c - array with space dimensions ! d - depth of the halo layer call ca_spalloc( space, c, 1 ) call ca_spalloc( space0, c, 1 ) ! allocate hx arrays, implicit sync all ! mask_array is set inside too. ! call ca_1D_halloc ! Init RND !call cgca_irs( debug = .false. ) ! Use a reproducible RND here for verification call random_seed( size = seed_size ) allocate( seed_array( seed_size ) ) seed_array = (/ (i, i=1,seed_size) /) call random_seed( put = seed_array ) ! Set space arrays allocate( rnd_arr( lbound(space, dim=1) : ubound(space, dim=1), & lbound(space, dim=2) : ubound(space, dim=2), & lbound(space, dim=3) : ubound(space, dim=3) ) ) call random_number( rnd_arr ) space = nint( rnd_arr, kind=iarr ) ! MPI ! init ! in ! test_mpi_ising_perf ! ! ! MPI types ! creation ! Calculate initial energy and magnetisation ! run=1 => ca_iter_tl ! run=2 => ca_iter_dc ! run=3 => ca_iter_omp do run=1,3 select case(run) case(1) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_tl, kernel = ca_kernel_ising_ener, & energy = energy0, magnet = magnet0 ) case(2) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_dc, kernel = ca_kernel_ising_ener, & energy = energy1, magnet = magnet1 ) case(3) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_omp, kernel = ca_kernel_ising_ener, & energy = energy2, magnet = magnet2 ) end select end do if (img .eq. 1 ) then write (*,*) "Initial energy and magnetisation" write (*,*) "ca_iter_tl :", energy0, magnet0 write (*,*) "ca_iter_dc :", energy1, magnet1 write (*,*) "ca_iter_omp:", energy2, magnet2 if ( energy0 .ne. energy1 .or. magnet0 .ne. magnet1 .or. & energy0 .ne. energy2 .or. magnet0 .ne. magnet2 ) then write (*,*) "FAIL: ca_iter_tl, ca_iter_dc, ca_iter_omp differ" error stop else write (*,*) "PASS: ca_iter_tl, ca_iter_dc, ca_iter_omp agree" end if end if ! save old space as space0 space0 = space ! run=1 => ca_iter_tl ! run=2 => ca_iter_dc ! run=3 => ca_iter_omp main: do run=1,1 ! Reset space to space0 space = space0 ! Start counter if ( img .eq. 1 ) call cpu_time( time1 ) ! CA iterations loop: do iter = 1,100 ! Check energy after every iter ! 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_ising, niter = 1 ) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_tl, kernel = ca_kernel_ising_ener, & energy = energy1, magnet = magnet1 ) case(2) call ca_run( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_dc, kernel = ca_kernel_ising, niter = 1 ) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_dc, kernel = ca_kernel_ising_ener, & energy = energy1, magnet = magnet1 ) case(3) call ca_run( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_omp, kernel = ca_kernel_ising, niter = 1 ) call ca_ising_energy_col( space = space, hx_sub = ca_1D_hx_sall, & iter_sub = ca_iter_omp, kernel = ca_kernel_ising_ener, & energy = energy1, magnet = magnet1 ) end select if ( img .eq. 1 ) then if ( energy1 .ne. energy0 ) then write (*,*) "FAIL: energy0:", energy0, "energy1:", energy1 error stop else ! write (*,"(a,i0,a,es18.6)") "Magnetisation_after_iter_", & ! iter, ":", real(magnet1) / real(mcells) end if end if end do loop if ( img .eq. 1 ) then ! Stop counter call cpu_time( time2 ) select case(run) case(1) write (*,*) "TIME ca_iter_tl, s:", time2-time1 case(2) write (*,*) "TIME ca_iter_dc, s:", time2-time1 case(3) write (*,*) "TIME ca_iter_omp,s:", time2-time1 end select end if end do main ! deallocate halos, implicit sync all call ca_1D_hdalloc ! No need to free ! MPI types ! deallocate space deallocate( space ) deallocate( space0 ) end program ising_perf_1D