1!     This is part of the netCDF package. Copyright 2006-2019
2!     University Corporation for Atmospheric Research/Unidata. See
3!     COPYRIGHT file for conditions of use.
4
5!     This program tests netCDF-4 variable functions from fortran.
6
7!     $Id: f90tst_vars2.f90,v 1.7 2010/01/25 21:01:07 ed Exp $
8
9program f90tst_vars2
10  use typeSizes
11  use netcdf
12  implicit none
13
14  ! This is the name of the data file we will create.
15  character (len = *), parameter :: FILE_NAME = "f90tst_vars2.nc"
16
17  ! We are writing 2D data, a 6 x 12 grid.
18  integer, parameter :: MAX_DIMS = 2
19  integer, parameter :: NX = 6, NY = 12
20  integer :: data_out(NY, NX), data_in(NY, NX)
21  integer :: data_out_1d(NX), data_in_1d(NX)
22
23  ! We need these ids and other gunk for netcdf.
24  integer :: ncid, varid1, varid2, varid3, varid4, varid5, dimids(MAX_DIMS)
25  integer :: chunksizes(MAX_DIMS), chunksizes_in(MAX_DIMS)
26  integer :: x_dimid, y_dimid
27  integer :: nvars, ngatts, ndims, unlimdimid, file_format
28  integer :: x, y
29  integer, parameter :: DEFLATE_LEVEL = 4
30  integer (kind = EightByteInt), parameter :: TOE_SAN_VALUE = 2147483648_EightByteInt
31  character (len = *), parameter :: VAR1_NAME = "Chon-Ji"
32  character (len = *), parameter :: VAR2_NAME = "Tan-Gun"
33  character (len = *), parameter :: VAR3_NAME = "Toe-San"
34  character (len = *), parameter :: VAR4_NAME = "Won-Hyo"
35  character (len = *), parameter :: VAR5_NAME = "Yul-Guk"
36  integer, parameter :: CACHE_SIZE = 8, CACHE_NELEMS = 571
37  integer, parameter :: CACHE_PREEMPTION = 66
38
39  ! Information read back from the file to check correctness.
40  integer :: varid1_in, varid2_in, varid3_in, varid4_in, varid5_in
41  integer :: xtype_in, ndims_in, natts_in, dimids_in(MAX_DIMS)
42  character (len = nf90_max_name) :: name_in
43  integer :: endianness_in, deflate_level_in
44  logical :: shuffle_in, fletcher32_in, contiguous_in
45  integer (kind = EightByteInt) :: toe_san_in
46  integer :: cache_size_in, cache_nelems_in, cache_preemption_in
47
48  print *, ''
49  print *,'*** Testing definition of netCDF-4 vars from Fortran 90.'
50
51  ! Create some pretend data.
52  do x = 1, NX
53     do y = 1, NY
54        data_out(y, x) = (x - 1) * NY + (y - 1)
55     end do
56  end do
57  do x = 1, NX
58        data_out_1d(x) = x
59  end do
60
61  ! Create the netCDF file.
62  call check(nf90_create(FILE_NAME, nf90_netcdf4, ncid, cache_nelems = CACHE_NELEMS, &
63       cache_size = CACHE_SIZE))
64
65  ! Define the dimensions.
66  call check(nf90_def_dim(ncid, "x", NX, x_dimid))
67  call check(nf90_def_dim(ncid, "y", NY, y_dimid))
68  dimids =  (/ y_dimid, x_dimid /)
69
70  ! Define some variables.
71  chunksizes = (/ NY, NX /)
72  call check(nf90_def_var(ncid, VAR1_NAME, NF90_INT, dimids, varid1, chunksizes = chunksizes, &
73       shuffle = .TRUE., fletcher32 = .TRUE., endianness = nf90_endian_big, deflate_level = DEFLATE_LEVEL))
74  call check(nf90_def_var(ncid, VAR2_NAME, NF90_INT, dimids, varid2, contiguous = .TRUE.))
75  call check(nf90_def_var(ncid, VAR3_NAME, NF90_INT64, varid3))
76  call check(nf90_def_var(ncid, VAR4_NAME, NF90_INT, x_dimid, varid4, contiguous = .TRUE.))
77  call check(nf90_def_var(ncid, VAR5_NAME, NF90_INT, dimids, varid5, chunksizes = chunksizes))
78
79  ! Write the pretend data to the file.
80  call check(nf90_put_var(ncid, varid1, data_out))
81  call check(nf90_put_var(ncid, varid2, data_out))
82  call check(nf90_put_var(ncid, varid3, TOE_SAN_VALUE))
83  call check(nf90_put_var(ncid, varid4, data_out_1d))
84  call check(nf90_put_var(ncid, varid5, data_out))
85
86  ! Close the file.
87  call check(nf90_close(ncid))
88
89  ! Reopen the file.
90  call check(nf90_open(FILE_NAME, nf90_nowrite, ncid))
91
92  ! Check some stuff out.
93  call check(nf90_inquire(ncid, ndims, nvars, ngatts, unlimdimid, file_format))
94  if (ndims /= 2 .or. nvars /= 5 .or. ngatts /= 0 .or. unlimdimid /= -1 .or. &
95       file_format /= nf90_format_netcdf4) stop 2
96
97  ! Get varids.
98  call check(nf90_inq_varid(ncid, VAR1_NAME, varid1_in))
99  call check(nf90_inq_varid(ncid, VAR2_NAME, varid2_in))
100  call check(nf90_inq_varid(ncid, VAR3_NAME, varid3_in))
101  call check(nf90_inq_varid(ncid, VAR4_NAME, varid4_in))
102  call check(nf90_inq_varid(ncid, VAR5_NAME, varid5_in))
103
104  ! Check variable 1.
105  call check(nf90_inquire_variable(ncid, varid1_in, name_in, xtype_in, ndims_in, dimids_in, &
106       natts_in, chunksizes = chunksizes_in, endianness = endianness_in, fletcher32 = fletcher32_in, &
107       deflate_level = deflate_level_in, shuffle = shuffle_in, cache_size = cache_size_in, &
108       cache_nelems = cache_nelems_in, cache_preemption = cache_preemption_in))
109  if (name_in .ne. VAR1_NAME .or. xtype_in .ne. NF90_INT .or. ndims_in .ne. MAX_DIMS .or. &
110       natts_in .ne. 0 .or. dimids_in(1) .ne. dimids(1) .or. dimids_in(2) .ne. dimids(2)) stop 3
111  if (chunksizes_in(1) /= chunksizes(1) .or. chunksizes_in(2) /= chunksizes(2)) &
112       stop 4
113  if (endianness_in .ne. nf90_endian_big) stop 5
114
115  ! This test code commented out because it fails parallel builds,
116  ! which don't use the cache, so don't get the same size
117  ! settings. Since we are not (yet) using a preprocessor on the
118  ! fortran code, there's no way to ifdef out these tests.
119  ! print *, cache_size_in, cache_nelems_in, cache_preemption
120  ! if (cache_size_in .ne. 16 .or. cache_nelems_in .ne. 4133 .or. &
121  !      cache_preemption .ne. CACHE_PREEMPTION) stop 555
122
123  ! Check variable 2.
124  call check(nf90_inquire_variable(ncid, varid2_in, name_in, xtype_in, ndims_in, dimids_in, &
125       natts_in, contiguous = contiguous_in, endianness = endianness_in, fletcher32 = fletcher32_in, &
126       deflate_level = deflate_level_in, shuffle = shuffle_in))
127  if (name_in .ne. VAR2_NAME .or. xtype_in .ne. NF90_INT .or. ndims_in .ne. MAX_DIMS .or. &
128       natts_in .ne. 0 .or. dimids_in(1) .ne. dimids(1) .or. dimids_in(2) .ne. dimids(2)) stop 6
129  if (deflate_level_in .ne. 0 .or. .not. contiguous_in .or. fletcher32_in .or. shuffle_in) stop 7
130
131  ! Check variable 3.
132  call check(nf90_inquire_variable(ncid, varid3_in, name_in, xtype_in, ndims_in, dimids_in, &
133       natts_in, contiguous = contiguous_in, endianness = endianness_in, fletcher32 = fletcher32_in, &
134       deflate_level = deflate_level_in, shuffle = shuffle_in))
135  if (name_in .ne. VAR3_NAME .or. xtype_in .ne. NF90_INT64 .or. ndims_in .ne. 0 .or. &
136       natts_in .ne. 0) stop 8
137  if (deflate_level_in .ne. 0 .or. .not. contiguous_in .or. fletcher32_in .or. shuffle_in) stop 9
138
139  ! Check variable 4.
140  call check(nf90_inquire_variable(ncid, varid4_in, name_in, xtype_in, ndims_in, dimids_in, &
141       natts_in, contiguous = contiguous_in, endianness = endianness_in, fletcher32 = fletcher32_in, &
142       deflate_level = deflate_level_in, shuffle = shuffle_in))
143  if (name_in .ne. VAR4_NAME .or. xtype_in .ne. NF90_INT .or. ndims_in .ne. 1 .or. &
144       natts_in .ne. 0 .or. dimids_in(1) .ne. x_dimid) stop 10
145  if (deflate_level_in .ne. 0 .or. .not. contiguous_in .or. fletcher32_in .or. shuffle_in) stop 11
146
147  ! Check the data.
148  call check(nf90_get_var(ncid, varid1_in, data_in))
149  do x = 1, NX
150     do y = 1, NY
151        if (data_out(y, x) .ne. data_in(y, x)) stop 12
152     end do
153  end do
154  call check(nf90_get_var(ncid, varid2_in, data_in))
155  do x = 1, NX
156     do y = 1, NY
157        if (data_out(y, x) .ne. data_in(y, x)) stop 13
158     end do
159  end do
160  call check(nf90_get_var(ncid, varid3_in, toe_san_in))
161  if (toe_san_in .ne. TOE_SAN_VALUE) stop 14
162  call check(nf90_get_var(ncid, varid4_in, data_in_1d))
163  do x = 1, NX
164     if (data_out_1d(x) .ne. data_in_1d(x)) stop 15
165  end do
166  call check(nf90_get_var(ncid, varid5_in, data_in))
167  do x = 1, NX
168     do y = 1, NY
169        if (data_out(y, x) .ne. data_in(y, x)) stop 12
170     end do
171  end do
172
173  ! Close the file.
174  call check(nf90_close(ncid))
175
176  print *,'*** SUCCESS!'
177
178contains
179!     This subroutine handles errors by printing an error message and
180!     exiting with a non-zero status.
181  subroutine check(errcode)
182    use netcdf
183    implicit none
184    integer, intent(in) :: errcode
185
186    if(errcode /= nf90_noerr) then
187       print *, 'Error: ', trim(nf90_strerror(errcode))
188       stop 2
189    endif
190  end subroutine check
191end program f90tst_vars2
192
193