File: genarray.F90

package info (click to toggle)
adios 1.13.1-31
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 23,692 kB
  • sloc: ansic: 133,236; f90: 8,791; sh: 7,779; python: 7,648; xml: 3,793; makefile: 2,996; cpp: 2,340; java: 626; sed: 16; perl: 8
file content (313 lines) | stat: -rw-r--r-- 10,760 bytes parent folder | download | duplicates (3)
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
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313

!  ADIOS is freely available under the terms of the BSD license described
!  in the COPYING file in the top level directory of this source distribution.
!
!  Copyright (c) 2008 - 2009.  UT-BATTELLE, LLC. All rights reserved.
!

!
!  GENARRAY
!
!  Write an ADIOS BP file from many processor for test purposes.
!
!  nx * ny * nz     processes write a 3D array, where each process writes an
!  ndx * ndy * ndz  piece with filling with its rank as integer (4 bytes) value
!
!
! (c) Oak Ridge National Laboratory, 2009
! Author: Norbert Podhorszki
!
module genarray_comm
    ! arguments
    character(len=256) :: outputfile, inputfile
    integer :: npx, npy, npz  ! # of processors in x-y-z direction
    integer :: ndx, ndy, ndz  ! size of array per processor
    integer :: timesteps      ! number of timesteps to write
    integer :: sleeptime      ! time to sleep between time steps
    logical :: common_size    ! .true.  if common local sizes are given as argument
                              ! .false. if we have to read sizes from a file

    integer :: gndx, gndy, gndz  ! size of the global array
    integer :: offx,offy,offz    ! offsets of local array in the global array

    real*8, dimension(:,:,:), allocatable :: double_xyz

    ! MPI variables
    integer :: group_comm
    integer :: rank, nproc
    integer :: ierr

    ! ADIOS variables
    character (len=200) :: group
    character (len=200) :: filename
    !character (len=6)   :: nprocstr
    integer*8 :: handle, total_size, group_size, adios_totalsize

    real*8 :: start_time, end_time, total_time,gbs,sz
    real*8 :: io_start_time, io_end_time, io_total_time


end module genarray_comm


program genarray
    use genarray_comm
    use adios_write_mod
    implicit none
    include 'mpif.h'

    call MPI_Init (ierr)
    call MPI_Comm_dup (MPI_COMM_WORLD, group_comm, ierr)
    call MPI_Comm_rank (MPI_COMM_WORLD, rank, ierr)
    call MPI_Comm_size (group_comm, nproc , ierr)

    call adios_init ("genarray3d.xml", group_comm, ierr)
    !call MPI_Barrier (group_comm, ierr)

    call processArgs()
    if (rank == 0) then
        print *,"Output file(s): "//trim(outputfile)//".<step>.bp"
        print '(" Process number        : ",i0," x ",i0," x ",i0)', npx,npy,npz
        if (common_size) then
            print '(" Array size per process: ",i0," x ",i0," x ",i0)', ndx,ndy,ndz
        else
            print *," Array sizes per processes taken from file: "//trim(inputfile)
        endif

        if (nproc .ne. npx*npy*npz) then
            print '(" Error: Number of processors ",i0,"does not match ndx*ndy*ndz=",i0)', nproc, npx*npy*npz
            call exit(1)
        endif
    endif

    !write (*,*) 'rank ', rank, "init completed"
    !write (nprocstr,'(i0)') nproc

    call determineLocalSize()
    call determineGlobalSize()
    call determineOffsets()
    call generateLocalArray()

    call writeArray()
    ! Terminate
    call MPI_Barrier (MPI_COMM_WORLD, ierr)
    call adios_finalize (rank, ierr)
    call MPI_Finalize (ierr)
end program genarray


!!***************************
subroutine determineLocalSize()
    use genarray_comm
    implicit none
    if (common_size) then
       ! we are done since we know them from argument
    else
       ! have to read from file
       print *, "To be implemented: read sizes from file 1"
       call exit(2)
    endif
end subroutine determineLocalSize

!!***************************
subroutine determineGlobalSize()
    use genarray_comm
    implicit none
    if (common_size) then
        gndx = npx * ndx
        gndy = npy * ndy
        gndz = npz * ndz
    else
       ! have to read from file
       print *, "To be implemented: read sizes from file 2"
       call exit(2)
    endif
end subroutine determineGlobalSize

!!***************************
subroutine determineOffsets()
    use genarray_comm
    implicit none
    integer :: posx, posy, posz ! position index in the array
    if (common_size) then
        posx = mod(rank, npx)     ! 1st dim easy: 0, npx, 2npx... are in the same X position
        posy = mod(rank/npx, npy) ! 2nd dim: (0, npx-1) have the same dim (so divide with npx first)
        posz = rank/(npx*npy)     ! 3rd dim: npx*npy processes belong into one dim
        offx = posx * ndx
        offy = posy * ndy
        offz = posz * ndz
    else
       ! have to read from file
       print *, "To be implemented: read sizes from file 3"
       call exit(2)
    endif
end subroutine determineOffsets


!!***************************
subroutine generateLocalArray()
    use genarray_comm
    implicit none
    integer :: i,j,k
    allocate( double_xyz(1:ndx, 1:ndy, 1:ndz) )
    do k=1,ndz
        do j=1,ndy
            do i=1,ndx
                double_xyz(i,j,k) = 1.0d0*rank
            enddo
        enddo
    enddo
