1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
|
! This file created from f77/io/atomicityf.f with f77tof90
!
! Copyright (C) by Argonne National Laboratory
! See COPYRIGHT in top-level directory
!
program main
use mpi
integer (kind=MPI_OFFSET_KIND) disp
! tests whether atomicity semantics are satisfied for overlapping accesses
! in atomic mode. The probability of detecting errors is higher if you run
! it on 8 or more processes.
! This is a version of the test in romio/test/atomicity.c .
integer BUFSIZE
parameter (BUFSIZE=10000)
integer writebuf(BUFSIZE), readbuf(BUFSIZE)
integer i, mynod, nprocs, len, ierr, errs, toterrs
character*50 filename
integer newtype, fh, info, status(MPI_STATUS_SIZE)
errs = 0
call MTest_Init(ierr)
call MPI_Comm_rank(MPI_COMM_WORLD, mynod, ierr )
call MPI_Comm_size(MPI_COMM_WORLD, nprocs, ierr )
! Unlike the C version, we fix the filename because of the difficulties
! in accessing the command line from different Fortran environments
filename = "testfile.txt"
! test atomicity of contiguous accesses
! initialize file to all zeros
if (mynod .eq. 0) then
call MPI_File_delete(filename, MPI_INFO_NULL, ierr )
call MPI_File_open(MPI_COMM_SELF, filename, MPI_MODE_CREATE + &
& MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr )
do i=1, BUFSIZE
writebuf(i) = 0
enddo
call MPI_File_write(fh, writebuf, BUFSIZE, MPI_INTEGER, status, &
& ierr)
call MPI_File_close(fh, ierr )
endif
call MPI_Barrier(MPI_COMM_WORLD, ierr )
do i=1, BUFSIZE
writebuf(i) = 10
readbuf(i) = 20
enddo
call MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_CREATE + &
& MPI_MODE_RDWR, MPI_INFO_NULL, fh, ierr )
! set atomicity to true
call MPI_File_set_atomicity(fh, .true., ierr)
if (ierr .ne. MPI_SUCCESS) then
print *, "Atomic mode not supported on this file system."
call MPI_Abort(MPI_COMM_WORLD, 1, ierr )
endif
call MPI_Barrier(MPI_COMM_WORLD, ierr )
! process 0 writes and others concurrently read. In atomic mode,
! the data read must be either all old values or all new values; nothing
! in between.
if (mynod .eq. 0) then
call MPI_File_write(fh, writebuf, BUFSIZE, MPI_INTEGER, status, &
& ierr)
else
call MPI_File_read(fh, readbuf, BUFSIZE, MPI_INTEGER, status, &
& ierr )
if (ierr .eq. MPI_SUCCESS) then
if (readbuf(1) .eq. 0) then
! the rest must also be 0
do i=2, BUFSIZE
if (readbuf(i) .ne. 0) then
errs = errs + 1
print *, "(contig)Process ", mynod, ": readbuf(", i &
& ,") is ", readbuf(i), ", should be 0"
call MPI_Abort(MPI_COMM_WORLD, 1, ierr )
endif
enddo
else if (readbuf(1) .eq. 10) then
! the rest must also be 10
do i=2, BUFSIZE
if (readbuf(i) .ne. 10) then
errs = errs + 1
print *, "(contig)Process ", mynod, ": readbuf(", i &
& ,") is ", readbuf(i), ", should be 10"
call MPI_Abort(MPI_COMM_WORLD, 1, ierr )
endif
enddo
else
errs = errs + 1
print *, "(contig)Process ", mynod, ": readbuf(1) is ", &
& readbuf(1), ", should be either 0 or 10"
endif
endif
endif
call MPI_File_close( fh, ierr )
call MPI_Barrier( MPI_COMM_WORLD, ierr )
! repeat the same test with a noncontiguous filetype
call MPI_Type_vector(BUFSIZE, 1, 2, MPI_INTEGER, newtype, ierr)
call MPI_Type_commit(newtype, ierr )
call MPI_Info_create(info, ierr )
! I am setting these info values for testing purposes only. It is
! better to use the default values in practice. */
call MPI_Info_set(info, "ind_rd_buffer_size", "1209", ierr )
call MPI_Info_set(info, "ind_wr_buffer_size", "1107", ierr )
if (mynod .eq. 0) then
call MPI_File_delete(filename, MPI_INFO_NULL, ierr )
call MPI_File_open(MPI_COMM_SELF, filename, MPI_MODE_CREATE + &
& MPI_MODE_RDWR, info, fh, ierr )
do i=1, BUFSIZE
writebuf(i) = 0
enddo
disp = 0
call MPI_File_set_view(fh, disp, MPI_INTEGER, newtype, "native" &
& ,info, ierr)
call MPI_File_write(fh, writebuf, BUFSIZE, MPI_INTEGER, status, &
& ierr )
call MPI_File_close( fh, ierr )
endif
call MPI_Barrier( MPI_COMM_WORLD, ierr )
do i=1, BUFSIZE
writebuf(i) = 10
readbuf(i) = 20
enddo
call MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_CREATE + &
& MPI_MODE_RDWR, info, fh, ierr )
call MPI_File_set_atomicity(fh, .true., ierr)
disp = 0
call MPI_File_set_view(fh, disp, MPI_INTEGER, newtype, "native", &
& info, ierr )
call MPI_Barrier(MPI_COMM_WORLD, ierr )
if (mynod .eq. 0) then
call MPI_File_write(fh, writebuf, BUFSIZE, MPI_INTEGER, status, &
& ierr )
else
call MPI_File_read(fh, readbuf, BUFSIZE, MPI_INTEGER, status, &
& ierr )
if (ierr .eq. MPI_SUCCESS) then
if (readbuf(1) .eq. 0) then
do i=2, BUFSIZE
if (readbuf(i) .ne. 0) then
errs = errs + 1
print *, "(noncontig)Process ", mynod, ": readbuf(" &
& , i,") is ", readbuf(i), ", should be 0"
call MPI_Abort(MPI_COMM_WORLD, 1, ierr )
endif
enddo
else if (readbuf(1) .eq. 10) then
do i=2, BUFSIZE
if (readbuf(i) .ne. 10) then
errs = errs + 1
print *, "(noncontig)Process ", mynod, ": readbuf(" &
& , i,") is ", readbuf(i), ", should be 10"
call MPI_Abort(MPI_COMM_WORLD, 1, ierr )
endif
enddo
else
errs = errs + 1
print *, "(noncontig)Process ", mynod, ": readbuf(1) is " &
& ,readbuf(1), ", should be either 0 or 10"
endif
endif
endif
call MPI_File_close( fh, ierr )
call MPI_Barrier(MPI_COMM_WORLD, ierr )
call MPI_Type_free(newtype, ierr )
call MPI_Info_free(info, ierr )
call MTest_Finalize(errs)
end
|