summaryrefslogtreecommitdiffstats
path: root/fortran/examples/rwdsetexample.f90
blob: b3cc42430bd10abc57a8fe7874a10bd55e775738 (plain)
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
!
! The following example shows how to write and read to/from an existing dataset. 
! It opens the file created in the previous example, obtains the dataset 
! identifier, writes the data to the dataset in the file, 
! then reads the dataset  to memory. 
!


     PROGRAM RWDSETEXAMPLE

     USE HDF5 ! This module contains all necessary modules 
        
     IMPLICIT NONE

     CHARACTER(LEN=8), PARAMETER :: filename = "dsetf.h5" ! File name
     CHARACTER(LEN=4), PARAMETER :: dsetname = "dset"     ! Dataset name

     INTEGER(HID_T) :: file_id       ! File identifier 
     INTEGER(HID_T) :: dset_id       ! Dataset identifier 

     INTEGER     ::   error ! Error flag
     INTEGER     ::  i, j

     INTEGER, DIMENSION(4,6) :: dset_data, data_out ! Data buffers
     INTEGER, DIMENSION(7) :: data_dims
     
     !
     ! Initialize the dset_data array.
     !
     do i = 1, 4
          do j = 1, 6
               dset_data(i,j) = (i-1)*6 + j;
          end do
     end do

     !
     ! Initialize FORTRAN interface.
     !
     CALL h5open_f(error) 

     !
     ! Open an existing file.
     !
     CALL h5fopen_f (filename, H5F_ACC_RDWR_F, file_id, error)

     !
     ! Open an existing dataset. 
     !
     CALL h5dopen_f(file_id, dsetname, dset_id, error)

     !
     ! Write the dataset.
     !
     data_dims(1) = 4
     data_dims(2) = 6 
     CALL h5dwrite_f(dset_id, H5T_NATIVE_INTEGER, dset_data, data_dims, error)

     !
     ! Read the dataset.
     !
     CALL h5dread_f(dset_id, H5T_NATIVE_INTEGER, data_out, data_dims, error)

     !
     ! Close the dataset.
     !
     CALL h5dclose_f(dset_id, error)

     !
     ! Close the file.
     !
     CALL h5fclose_f(file_id, error)
     
     !
     ! Close FORTRAN interface.
     !
     CALL h5close_f(error)

     END PROGRAM RWDSETEXAMPLE