end subroutine generateLocalArray


!!***************************
subroutine writeArray()
    use genarray_comm
    use adios_write_mod
    implicit none
    integer*8 adios_handle, adios_groupsize
    integer adios_err
    integer :: tstep
    character(2) :: mode = "w"
    character(len=256) :: outfilename
    include 'mpif.h'


    if (rank==0) print '("Writing: "," filename ",14x,"size(GB)",4x,"io_time(sec)",6x,"GB/s")'
    do tstep=1,timesteps
        !if (tstep > 1) mode = "a"
        double_xyz = tstep + double_xyz
        call MPI_BARRIER(MPI_COMM_WORLD,ierr)
        io_start_time = MPI_WTIME()
        group = "genarray"
        write (outfilename,'(a,".",i3.3,".bp")') trim(outputfile),tstep
        call adios_open (adios_handle, group, outfilename, mode, group_comm, adios_err)
        adios_groupsize = 0
        adios_groupsize = adios_groupsize + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 4_8 &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz) &
                + 8_8 * (ndx) * (ndy) * (ndz)
        call adios_group_size (adios_handle, adios_groupsize, adios_totalsize, adios_err)
        call adios_write (adios_handle, "/dimensions/gndx", gndx, adios_err)
        call adios_write (adios_handle, "/dimensions/gndy", gndy, adios_err)
        call adios_write (adios_handle, "/dimensions/gndz", gndz, adios_err)
        call adios_write (adios_handle, "/info/nproc", nproc, adios_err)
        call adios_write (adios_handle, "/info/npx", npx, adios_err)
        call adios_write (adios_handle, "/info/npy", npy, adios_err)
        call adios_write (adios_handle, "/info/npz", npz, adios_err)
        call adios_write (adios_handle, "/aux/offx", offx, adios_err)
        call adios_write (adios_handle, "/aux/offy", offy, adios_err)
        call adios_write (adios_handle, "/aux/offz", offz, adios_err)
        call adios_write (adios_handle, "/aux/ndx", ndx, adios_err)
        call adios_write (adios_handle, "/aux/ndy", ndy, adios_err)
        call adios_write (adios_handle, "/aux/ndz", ndz, adios_err)
        call adios_write (adios_handle, "/var/var1", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var2", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var3", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var4", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var5", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var6", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var7", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var8", double_xyz, adios_err)
        call adios_write (adios_handle, "/var/var9", double_xyz, adios_err)
        call adios_close (adios_handle, adios_err)
        call MPI_BARRIER(MPI_COMM_WORLD,ierr)
        io_end_time = MPI_WTIME()
        io_total_time = io_end_time - io_start_time
        sz = adios_totalsize * nproc/1024.d0/1024.d0/1024.d0 !size in GB
        gbs = sz/io_total_time
        if (rank==0) print '("Writing: ",a20,d12.2,2x,d12.2,2x,d12.3)', outfilename,sz,io_total_time,gbs
        if (tstep<timesteps) call sleep(sleeptime)
     end do
end subroutine writeArray


!!***************************
subroutine usage()
    print *, "Usage: genarray  output N  M  K  [infile|nx  ny  nz timesteps sleeptime]"
    print *, "output: name of output file"
    print *, "N:      number of processes in X dimension"
    print *, "M:      number of processes in Y dimension"
    print *, "K:      number of processes in Z dimension"
    print *, "nx:     local array size in X dimension per processor"
    print *, "ny:     local array size in Y dimension per processor"
    print *, "nz:     local array size in Z dimension per processor"
    print *, "infile: file that describes nx ny nz for each processor"
    print *, "timesteps: the total number of timesteps to output" 
    print *, "sleeptime: the time to sleep (s)"
end subroutine usage

!!***************************
subroutine processArgs()
    use genarray_comm

#ifndef __GFORTRAN__
!#ifndef __GNUC__
    interface
         integer function iargc()
         end function iargc
    end interface
!#endif
#endif

    character(len=256) :: npx_str, npy_str, npz_str, ndx_str, ndy_str, ndz_str
    character(len=256) :: time_str,sleep_str
    integer :: numargs

    !! process arguments
    numargs = iargc()
    !print *,"Number of arguments:",numargs
    if ( numargs < 5 ) then
        call usage()
        call exit(1)
    endif
    call getarg(1, outputfile)
    call getarg(2, npx_str)
    call getarg(3, npy_str)
    call getarg(4, npz_str)
    read (npx_str,'(i5)') npx
    read (npy_str,'(i5)') npy
    read (npz_str,'(i5)') npz
    if ( numargs == 5 ) then
        call getarg(5, inputfile)
        ndx = 0
        ndy = 0
        ndz = 0
        common_size = .false.
    else if (numargs == 9) then
        call getarg(5, ndx_str)
        call getarg(6, ndy_str)
        call getarg(7, ndz_str)
        read (ndx_str,'(i12)') ndx
        read (ndy_str,'(i12)') ndy
        read (ndz_str,'(i12)') ndz
        inputfile=char(0)
        common_size = .true.
        call getarg(8, time_str)
        call getarg(9, sleep_str)
        read (time_str,'(i6)') timesteps
        read (sleep_str,'(i6)') sleeptime
    else
        call usage()
        call exit(1)
    endif

end subroutine processArgs