! ! Copyright (C) by Argonne National Laboratory ! See COPYRIGHT in top-level directory ! program main implicit none include 'mpif.h' ! Fortran equivalent of perf.c integer SIZE parameter (SIZE=1048576*4) ! read/write size per node in bytes integer buf(SIZE/4), j, mynod, nprocs, ntimes, flag double precision stim, read_tim, write_tim, new_read_tim double precision new_write_tim, min_read_tim, min_write_tim double precision read_bw, write_bw integer fh, status(MPI_STATUS_SIZE), ierr, argc, iargc, i character*1024 str ! used to store the filename integer*8 offset ntimes = 5 min_read_tim = 10000000.0D0 min_write_tim = 10000000.0D0 call MPI_INIT(ierr) call MPI_COMM_SIZE(MPI_COMM_WORLD, nprocs, ierr) call MPI_COMM_RANK(MPI_COMM_WORLD, mynod, ierr) ! process 0 takes the file name as a command-line argument and ! broadcasts it to other processes if (mynod .eq. 0) then argc = iargc() i = 0 call getarg(i,str) do while ((i .lt. argc) .and. (str .ne. '-fname')) i = i + 1 call getarg(i,str) end do if (i .ge. argc) then print * print *, '*# Usage: fperf -fname filename' print * call MPI_ABORT(MPI_COMM_WORLD, 1, ierr) end if i = i + 1 call getarg(i,str) call MPI_BCAST(str, 1024, MPI_CHARACTER, 0, & & MPI_COMM_WORLD, ierr) print *, 'Access size per process = ', SIZE, ' bytes', & & ', ntimes = ', ntimes else call MPI_BCAST(str, 1024, MPI_CHARACTER, 0, & & MPI_COMM_WORLD, ierr) end if offset = mynod*SIZE do j=1, ntimes call MPI_FILE_OPEN(MPI_COMM_WORLD, str, & & MPI_MODE_CREATE+MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr) call MPI_FILE_SEEK(fh, offset, MPI_SEEK_SET, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) stim = MPI_WTIME() call MPI_FILE_WRITE(fh, buf, SIZE, MPI_BYTE, status, ierr) write_tim = MPI_WTIME() - stim call MPI_FILE_CLOSE(fh, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) call MPI_FILE_OPEN(MPI_COMM_WORLD, str, & & MPI_MODE_CREATE+MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr) call MPI_FILE_SEEK(fh, offset, MPI_SEEK_SET, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) stim = MPI_WTIME() call MPI_FILE_READ(fh, buf, SIZE, MPI_BYTE, status, ierr) read_tim = MPI_WTIME() - stim call MPI_FILE_CLOSE(fh, ierr) call MPI_ALLREDUCE(write_tim, new_write_tim, 1, & & MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr) call MPI_ALLREDUCE(read_tim, new_read_tim, 1, & & MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr) if (new_read_tim .lt. min_read_tim) then min_read_tim = new_read_tim end if if (new_write_tim .lt. min_write_tim) then min_write_tim = new_write_tim end if end do if (mynod .eq. 0) then read_bw = (SIZE*nprocs*1.0D0)/(min_read_tim*1000000.0D0) write_bw = (SIZE*nprocs*1.0D0)/(min_write_tim*1000000.0D0) print *, 'Write bandwidth without file sync = ', & & write_bw, ' Mbytes/sec' print *, 'Read bandwidth without prior file sync = ', & & read_bw, ' Mbytes/sec' end if min_read_tim = 10000000.0D0 min_write_tim = 10000000.0D0 flag = 0 do j=1, ntimes call MPI_FILE_OPEN(MPI_COMM_WORLD, str, & & MPI_MODE_CREATE+MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr) call MPI_FILE_SEEK(fh, offset, MPI_SEEK_SET, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) stim = MPI_WTIME() call MPI_FILE_WRITE(fh, buf, SIZE, MPI_BYTE, status, ierr) call MPI_FILE_SYNC(fh, ierr) write_tim = MPI_WTIME() - stim if (ierr .eq. MPI_ERR_UNKNOWN) then flag = 1 end if call MPI_FILE_CLOSE(fh, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) call MPI_FILE_OPEN(MPI_COMM_WORLD, str, & & MPI_MODE_CREATE+MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr) call MPI_FILE_SEEK(fh, offset, MPI_SEEK_SET, ierr) call MPI_BARRIER(MPI_COMM_WORLD, ierr) stim = MPI_WTIME() call MPI_FILE_READ(fh, buf, SIZE, MPI_BYTE, status, ierr) read_tim = MPI_WTIME() - stim call MPI_FILE_CLOSE(fh, ierr) call MPI_ALLREDUCE(write_tim, new_write_tim, 1, & & MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr) call MPI_ALLREDUCE(read_tim, new_read_tim, 1, & & MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr) if (new_read_tim .lt. min_read_tim) then min_read_tim = new_read_tim end if if (new_write_tim .lt. min_write_tim) then min_write_tim = new_write_tim end if end do if (mynod .eq. 0) then if (flag .eq. 1) then print *, 'MPI_FILE_SYNC returns error.' else read_bw = (SIZE*nprocs*1.0D0)/(min_read_tim*1000000.0D0) write_bw = (SIZE*nprocs*1.0D0)/(min_write_tim*1000000.0D0) print *, 'Write bandwidth including file sync = ', & & write_bw, ' Mbytes/sec' print *, 'Read bandwidth after file sync = ', & & read_bw, ' Mbytes/sec' end if end if call MPI_FINALIZE(ierr) end