summaryrefslogtreecommitdiffstats
path: root/HDF5Examples/FORTRAN/H5D/h5ex_d_compact.F90
blob: c2e5aafd87f48d4d8d6177b5fddc2f7c6476ab38 (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
! ************************************************************
!
!  This example shows how to read and write data to a compact
!  dataset.  The program first writes integers to a compact
!  dataset with dataspace dimensions of DIM0xDIM1, then
!  closes the file.  Next, it reopens the file, reads back
!  the data, and outputs it to the screen.
!
!  This file is intended for use with HDF5 Library version 1.8
!
! ************************************************************

! An optional example for determining the correct HDF5 version 
! for picking the correct HDF5 API parameters. This is not 
! part of the HDF5 library.
#include "h5_version.h"

PROGRAM main

  USE HDF5

  IMPLICIT NONE

  CHARACTER(LEN=17), PARAMETER :: filename = "h5ex_d_compact.h5"
  CHARACTER(LEN=3) , PARAMETER :: dataset  = "DS1"
  INTEGER          , PARAMETER :: dim0     = 4
  INTEGER          , PARAMETER :: dim1     = 7

  INTEGER :: hdferr
  INTEGER :: layout
  INTEGER(HID_T)  :: file, space, dset, dcpl ! Handles
  INTEGER(HSIZE_T), DIMENSION(1:2) :: dims = (/dim0, dim1/)

  INTEGER, DIMENSION(1:dim0, 1:dim1) :: wdata, & ! Write buffer
                                        rdata    ! Read buffer
  INTEGER :: i, j
  !
  ! Initialize FORTRAN interface.
  !
  CALL h5open_f(hdferr)
  !
  ! Initialize data.
  !
  DO i = 1, dim0
     DO j = 1, dim1
        wdata(i,j) = (i-1)*(j-1)-(j-1)
     ENDDO
  ENDDO
  !
  ! Create a new file using the default properties.
  !
  CALL h5fcreate_f(filename, H5F_ACC_TRUNC_F, file, hdferr)
  !
  ! Create dataspace.  Setting maximum size to NULL sets the maximum
  ! size to be the current size.
  !
  CALL h5screate_simple_f(2, dims, space, hdferr)
  !
  ! Create the dataset creation property list, set the layout to
  ! compact.
  !
  CALL h5pcreate_f(H5P_DATASET_CREATE_F, dcpl, hdferr)
  CALL h5pset_layout_f(dcpl, H5D_COMPACT_F, hdferr)
  !
  ! Create the dataset.  We will use all default properties for this
  ! example.
  !
  CALL h5dcreate_f(file, dataset, H5T_STD_I32LE, space, dset, hdferr, dcpl)
  !
  ! Write the data to the dataset.
  !
  CALL h5dwrite_f(dset, H5T_NATIVE_INTEGER, wdata, dims, hdferr)
  !
  ! Close and release resources.
  !
  CALL h5pclose_f(dcpl , hdferr)
  CALL h5dclose_f(dset , hdferr)
  CALL h5sclose_f(space, hdferr)
  CALL h5fclose_f(file , hdferr)
  !
  ! Now we begin the READ section of this example.
  !
  !
  ! Open file and dataset using the default properties.
  !
  CALL h5fopen_f(filename, H5F_ACC_RDONLY_F, file, hdferr)
  CALL h5dopen_f (file, dataset, dset, hdferr)
  !
  ! Retrieve the dataset creation property list, and print the
  ! storage layout.
  !
  CALL h5dget_create_plist_f(dset, dcpl, hdferr)
  CALL h5pget_layout_f(dcpl, layout, hdferr)
  WRITE(*,'(/,"Storage layout for ", A," is: ")', ADVANCE='NO') dataset
  IF(layout.EQ.H5D_COMPACT_F)THEN
     WRITE(*,'("H5D_COMPACT_F",/)')
  ELSE IF (layout.EQ.H5D_CONTIGUOUS_F)THEN
     WRITE(*,'("H5D_CONTIGUOUS_F",/)')
  ELSE IF (layout.EQ.H5D_CHUNKED_F)THEN
     WRITE(*,'("H5D_CHUNKED_F",/)')
#if H5_VERSION_GE(1,12,0)
  ELSE IF (layout.EQ.H5D_VIRTUAL_F)THEN
     WRITE(*,'("H5D_VIRTUAL_F",/)')
#endif
  ELSE
     WRITE(*,'("Layout Error",/)')
  ENDIF
  !
  ! Read the data using the default properties.
  !
  CALL h5dread_f(dset, H5T_NATIVE_INTEGER, rdata, dims, hdferr)
  !
  ! Output the data to the screen.
  !
  WRITE(*, '(A,":")') dataset
  DO i=1, dim0
     WRITE(*,'(" [")', ADVANCE='NO')
     WRITE(*,'(80i3)', ADVANCE='NO') rdata(i,:)
     WRITE(*,'(" ]")')
  ENDDO
  !
  ! Close and release resources.
  !
  CALL h5pclose_f(dcpl , hdferr)
  CALL h5dclose_f(dset , hdferr)
  CALL h5fclose_f(file , hdferr)

END PROGRAM main