summaryrefslogtreecommitdiffstats
path: root/HDF5Examples/FORTRAN/H5PAR/ph5_f90_hyperslab_by_chunk.F90
blob: c74e55dbdcf392302fbde9e6b5ff66ed5953972b (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
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
!
! Number of processes is assumed to be 4
!
     PROGRAM DATASET_BY_CHUNK

     USE HDF5 ! This module contains all necessary modules
!     USE MPI

     IMPLICIT NONE

     include 'mpif.h'
     CHARACTER(LEN=11), PARAMETER :: filename = "sds_chnk.h5"  ! File name
     CHARACTER(LEN=8), PARAMETER :: dsetname = "IntArray" ! Dataset name

     INTEGER(HID_T) :: file_id       ! File identifier
     INTEGER(HID_T) :: dset_id       ! Dataset identifier
     INTEGER(HID_T) :: filespace     ! Dataspace identifier in file
     INTEGER(HID_T) :: memspace      ! Dataspace identifier in memory
     INTEGER(HID_T) :: plist_id      ! Property list identifier

     INTEGER(HSIZE_T), DIMENSION(2) :: dimsf = (/4,8/) ! Dataset dimensions
                                                       ! in the file.
     INTEGER(HSIZE_T), DIMENSION(2) :: dimsfi = (/4,8/)
     INTEGER(HSIZE_T), DIMENSION(2) :: chunk_dims = (/2,4/) ! Chunks dimensions

     INTEGER(HSIZE_T),  DIMENSION(2) :: count
     INTEGER(HSSIZE_T), DIMENSION(2) :: offset
     INTEGER(HSIZE_T),  DIMENSION(2) :: stride
     INTEGER(HSIZE_T),  DIMENSION(2) :: block

     INTEGER, ALLOCATABLE :: data (:,:)  ! Data to write
     INTEGER :: rank = 2 ! Dataset rank

     INTEGER :: error, error_n  ! Error flags
     !
     ! MPI definitions and calls.
     !
     INTEGER :: mpierror       ! MPI error flag
     INTEGER :: comm, info
     INTEGER :: mpi_size, mpi_rank

     comm = MPI_COMM_WORLD
     info = MPI_INFO_NULL

     CALL MPI_INIT(mpierror)
     CALL MPI_COMM_SIZE(comm, mpi_size, mpierror)
     CALL MPI_COMM_RANK(comm, mpi_rank, mpierror)
     ! Quit if mpi_size is not 4
     if (mpi_size .NE. 4) then
        write(*,*) 'This example is set up to use only 4 processes'
        write(*,*) 'Quitting....'
        goto 100
     endif

     !
     ! Initialize HDF5 library and Fortran interfaces.
     !
     CALL h5open_f(error)

     !
     ! Setup file access property list with parallel I/O access.
     !
     CALL h5pcreate_f(H5P_FILE_ACCESS_F, plist_id, error)
     CALL h5pset_fapl_mpio_f(plist_id, comm, info, error)

     !
     ! Create the file collectively.
     !
     CALL h5fcreate_f(filename, H5F_ACC_TRUNC_F, file_id, error, access_prp = plist_id)
     CALL h5pclose_f(plist_id, error)
     !
     ! Create the data space for the  dataset.
     !
     CALL h5screate_simple_f(rank, dimsf, filespace, error)
     CALL h5screate_simple_f(rank, chunk_dims, memspace, error)

     !
     ! Create chunked dataset.
     !
     CALL h5pcreate_f(H5P_DATASET_CREATE_F, plist_id, error)
     CALL h5pset_chunk_f(plist_id, rank, chunk_dims, error)
     CALL h5dcreate_f(file_id, dsetname, H5T_NATIVE_INTEGER, filespace, &
                      dset_id, error, plist_id)
     CALL h5sclose_f(filespace, error)
     !
     ! Each process defines dataset in memory and writes it to the hyperslab
     ! in the file.
     !
     stride(1) = 1
     stride(2) = 1
     count(1) =  1
     count(2) =  1
     block(1) = chunk_dims(1)
     block(2) = chunk_dims(2)
     if (mpi_rank .EQ. 0) then
         offset(1) = 0
         offset(2) = 0
     endif
     if (mpi_rank .EQ. 1) then
         offset(1) = chunk_dims(1)
         offset(2) = 0
     endif
     if (mpi_rank .EQ. 2) then
         offset(1) = 0
         offset(2) = chunk_dims(2)
     endif
     if (mpi_rank .EQ. 3) then
         offset(1) = chunk_dims(1)
         offset(2) = chunk_dims(2)
     endif
     !
     ! Select hyperslab in the file.
     !
     CALL h5dget_space_f(dset_id, filespace, error)
     CALL h5sselect_hyperslab_f (filespace, H5S_SELECT_SET_F, offset, count, error, &
                                 stride, block)
     !
     ! Initialize data buffer with trivial data.
     !
     ALLOCATE (data(chunk_dims(1),chunk_dims(2)))
     data = mpi_rank + 1
     !
     ! Create property list for collective dataset write
     !
     CALL h5pcreate_f(H5P_DATASET_XFER_F, plist_id, error)
     CALL h5pset_dxpl_mpio_f(plist_id, H5FD_MPIO_COLLECTIVE_F, error)

     !
     ! Write the dataset collectively.
     !
     CALL h5dwrite_f(dset_id, H5T_NATIVE_INTEGER, data, dimsfi, error, &
                     file_space_id = filespace, mem_space_id = memspace, xfer_prp = plist_id)
     !
     ! Write the dataset independently.
     !
!    CALL h5dwrite_f(dset_id, H5T_NATIVE_INTEGER, data, dimsfi,error, &
!                     file_space_id = filespace, mem_space_id = memspace)
     !
     ! Deallocate data buffer.
     !
     DEALLOCATE(data)

     !
     ! Close dataspaces.
     !
     CALL h5sclose_f(filespace, error)
     CALL h5sclose_f(memspace, error)
     !
     ! Close the dataset.
     !
     CALL h5dclose_f(dset_id, error)
     !
     ! Close the property list.
     !
     CALL h5pclose_f(plist_id, error)
     !
     ! Close the file.
     !
     CALL h5fclose_f(file_id, error)

     !
     ! Close FORTRAN interfaces and HDF5 library.
     !
     CALL h5close_f(error)
     IF(mpi_rank.EQ.0) WRITE(*,'(A)') "PHDF5 example finished with no errors"
100  continue
     CALL MPI_FINALIZE(mpierror)

     END PROGRAM DATASET_BY_CHUNK