forked from CODARcode/Example-Heat_Transfer
-
Notifications
You must be signed in to change notification settings - Fork 0
/
io_phdf5.F90
167 lines (122 loc) · 4.39 KB
/
io_phdf5.F90
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
! 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.
!
!
! Parallel HDF5 based I/O for the heat_transfer example
!
! (c) Oak Ridge National Laboratory, 2014
! Authors: Jeremy Logan, Norbert Podhorszki
!
module heat_io
contains
subroutine io_init()
use heat_vars
use HDF5
! Initialize FORTRAN interface
call h5open_f(err)
end subroutine io_init
subroutine io_finalize()
use heat_vars
use HDF5
! Close FORTRAN interface.
call h5close_f(err)
end subroutine io_finalize
subroutine io_write(tstep,curr)
use heat_vars
use HDF5
implicit none
include 'mpif.h'
integer, intent(in) :: tstep
integer, intent(in) :: curr
integer :: ndims
!integer*8, dimension(1:2) :: dims !local chunk
!integer*8, dimension(1:2) :: global_dims
INTEGER(HSSIZE_T), DIMENSION(3) :: offset
INTEGER(HSIZE_T), DIMENSION(3) :: dims
INTEGER(HSIZE_T), DIMENSION(3) :: global_dims
INTEGER(HSIZE_T), DIMENSION(3) :: max_dims
integer*8 io_size
INTEGER(HID_T) :: file_id
INTEGER(HID_T) :: dset_id
INTEGER(HID_T) :: dspace_id
INTEGER(HID_T) :: memspace
INTEGER(HID_T) :: plist_id
INTEGER :: comm, info
character (len=200) :: filename
character(2) :: mode = "w"
comm = MPI_COMM_WORLD
info = MPI_INFO_NULL
if (rank==0.and.tstep==1) then
print '("Writing: "," filename ",14x,"size(GB)",4x,"io_time(sec)",6x,"GB/s")'
endif
write(filename,'(a,".h5")') trim(outputfile)
!print '("rank ",i0," writes to: ",a)', rank, trim(filename)
call MPI_BARRIER(app_comm, err)
io_start_time = MPI_WTIME()
ndims = 3
dims(1) = ndx
dims(2) = ndy
dims(3) = 1
global_dims(1) = gndx
global_dims(2) = gndy
global_dims(3) = 1
max_dims(1) = gndx
max_dims(2) = gndy
max_dims(3) = H5S_UNLIMITED_F
offset(1) = offx
offset(2) = offy
offset(3) = 0
io_size = 11*4 + 2*8*ndx*ndy
call h5pcreate_f(H5P_FILE_ACCESS_F, plist_id, err)
call h5pset_fapl_mpio_f(plist_id, comm, info, err)
IF (tstep == 1) THEN
call h5fcreate_f (filename, H5F_ACC_TRUNC_F, file_id, err, access_prp = plist_id)
call h5pclose_f(plist_id, err)
call h5screate_simple_f(ndims, global_dims, dspace_id, err, max_dims)
call h5pcreate_f(H5P_DATASET_CREATE_F, plist_id, err)
call h5pset_chunk_f(plist_id, ndims, dims, err)
call h5dcreate_f(file_id, "T", H5T_NATIVE_DOUBLE, dspace_id, &
dset_id, err, plist_id)
call h5pclose_f(plist_id, err)
!!call h5sclose_f(dspace_id, err)
!!call h5dget_space_f(dset_id, dspace_id, err)
ELSE
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, err, &
access_prp = plist_id)
call h5pclose_f(plist_id, err)
call h5dopen_f(file_id, "T", dset_id, err)
global_dims(3) = tstep
call h5dset_extent_f(dset_id, global_dims, err)
offset(3) = tstep-1
call h5dget_space_f(dset_id, dspace_id, err)
END IF
call h5screate_simple_f(ndims, dims, memspace, err)
call h5sselect_hyperslab_f (dspace_id, H5S_SELECT_SET_F, &
offset, dims, err)
! For collective writes
call h5pcreate_f(H5P_DATASET_XFER_F, plist_id, err)
call h5pset_dxpl_mpio_f(plist_id, H5FD_MPIO_COLLECTIVE_F, err)
call h5dwrite_f(dset_id, H5T_NATIVE_DOUBLE, T(1:dims(1),1:dims(2),curr), &
dims, err, &
memspace, dspace_id, plist_id)
! For independent writes
!call h5dwrite_f(dset_id, H5T_NATIVE_DOUBLE, T(1:dims(1),1:dims(2),curr), &
! dims, err, &
! memspace, dspace_id)
! close hdf5 objects
call h5pclose_f(plist_id, err)
call h5dclose_f(dset_id, err)
call h5sclose_f(dspace_id, err)
call h5sclose_f(memspace, err)
call h5fclose_f(file_id, err)
call MPI_BARRIER(app_comm ,err)
io_end_time = MPI_WTIME()
io_total_time = io_end_time - io_start_time
sz = io_size * nproc/1024.d0/1024.d0/1024.d0 !size in GB
gbs = sz/io_total_time
if (rank==0) print '("Step ",i3,": ",a20,f12.4,2x,f12.3,2x,f12.3)', &
tstep,filename,sz,io_total_time,gbs
end subroutine io_write
end module heat_io