[496c607] | 1 | MODULE fp2nc4io_mod |
---|
| 2 | |
---|
| 3 | !**************************************************************** |
---|
| 4 | ! * |
---|
| 5 | ! Contains data and routines for dumping selected FLEXPART * |
---|
| 6 | ! array variables to a NetCDF4 file. * |
---|
| 7 | ! * |
---|
| 8 | ! Don Morton (Boreal Scientific Computing LLC) * |
---|
| 9 | ! * |
---|
| 10 | ! May 2016 * |
---|
| 11 | ! * |
---|
| 12 | !**************************************************************** |
---|
| 13 | |
---|
| 14 | USE par_mod |
---|
| 15 | USE com_mod |
---|
| 16 | |
---|
| 17 | USE netcdf |
---|
| 18 | |
---|
| 19 | IMPLICIT NONE |
---|
| 20 | |
---|
| 21 | ! This variable should be in the range [1,9]. It has been suggested |
---|
| 22 | ! that 2 offers reasonable compression in a reasonable time. |
---|
| 23 | ! Higher values will offer more compression, but will take more time |
---|
| 24 | INTEGER, PARAMETER :: DEFAULT_DEFLATE_LEVEL = 2 |
---|
| 25 | PRIVATE DEFAULT_DEFLATE_LEVEL |
---|
| 26 | |
---|
| 27 | ! These are valid variable names for the user of this module to reference |
---|
[87d9684] | 28 | !!! DJM - 2016-06-13 -- added specific value in DIMENSION statement. |
---|
| 29 | !!! can't be "*" in some Fortran versions |
---|
| 30 | CHARACTER, DIMENSION(10), PARAMETER :: VALID_VARS = & |
---|
[496c607] | 31 | & (/ 't', 'u', 'v', 'w', 'q', & |
---|
| 32 | & 'T', 'U', 'V', 'W', 'Q' /) |
---|
| 33 | PRIVATE VALID_VARS |
---|
| 34 | |
---|
| 35 | ! Private routines in this module |
---|
| 36 | PRIVATE private_dump_3dfield |
---|
| 37 | PRIVATE private_read_3dfield |
---|
| 38 | PRIVATE to_upper |
---|
| 39 | |
---|
| 40 | |
---|
| 41 | CONTAINS |
---|
| 42 | |
---|
| 43 | SUBROUTINE fp2nc4io_print_valid_vars |
---|
| 44 | |
---|
| 45 | ! Prints the list of met variables that are considered valid in this |
---|
| 46 | ! module |
---|
| 47 | |
---|
| 48 | IMPLICIT NONE |
---|
| 49 | INTEGER :: i |
---|
| 50 | |
---|
| 51 | DO i=1,SIZE(VALID_VARS) |
---|
| 52 | PRINT *, VALID_VARS(i) |
---|
| 53 | ENDDO |
---|
| 54 | |
---|
| 55 | END SUBROUTINE fp2nc4io_print_valid_vars |
---|
| 56 | |
---|
| 57 | |
---|
| 58 | |
---|
| 59 | LOGICAL FUNCTION fp2nc4io_vars_are_valid(num_vars, dump_vars) |
---|
| 60 | |
---|
| 61 | ! Returns True or False depending on whether all of the variables |
---|
| 62 | ! in dump_vars are valid names (according to VALID_VARS) |
---|
| 63 | |
---|
| 64 | IMPLICIT NONE |
---|
| 65 | |
---|
| 66 | INTEGER, INTENT(IN) :: num_vars |
---|
| 67 | CHARACTER, DIMENSION(num_vars), INTENT(IN) :: dump_vars ! var list |
---|
| 68 | |
---|
| 69 | LOGICAL :: all_good = .TRUE. |
---|
| 70 | INTEGER :: i |
---|
| 71 | |
---|
| 72 | DO i=1,num_vars |
---|
| 73 | IF( .NOT. ANY(VALID_VARS == dump_vars(i)) ) THEN |
---|
| 74 | all_good = .FALSE. |
---|
| 75 | ENDIF |
---|
| 76 | ENDDO |
---|
| 77 | |
---|
| 78 | fp2nc4io_vars_are_valid = all_good |
---|
| 79 | |
---|
| 80 | END FUNCTION fp2nc4io_vars_are_valid |
---|
| 81 | |
---|
| 82 | |
---|
| 83 | |
---|
| 84 | SUBROUTINE fp2nc4io_dump(nc4_filepath, num_vars, dump_vars, deflate_level) |
---|
| 85 | |
---|
| 86 | ! Writes metadata plus variables in dump_vars to NetCDF4 file |
---|
| 87 | ! All of the dumped variables come from FLEXPART modules |
---|
| 88 | ! par_mod and com_mod |
---|
| 89 | |
---|
| 90 | IMPLICIT NONE |
---|
| 91 | |
---|
| 92 | CHARACTER(LEN=*), INTENT(IN) :: nc4_filepath ! Full path to dump file |
---|
| 93 | INTEGER, INTENT(IN) :: num_vars ! Num variables in dump_vars |
---|
| 94 | CHARACTER, DIMENSION(num_vars), INTENT(IN) :: dump_vars ! var list |
---|
| 95 | INTEGER, OPTIONAL, INTENT(IN) :: deflate_level ! (should be 0-9) |
---|
| 96 | |
---|
| 97 | |
---|
| 98 | INTEGER :: i, j, k |
---|
| 99 | INTEGER :: ncfunc_retval ! NetCDF function call return values |
---|
| 100 | INTEGER :: ncid ! NetCDF file id |
---|
| 101 | |
---|
| 102 | |
---|
| 103 | ! Variables used by NetCDF routines |
---|
| 104 | INTEGER :: x_dimid, y_dimid, z_dimid, dimids(3) |
---|
| 105 | INTEGER :: varid |
---|
| 106 | INTEGER :: deflevel ! Deflate level |
---|
| 107 | |
---|
| 108 | #ifdef TESTING |
---|
| 109 | INTEGER :: nx_test, ny_test, nz_test |
---|
| 110 | REAL, ALLOCATABLE, DIMENSION(:,:,:) :: testvar_array |
---|
| 111 | CHARACTER(LEN=NF90_MAX_NAME) :: x_dimname, y_dimname, z_dimname |
---|
| 112 | REAL :: orig_array_sum, test_array_sum |
---|
| 113 | #endif |
---|
| 114 | |
---|
| 115 | #ifdef TESTING |
---|
| 116 | PRINT *, |
---|
| 117 | PRINT *, '*** Running in testing mode ***' |
---|
| 118 | PRINT *, |
---|
| 119 | #endif |
---|
| 120 | |
---|
| 121 | ! Use default deflate level if it wasn't passed in, or if a bad |
---|
| 122 | ! value was passed in. |
---|
| 123 | IF (PRESENT(deflate_level)) THEN |
---|
| 124 | IF (deflate_level < 0 .OR. deflate_level > 9) THEN |
---|
| 125 | deflevel = DEFAULT_DEFLATE_LEVEL |
---|
| 126 | ELSE |
---|
| 127 | deflevel = deflate_level |
---|
| 128 | ENDIF |
---|
| 129 | ELSE |
---|
| 130 | deflevel = DEFAULT_DEFLATE_LEVEL |
---|
| 131 | ENDIF |
---|
| 132 | |
---|
| 133 | PRINT *, 'Using deflate level: ', deflevel |
---|
| 134 | |
---|
| 135 | !!!!!!--------------------------------------------------- |
---|
| 136 | !!!!!! Now we are ready to dump to NetCDF4 file |
---|
| 137 | !!!!!!--------------------------------------------------- |
---|
| 138 | |
---|
| 139 | ncfunc_retval = nf90_create(nc4_filepath, & |
---|
| 140 | & OR(NF90_CLOBBER, NF90_HDF5), ncid) |
---|
| 141 | PRINT *, 'Created file: ', TRIM(nc4_filepath) |
---|
| 142 | |
---|
| 143 | ! Define the dimensions, and get dimension ids passed back |
---|
| 144 | ! The values nx, ny and nz come from FP par_mod |
---|
| 145 | ncfunc_retval = nf90_def_dim(ncid, 'x', nx, x_dimid) |
---|
| 146 | ncfunc_retval = nf90_def_dim(ncid, 'y', ny, y_dimid) |
---|
| 147 | ncfunc_retval = nf90_def_dim(ncid, 'z', nz, z_dimid) |
---|
| 148 | dimids = (/ x_dimid, y_dimid, z_dimid /) |
---|
| 149 | |
---|
| 150 | ! Write each of the 3d variables to the NetCDF file |
---|
| 151 | DO i=1,num_vars |
---|
| 152 | CALL private_dump_3dfield(ncid, dump_vars(i), dimids, deflevel) |
---|
| 153 | PRINT *, 'Dumped 3d field: ', dump_vars(i) |
---|
| 154 | ENDDO |
---|
| 155 | |
---|
| 156 | ! Write the height field - variable 'height' is defined in com_mod |
---|
| 157 | ncfunc_retval = nf90_def_var(ncid, 'height', NF90_DOUBLE, & |
---|
| 158 | & z_dimid, varid) |
---|
[8624a75] | 159 | ! attributes |
---|
| 160 | ncfunc_retval = nf90_put_att(ncid, varid, "description","height of the FLEXPART model levels") |
---|
| 161 | ncfunc_retval = nf90_put_att(ncid, varid, "units","m a.g.l") |
---|
[496c607] | 162 | |
---|
| 163 | ncfunc_retval = nf90_def_var_deflate(ncid, varid, & |
---|
| 164 | & shuffle=0, & |
---|
| 165 | & deflate=1, & |
---|
| 166 | & deflate_level=deflevel) |
---|
| 167 | |
---|
| 168 | ncfunc_retval = nf90_put_var(ncid, varid, height(1:nz)) |
---|
| 169 | |
---|
| 170 | ! Write some of the scalar metadata variables |
---|
| 171 | ! dx, dy, xlon0, xlat0 are all defined in com_mod |
---|
| 172 | ncfunc_retval = nf90_def_var(ncid, 'dx', NF90_DOUBLE, varid) |
---|
| 173 | ncfunc_retval = nf90_put_var(ncid, varid, dx) |
---|
[8624a75] | 174 | ! attributes |
---|
| 175 | ncfunc_retval = nf90_put_att(ncid, varid, "description","grid distance in x direction") |
---|
| 176 | ncfunc_retval = nf90_put_att(ncid, varid, "units","degrees") |
---|
| 177 | |
---|
[496c607] | 178 | |
---|
| 179 | ncfunc_retval = nf90_def_var(ncid, 'dy', NF90_DOUBLE, varid) |
---|
| 180 | ncfunc_retval = nf90_put_var(ncid, varid, dy) |
---|
[8624a75] | 181 | ! attributes |
---|
| 182 | ncfunc_retval = nf90_put_att(ncid, varid, "description","grid distance in y direction") |
---|
| 183 | ncfunc_retval = nf90_put_att(ncid, varid, "units","degrees") |
---|
| 184 | |
---|
[496c607] | 185 | |
---|
| 186 | ncfunc_retval = nf90_def_var(ncid, 'xlon0', NF90_DOUBLE, varid) |
---|
| 187 | ncfunc_retval = nf90_put_var(ncid, varid, xlon0) |
---|
[8624a75] | 188 | ! attributes |
---|
| 189 | ncfunc_retval = nf90_put_att(ncid, varid, "description","longitude of the lowest left corner") |
---|
| 190 | ncfunc_retval = nf90_put_att(ncid, varid, "units","degrees") |
---|
| 191 | |
---|
[496c607] | 192 | |
---|
| 193 | ncfunc_retval = nf90_def_var(ncid, 'ylat0', NF90_DOUBLE, varid) |
---|
| 194 | ncfunc_retval = nf90_put_var(ncid, varid, ylat0) |
---|
[8624a75] | 195 | ! attributes |
---|
| 196 | ncfunc_retval = nf90_put_att(ncid, varid, "description","latitude of the lowest left corner") |
---|
| 197 | ncfunc_retval = nf90_put_att(ncid, varid, "units","degrees") |
---|
| 198 | |
---|
[496c607] | 199 | |
---|
| 200 | ! All done, close the NetCDF file |
---|
| 201 | ncfunc_retval = nf90_close(ncid) |
---|
| 202 | |
---|
| 203 | #ifdef TESTING |
---|
| 204 | !!!!!!!!!!!!!! Reading !!!!!!!!!!!!!!!!!!!! |
---|
| 205 | print *, "Opening nc file for reading" |
---|
| 206 | ncfunc_retval = nf90_open(nc4_filepath, NF90_NOWRITE, ncid) |
---|
| 207 | |
---|
| 208 | ! Get dimensions |
---|
| 209 | ncfunc_retval = nf90_inq_dimid(ncid, 'x', x_dimid) |
---|
| 210 | ncfunc_retval = nf90_inquire_dimension(ncid, x_dimid, x_dimname, & |
---|
| 211 | & nx_test) |
---|
| 212 | PRINT *, 'nx_test: ', nx_test |
---|
| 213 | |
---|
| 214 | ncfunc_retval = nf90_inq_dimid(ncid, 'y', y_dimid) |
---|
| 215 | ncfunc_retval = nf90_inquire_dimension(ncid, y_dimid, y_dimname, & |
---|
| 216 | & ny_test) |
---|
| 217 | PRINT *, 'ny_test: ', ny_test |
---|
| 218 | |
---|
| 219 | ncfunc_retval = nf90_inq_dimid(ncid, 'z', z_dimid) |
---|
| 220 | ncfunc_retval = nf90_inquire_dimension(ncid, z_dimid, z_dimname, & |
---|
| 221 | & nz_test) |
---|
| 222 | PRINT *, 'nz_test: ', nz_test |
---|
| 223 | |
---|
| 224 | ALLOCATE( testvar_array(0:nx_test-1, 0:ny_test-1, nz_test) ) |
---|
| 225 | |
---|
| 226 | ! Read each variable and compare with original data |
---|
| 227 | DO i=1,num_vars |
---|
| 228 | CALL private_read_3dfield(ncid, dump_vars(i), & |
---|
| 229 | & nx_test, ny_test, nz_test, & |
---|
| 230 | & testvar_array) |
---|
| 231 | |
---|
| 232 | |
---|
| 233 | IF (to_upper(dump_vars(i)) == 'U') THEN |
---|
| 234 | orig_array_sum = SUM( uu(0:nx_test-1, 0:ny_test-1, & |
---|
| 235 | & 1:nz_test, 1) ) |
---|
| 236 | ELSEIF (to_upper(dump_vars(i)) == 'V') THEN |
---|
| 237 | orig_array_sum = SUM( vv(0:nx_test-1, 0:ny_test-1, & |
---|
| 238 | & 1:nz_test, 1) ) |
---|
| 239 | ELSEIF (to_upper(dump_vars(i)) == 'T') THEN |
---|
| 240 | orig_array_sum = SUM( tt(0:nx_test-1, 0:ny_test-1, & |
---|
| 241 | & 1:nz_test, 1) ) |
---|
| 242 | ELSEIF (to_upper(dump_vars(i)) == 'W') THEN |
---|
| 243 | orig_array_sum = SUM( ww(0:nx_test-1, 0:ny_test-1, & |
---|
| 244 | & 1:nz_test, 1) ) |
---|
| 245 | ELSEIF (to_upper(dump_vars(i)) == 'Q') THEN |
---|
| 246 | orig_array_sum = SUM( qv(0:nx_test-1, 0:ny_test-1, & |
---|
| 247 | & 1:nz_test, 1) ) |
---|
| 248 | ENDIF |
---|
| 249 | |
---|
| 250 | test_array_sum = SUM( testvar_array(0:nx_test-1, 0:ny_test-1, & |
---|
| 251 | & 1:nz_test) ) |
---|
| 252 | |
---|
| 253 | PRINT *, dump_vars(i), ': ', 'SUM of differences = ', & |
---|
| 254 | & test_array_sum - orig_array_sum |
---|
| 255 | IF ( ABS(test_array_sum - orig_array_sum) .GT. 1.0E-3 ) THEN |
---|
| 256 | PRINT *, & |
---|
| 257 | & 'WARNING WILL ROBINSON!: Sum of differences exceeds 1.0E-3' |
---|
| 258 | ENDIF |
---|
| 259 | ENDDO |
---|
| 260 | |
---|
| 261 | ncfunc_retval = nf90_close(ncid) |
---|
| 262 | PRINT *, 'Closed file: ', ncfunc_retval |
---|
| 263 | |
---|
| 264 | #endif |
---|
| 265 | |
---|
| 266 | END SUBROUTINE fp2nc4io_dump |
---|
| 267 | |
---|
| 268 | |
---|
| 269 | SUBROUTINE private_dump_3dfield(ncid, varname, dimids, deflevel) |
---|
| 270 | |
---|
| 271 | ! Private routine meant to provide low level access for writing |
---|
| 272 | ! specified varname to NetCDF4 file. It is assumed that the |
---|
| 273 | ! NC4 file has already been opened, and that dimension id's have |
---|
| 274 | ! already been obtained |
---|
| 275 | |
---|
| 276 | IMPLICIT NONE |
---|
| 277 | |
---|
| 278 | INTEGER, INTENT(IN) :: ncid ! NC4 file id |
---|
| 279 | CHARACTER, INTENT(IN) :: varname |
---|
| 280 | INTEGER, INTENT(IN) :: dimids(3) ! NC4 dimension ids |
---|
| 281 | INTEGER, INTENT(IN) :: deflevel ! compression level |
---|
| 282 | |
---|
| 283 | ! NetCDF4 variables |
---|
| 284 | CHARACTER :: nc_varname |
---|
| 285 | INTEGER :: ncfunc_retval, varid |
---|
| 286 | |
---|
| 287 | ! Check that we have a valid varname. If not, buh-bye |
---|
| 288 | IF( .NOT. ANY(VALID_VARS == varname) ) THEN |
---|
| 289 | PRINT *, |
---|
| 290 | PRINT *, 'fp2nc4io:private_dump_3d_field() bad var: ', varname |
---|
| 291 | PRINT *, ' ABORTING...' |
---|
| 292 | PRINT *, |
---|
| 293 | STOP |
---|
| 294 | ENDIF |
---|
| 295 | |
---|
| 296 | ! Convert varname to upper case for use in NetCDF file |
---|
| 297 | nc_varname = to_upper(varname) |
---|
| 298 | |
---|
| 299 | ! Create the variable in the NetCDF file |
---|
| 300 | ncfunc_retval = nf90_def_var(ncid, nc_varname, NF90_DOUBLE, & |
---|
| 301 | & dimids, varid) |
---|
| 302 | |
---|
| 303 | ncfunc_retval = nf90_def_var_deflate(ncid, varid, & |
---|
| 304 | & shuffle=0, & |
---|
| 305 | & deflate=1, & |
---|
| 306 | & deflate_level=deflevel) |
---|
| 307 | |
---|
| 308 | ! Write the data arrays |
---|
| 309 | ! The values nx, ny and nz come from module com_mod |
---|
| 310 | ! Likewise, the arrays uu, vv, tt, ww, qv are also from the |
---|
| 311 | ! same module, and we assume they all have the same dimensions |
---|
| 312 | ! (currently they do) |
---|
| 313 | PRINT *, 'Writing: ', nc_varname |
---|
| 314 | IF (nc_varname == 'U') THEN |
---|
| 315 | ncfunc_retval = nf90_put_var(ncid, varid, & |
---|
| 316 | & uu(0:nx-1, 0:ny-1, 1:nz, 1)) |
---|
[8624a75] | 317 | ! attributes |
---|
| 318 | ncfunc_retval = nf90_put_att(ncid, varid, "description","U component of wind in the X[horizontal] direction") |
---|
| 319 | ncfunc_retval = nf90_put_att(ncid, varid, "units","m s**-1") |
---|
| 320 | |
---|
| 321 | |
---|
[496c607] | 322 | ELSEIF (nc_varname == 'V') THEN |
---|
| 323 | ncfunc_retval = nf90_put_var(ncid, varid, & |
---|
| 324 | & vv(0:nx-1, 0:ny-1, 1:nz, 1)) |
---|
[8624a75] | 325 | ! attributes |
---|
| 326 | ncfunc_retval = nf90_put_att(ncid, varid, "description","V component of wind in the Y[horizontal] direction") |
---|
| 327 | ncfunc_retval = nf90_put_att(ncid, varid, "units","m s**-1") |
---|
| 328 | |
---|
| 329 | |
---|
[496c607] | 330 | ELSEIF (nc_varname == 'T') THEN |
---|
| 331 | ncfunc_retval = nf90_put_var(ncid, varid, & |
---|
| 332 | & tt(0:nx-1, 0:ny-1, 1:nz, 1)) |
---|
[8624a75] | 333 | ! attributes |
---|
| 334 | ncfunc_retval = nf90_put_att(ncid, varid, "description","temperature") |
---|
| 335 | ncfunc_retval = nf90_put_att(ncid, varid, "units","k") |
---|
| 336 | |
---|
[496c607] | 337 | ELSEIF (nc_varname == 'W') THEN |
---|
| 338 | ncfunc_retval = nf90_put_var(ncid, varid, & |
---|
| 339 | & ww(0:nx-1, 0:ny-1, 1:nz, 1)) |
---|
[8624a75] | 340 | ! attributes |
---|
| 341 | ncfunc_retval = nf90_put_att(ncid, varid, "description","wind component in the Z[vertical] direction") |
---|
| 342 | ncfunc_retval = nf90_put_att(ncid, varid, "units","m s**-1") |
---|
| 343 | |
---|
| 344 | |
---|
[496c607] | 345 | ELSEIF (nc_varname == 'Q') THEN |
---|
| 346 | ncfunc_retval = nf90_put_var(ncid, varid, & |
---|
| 347 | & qv(0:nx-1, 0:ny-1, 1:nz, 1)) |
---|
[8624a75] | 348 | ! attributes |
---|
| 349 | ncfunc_retval = nf90_put_att(ncid, varid, "description","specific humidity") |
---|
| 350 | ncfunc_retval = nf90_put_att(ncid, varid, "units"," ") |
---|
| 351 | |
---|
| 352 | |
---|
[496c607] | 353 | ELSE |
---|
| 354 | PRINT *, |
---|
| 355 | PRINT *, 'fp2nc4io:private_dump_3d_field() bad var: ', nc_varname |
---|
| 356 | PRINT *, ' ABORTING...' |
---|
| 357 | PRINT *, |
---|
| 358 | ENDIF |
---|
| 359 | |
---|
| 360 | IF (ncfunc_retval /= 0) THEN |
---|
| 361 | PRINT *, |
---|
| 362 | PRINT *, '*** WARNING ***' |
---|
| 363 | PRINT *, ' fp2nc4io:private_dump_3d_field()' |
---|
| 364 | PRINT *, ' nf90_put_var returned error for var: ', nc_varname |
---|
| 365 | PRINT *, |
---|
| 366 | |
---|
| 367 | ENDIF |
---|
| 368 | |
---|
| 369 | |
---|
| 370 | END SUBROUTINE private_dump_3dfield |
---|
| 371 | |
---|
| 372 | |
---|
| 373 | |
---|
| 374 | SUBROUTINE private_read_3dfield(ncid, varname, xdim, ydim, zdim, var_array) |
---|
| 375 | |
---|
| 376 | ! Private routine for reading full 3D array, specified by varname, |
---|
| 377 | ! from NC4 file. Reads into preallocated array of size |
---|
| 378 | ! xdim x ydim x zdim |
---|
| 379 | IMPLICIT NONE |
---|
| 380 | |
---|
| 381 | INTEGER, INTENT(IN) :: ncid ! NC4 file id |
---|
| 382 | CHARACTER, INTENT(IN) :: varname |
---|
| 383 | INTEGER, INTENT(IN) :: xdim, ydim, zdim ! NC4 dimension ids |
---|
| 384 | REAL, DIMENSION(xdim, ydim, zdim) :: var_array |
---|
| 385 | |
---|
| 386 | CHARACTER :: nc_varname |
---|
| 387 | INTEGER :: ncfunc_retval, varid |
---|
| 388 | |
---|
| 389 | ! Check that we have a valid varname. If not, buh-bye |
---|
| 390 | IF( .NOT. ANY(VALID_VARS == varname) ) THEN |
---|
| 391 | PRINT *, |
---|
| 392 | PRINT *, 'fp2nc4io:private_dump_3d_field() bad var: ', varname |
---|
| 393 | PRINT *, ' ABORTING...' |
---|
| 394 | PRINT *, |
---|
| 395 | STOP |
---|
| 396 | ENDIF |
---|
| 397 | |
---|
| 398 | ! Convert varname to upper case for use in NetCDF file |
---|
| 399 | nc_varname = to_upper(varname) |
---|
| 400 | |
---|
| 401 | ! Get the varid |
---|
| 402 | ncfunc_retval = nf90_inq_varid(ncid, nc_varname, varid) |
---|
| 403 | |
---|
| 404 | ! Read the variable into var_array |
---|
| 405 | ncfunc_retval = nf90_get_var(ncid, varid, var_array) |
---|
| 406 | |
---|
| 407 | END SUBROUTINE private_read_3dfield |
---|
| 408 | |
---|
| 409 | |
---|
| 410 | CHARACTER FUNCTION to_upper(c) |
---|
| 411 | |
---|
| 412 | ! Utility function to convert lower case char to upper case |
---|
| 413 | |
---|
| 414 | IMPLICIT NONE |
---|
| 415 | |
---|
| 416 | CHARACTER, INTENT(IN) :: c |
---|
| 417 | |
---|
| 418 | INTEGER :: c_ascii_code |
---|
| 419 | |
---|
| 420 | c_ascii_code = IACHAR(c) |
---|
| 421 | IF (c_ascii_code >= IACHAR("a") .AND. c_ascii_code <= IACHAR("z")) THEN |
---|
| 422 | to_upper = ACHAR(c_ascii_code - 32) |
---|
| 423 | ELSE |
---|
| 424 | to_upper = c |
---|
| 425 | ENDIF |
---|
| 426 | |
---|
| 427 | END FUNCTION to_upper |
---|
| 428 | |
---|
| 429 | |
---|
| 430 | |
---|
| 431 | END MODULE fp2nc4io_mod |
---|