-
Notifications
You must be signed in to change notification settings - Fork 39
/
Copy pathHdfRoutines.F90
498 lines (346 loc) · 14.1 KB
/
HdfRoutines.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
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
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: HdfRoutines.F90 !
! CONTAINS: subroutine hdf_read_serial_1d !
! !
! PURPOSE: I/O routines. !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine HdfCreateBlankFile(filename)
use hdf5
implicit none
character*30,intent(in) :: filename
integer(HID_T) :: file_id
integer :: hdf_error
call h5fcreate_f(filename, H5F_ACC_TRUNC_F, file_id, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfCreateBlankFile
!====================================================================
subroutine HdfSerialWriteRealScalar(dsetname,filename,n)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
real, intent(in) :: n
integer(HID_T) :: file_id
integer(HID_T) :: dset, filespace
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
logical :: fileexists
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=1
call h5lexists_f(file_id,dsetname,fileexists,hdf_error)
if(fileexists) call h5ldelete_f(file_id,dsetname,hdf_error)
call h5screate_simple_f(1, dims, &
& filespace, hdf_error)
call h5dcreate_f(file_id, dsetname, H5T_NATIVE_DOUBLE, &
& filespace, dset, hdf_error)
call h5dwrite_f(dset, H5T_NATIVE_DOUBLE, &
& n, dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5sclose_f(filespace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialWriteRealScalar
!====================================================================
subroutine HdfSerialWriteIntScalar(dsetname,filename,n)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
integer, intent(in) :: n
integer(HID_T) :: file_id
integer(HID_T) :: dset, filespace
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
logical :: fileexists
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=1
call h5lexists_f(file_id,dsetname,fileexists,hdf_error)
if(fileexists) call h5ldelete_f(file_id,dsetname,hdf_error)
call h5screate_simple_f(1, dims, &
& filespace, hdf_error)
call h5dcreate_f(file_id, dsetname, H5T_NATIVE_INTEGER, &
& filespace, dset, hdf_error)
call h5dwrite_f(dset, H5T_NATIVE_INTEGER, &
& n, dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5sclose_f(filespace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialWriteIntScalar
!====================================================================
subroutine HdfSerialWriteReal1D(dsetname,filename,var,sz)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
integer, intent(in) :: sz
real, dimension(sz), intent(in) :: var
integer(HID_T) :: file_id
integer(HID_T) :: dset, filespace
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
logical :: fileexists
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=sz
call h5screate_simple_f(1, dims, &
& filespace, hdf_error)
call h5lexists_f(file_id,dsetname,fileexists,hdf_error)
if(fileexists) call h5ldelete_f(file_id,dsetname,hdf_error)
call h5dcreate_f(file_id, dsetname, H5T_NATIVE_DOUBLE, &
& filespace, dset, hdf_error)
call h5dwrite_f(dset, H5T_NATIVE_DOUBLE, &
& var(1:sz), dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5sclose_f(filespace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialWriteReal1D
!====================================================================
subroutine HdfWriteReal2D(dsetname,filename,var)
use mpih
use param
use hdf5
use decomp_2d, only: xstart,xend
implicit none
character*30,intent(in) :: dsetname,filename
real, intent(in) :: var(xstart(2):xend(2) &
& ,xstart(3):xend(3))
integer(HID_T) :: file_id
integer(HID_T) :: filespace
integer(HID_T) :: slabspace
integer(HID_T) :: memspace
integer(HID_T) :: dset
integer(HID_T) :: plist_id
integer(HSIZE_T) :: dims(2)
integer(HSIZE_T), dimension(2) :: data_count
integer(HSSIZE_T), dimension(2) :: data_offset
integer :: hdf_error, ndims
integer :: comm, info
!RO Sort out MPI definitions
comm = MPI_COMM_WORLD
info = MPI_INFO_NULL
ndims = 2
dims(1)=nym
dims(2)=nzm
call h5screate_simple_f(ndims, dims, &
& filespace, hdf_error)
data_count(1) = xend(2)-xstart(2)+1
data_count(2) = xend(3)-xstart(3)+1
data_offset(1) = xstart(2)-1
data_offset(2) = xstart(3)-1
call h5pcreate_f(H5P_FILE_ACCESS_F, plist_id, &
& hdf_error)
call h5pset_fapl_mpio_f(plist_id, comm, info, &
& hdf_error)
call h5fcreate_f(filename, H5F_ACC_TRUNC_F, file_id, &
& hdf_error, access_prp=plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dcreate_f(file_id, dsetname, H5T_NATIVE_DOUBLE, &
& filespace, dset, hdf_error)
call h5screate_simple_f(ndims, data_count, memspace, hdf_error)
call h5dget_space_f(dset, slabspace, hdf_error)
call h5sselect_hyperslab_f (slabspace, H5S_SELECT_SET_F, &
& data_offset, data_count, hdf_error)
call h5pcreate_f(H5P_DATASET_XFER_F, plist_id, hdf_error)
call h5pset_dxpl_mpio_f(plist_id, H5FD_MPIO_COLLECTIVE_F, &
& hdf_error)
call h5dwrite_f(dset, H5T_NATIVE_DOUBLE, &
& var(xstart(2):xend(2),xstart(3):xend(3)), dims, &
& hdf_error, file_space_id = slabspace, mem_space_id = memspace, &
& xfer_prp = plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5sclose_f(memspace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfWriteReal2D
!====================================================================
subroutine HdfSerialReadRealScalar(dsetname,filename,n)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
real, intent(out) :: n
integer(HID_T) :: file_id
integer(HID_T) :: dset
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=1
call h5dopen_f(file_id, dsetname, dset, hdf_error)
call h5dread_f(dset, H5T_NATIVE_DOUBLE, &
& n, dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialReadRealScalar
!====================================================================
subroutine HdfSerialReadIntScalar(dsetname,filename,n)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
integer, intent(out) :: n
integer(HID_T) :: file_id
integer(HID_T) :: dset
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=1
call h5dopen_f(file_id, dsetname, dset, hdf_error)
call h5dread_f(dset, H5T_NATIVE_INTEGER, &
& n, dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialReadIntScalar
!====================================================================
subroutine HdfSerialReadReal1D(dsetname,filename,var,sz)
use hdf5
implicit none
character*30,intent(in) :: dsetname,filename
integer, intent(in) :: sz
real, dimension(sz), intent(out) :: var
integer(HID_T) :: file_id
integer(HID_T) :: dset
integer :: hdf_error
integer(HSIZE_T) :: dims(1)
call h5fopen_f(filename, H5F_ACC_RDWR_F, file_id, hdf_error)
dims(1)=sz
call h5dopen_f(file_id, dsetname, dset, hdf_error)
call h5dread_f(dset, H5T_NATIVE_DOUBLE, &
& var, dims, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfSerialReadReal1D
!====================================================================
subroutine HdfReadRealHalo3D(filnam1,qua,st2,en2,st3,en3,n1o,n2o,n3o)
use param
use mpih
use hdf5
implicit none
integer, intent(in) :: n1o,n2o,n3o,st2,en2,st3,en3
real, intent(out), dimension(1:n3o,st2-lvlhalo:en2+lvlhalo, &
st3-lvlhalo:en3+lvlhalo)::qua
integer hdf_error
integer(HID_T) :: file_id
integer(HID_T) :: slabspace
integer(HID_T) :: memspace
integer(HID_T) :: dset_qua
integer(HSIZE_T) :: dims(3)
integer(HID_T) :: plist_id
integer(HSIZE_T), dimension(3) :: data_count
integer(HSSIZE_T), dimension(3) :: data_offset
integer :: comm, info
integer :: ndims
! real, intent(in), dimension(1:n3,xstart(2)-1:xend(2)+1, &
! & xstart(3)-1:xend(3)+1)::qua
character*30,intent(in) :: filnam1
!RO Sort out MPI definitions
comm = MPI_COMM_WORLD
info = MPI_INFO_NULL
!RO Set offsets and element counts
ndims = 3
dims(1)=n3o
dims(2)=n2o-1
dims(3)=n1o-1
data_count(1) = n3o
data_count(2) = en2-st2+1
data_count(3) = en3-st3+1
data_offset(1) = 0
data_offset(2) = st2-1
data_offset(3) = st3-1
call h5pcreate_f(H5P_FILE_ACCESS_F, plist_id, &
hdf_error)
call h5pset_fapl_mpio_f(plist_id, comm, info, &
& hdf_error)
call h5fopen_f(filnam1, H5F_ACC_RDONLY_F, file_id, &
& hdf_error, access_prp=plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dopen_f(file_id, 'var', &
& dset_qua, hdf_error)
call h5screate_simple_f(ndims, data_count, memspace, hdf_error)
call h5dget_space_f(dset_qua, slabspace, hdf_error)
call h5sselect_hyperslab_f (slabspace, H5S_SELECT_SET_F, &
& data_offset, data_count, hdf_error)
call h5pcreate_f(H5P_DATASET_XFER_F, plist_id, hdf_error)
call h5pset_dxpl_mpio_f(plist_id, H5FD_MPIO_COLLECTIVE_F, &
& hdf_error)
call h5dread_f(dset_qua, H5T_NATIVE_DOUBLE, &
& qua(1:n3o,st2:en2,st3:en3), dims, &
& hdf_error, file_space_id = slabspace, mem_space_id = memspace, &
& xfer_prp = plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dclose_f(dset_qua, hdf_error)
call h5sclose_f(memspace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfReadRealHalo3D
!====================================================================
subroutine HdfWriteRealHalo3D(filnam1,qua)
use param
use mpih
use hdf5
use decomp_2d, only: xstart,xend
implicit none
integer hdf_error
integer(HID_T) :: file_id
integer(HID_T) :: filespace
integer(HID_T) :: slabspace
integer(HID_T) :: memspace
integer(HID_T) :: dset
integer(HSIZE_T) :: dims(3)
integer(HID_T) :: plist_id
integer(HSIZE_T), dimension(3) :: data_count
integer(HSSIZE_T), dimension(3) :: data_offset
integer :: comm, info
integer :: ndims
real, intent(in), dimension(1:nx, &
xstart(2)-lvlhalo:xend(2)+lvlhalo, &
xstart(3)-lvlhalo:xend(3)+lvlhalo)::qua
character*30,intent(in) :: filnam1
!RO Sort out MPI definitions
comm = MPI_COMM_WORLD
info = MPI_INFO_NULL
!RO Set offsets and element counts
ndims = 3
dims(1)=nx
dims(2)=nym
dims(3)=nzm
call h5screate_simple_f(ndims, dims, &
& filespace, hdf_error)
data_count(1) = nx
data_count(2) = xend(2)-xstart(2)+1
data_count(3) = xend(3)-xstart(3)+1
data_offset(1) = 0
data_offset(2) = xstart(2)-1
data_offset(3) = xstart(3)-1
call h5pcreate_f(H5P_FILE_ACCESS_F, plist_id, &
hdf_error)
call h5pset_fapl_mpio_f(plist_id, comm, info, &
& hdf_error)
call h5fcreate_f(filnam1, H5F_ACC_TRUNC_F, file_id, &
& hdf_error, access_prp=plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dcreate_f(file_id, 'var', H5T_NATIVE_DOUBLE, &
& filespace, dset, hdf_error)
call h5screate_simple_f(ndims, data_count, memspace, hdf_error)
call h5dget_space_f(dset, slabspace, hdf_error)
call h5sselect_hyperslab_f (slabspace, H5S_SELECT_SET_F, &
& data_offset, data_count, hdf_error)
call h5pcreate_f(H5P_DATASET_XFER_F, plist_id, hdf_error)
call h5pset_dxpl_mpio_f(plist_id, H5FD_MPIO_COLLECTIVE_F, &
& hdf_error)
call h5dwrite_f(dset, H5T_NATIVE_DOUBLE, &
& qua(1:nx,xstart(2):xend(2),xstart(3):xend(3)), dims, &
& hdf_error, file_space_id = slabspace, mem_space_id = memspace, &
& xfer_prp = plist_id)
call h5pclose_f(plist_id, hdf_error)
call h5dclose_f(dset, hdf_error)
call h5sclose_f(memspace, hdf_error)
call h5fclose_f(file_id, hdf_error)
end subroutine HdfWriteRealHalo3D
!====================================================================
subroutine HdfStart
use hdf5
implicit none
integer :: hdf_error
call h5open_f(hdf_error)
end subroutine HdfStart
!====================================================================
subroutine HdfClose
use hdf5
implicit none
integer :: hdf_error
call h5close_f(hdf_error)
end subroutine HdfClose