| 1 | MODULE io_netcdf |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! io_netcdf |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Provides input/output procedures and parameters for the PEM. |
|---|
| 8 | ! |
|---|
| 9 | ! AUTHORS & DATE |
|---|
| 10 | ! JB Clement, 12/2025 |
|---|
| 11 | ! |
|---|
| 12 | ! NOTES |
|---|
| 13 | ! |
|---|
| 14 | !----------------------------------------------------------------------- |
|---|
| 15 | |
|---|
| 16 | ! DEPEDENCIES |
|---|
| 17 | ! ----------- |
|---|
| 18 | use numerics, only: dp, di, k4, minieps |
|---|
| 19 | use netcdf, only: nf90_double, nf90_noerr, nf90_strerror, nf90_write, nf90_nowrite, & |
|---|
| 20 | nf90_open, nf90_close, nf90_redef, nf90_enddef, nf90_inquire, nf90_max_var_dims, & |
|---|
| 21 | nf90_inq_dimid, nf90_inquire_dimension, nf90_inq_varid, nf90_inquire_variable, & |
|---|
| 22 | nf90_def_var, nf90_get_var, nf90_put_var, nf90_get_att, nf90_put_att |
|---|
| 23 | use stoppage, only: stop_clean |
|---|
| 24 | |
|---|
| 25 | ! DECLARATION |
|---|
| 26 | ! ----------- |
|---|
| 27 | implicit none |
|---|
| 28 | |
|---|
| 29 | ! PARAMETERS |
|---|
| 30 | ! ---------- |
|---|
| 31 | character(8), parameter :: start_name = 'start.nc' |
|---|
| 32 | character(11), parameter :: start1D_name = 'start1D.txt' |
|---|
| 33 | character(10), parameter :: startfi_name = 'startfi.nc' |
|---|
| 34 | character(12), parameter :: startpem_name = 'startevol.nc' |
|---|
| 35 | character(19), parameter :: xios_day_name1 = 'Xoutdaily4pem_Y1.nc' ! XIOS daily output file, second to last PCM year |
|---|
| 36 | character(19), parameter :: xios_day_name2 = 'Xoutdaily4pem_Y2.nc' ! XIOS daily output file, last PCM year |
|---|
| 37 | character(20), parameter :: xios_yr_name1 = 'Xoutyearly4pem_Y1.nc' ! XIOS yearly output file, second to last PCM year |
|---|
| 38 | character(20), parameter :: xios_yr_name2 = 'Xoutyearly4pem_Y2.nc' ! XIOS yearly output file, last PCM year |
|---|
| 39 | character(11), parameter :: diagevol_name = 'diagevol.nc' |
|---|
| 40 | |
|---|
| 41 | ! VARIABLES |
|---|
| 42 | ! --------- |
|---|
| 43 | logical(k4), private :: open_ncfile = .false. ! Flag to true if a NetCDF file is already open |
|---|
| 44 | logical(k4), private :: open_diagevol = .false. ! Flag to true if "diagevol.nc" is already open |
|---|
| 45 | integer(di), private :: ncid_file ! File ID |
|---|
| 46 | integer(di), private :: varid ! Variable ID |
|---|
| 47 | integer(di) :: ncid_diagevol ! File ID specific to "diagevol.nc" |
|---|
| 48 | |
|---|
| 49 | ! INTERFACES |
|---|
| 50 | ! ---------- |
|---|
| 51 | interface get_var_nc |
|---|
| 52 | module procedure get_var0d_nc, get_var1d_nc, get_var2d_nc, get_var3d_nc, get_var4d_nc |
|---|
| 53 | end interface get_var_nc |
|---|
| 54 | |
|---|
| 55 | interface put_var_nc |
|---|
| 56 | module procedure put_var0d_nc, put_var1d_nc, put_var2d_nc, put_var3d_nc, put_var4d_nc |
|---|
| 57 | end interface put_var_nc |
|---|
| 58 | |
|---|
| 59 | interface check_valid_var_nc |
|---|
| 60 | module procedure check_valid_var0d_nc, check_valid_var1d_nc, check_valid_var2d_nc, check_valid_var3d_nc, check_valid_var4d_nc |
|---|
| 61 | end interface check_valid_var_nc |
|---|
| 62 | |
|---|
| 63 | contains |
|---|
| 64 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 65 | |
|---|
| 66 | !======================================================================= |
|---|
| 67 | SUBROUTINE check_nc(ierr,msg,found) |
|---|
| 68 | !----------------------------------------------------------------------- |
|---|
| 69 | ! NAME |
|---|
| 70 | ! check_nc |
|---|
| 71 | ! |
|---|
| 72 | ! DESCRIPTION |
|---|
| 73 | ! NetCDF error handler. |
|---|
| 74 | ! |
|---|
| 75 | ! AUTHORS & DATE |
|---|
| 76 | ! JB Clement, 12/2025 |
|---|
| 77 | ! |
|---|
| 78 | ! NOTES |
|---|
| 79 | ! |
|---|
| 80 | !----------------------------------------------------------------------- |
|---|
| 81 | |
|---|
| 82 | ! DEPENDENCIES |
|---|
| 83 | ! ------------ |
|---|
| 84 | use display, only: print_err |
|---|
| 85 | |
|---|
| 86 | ! ARGUMENTS |
|---|
| 87 | ! --------- |
|---|
| 88 | integer(di), intent(in) :: ierr |
|---|
| 89 | character(*), intent(in) :: msg |
|---|
| 90 | logical(k4), intent(out), optional :: found |
|---|
| 91 | |
|---|
| 92 | ! LOCAL VARIABLES |
|---|
| 93 | ! --------------- |
|---|
| 94 | logical(k4) :: tmp_found |
|---|
| 95 | |
|---|
| 96 | ! CODE |
|---|
| 97 | ! ---- |
|---|
| 98 | if (ierr == nf90_noerr) then |
|---|
| 99 | tmp_found = .true. |
|---|
| 100 | else |
|---|
| 101 | tmp_found = .false. |
|---|
| 102 | end if |
|---|
| 103 | |
|---|
| 104 | if (present(found)) then |
|---|
| 105 | found = tmp_found |
|---|
| 106 | else |
|---|
| 107 | if (.not. tmp_found) then |
|---|
| 108 | call print_err(trim(nf90_strerror(ierr))) |
|---|
| 109 | call stop_clean(__FILE__,__LINE__,'NetCDF error when '//trim(msg),ierr) |
|---|
| 110 | end if |
|---|
| 111 | end if |
|---|
| 112 | |
|---|
| 113 | END SUBROUTINE check_nc |
|---|
| 114 | !======================================================================= |
|---|
| 115 | |
|---|
| 116 | !======================================================================= |
|---|
| 117 | SUBROUTINE open_nc(filename,mode,itime) |
|---|
| 118 | !----------------------------------------------------------------------- |
|---|
| 119 | ! NAME |
|---|
| 120 | ! open_nc |
|---|
| 121 | ! |
|---|
| 122 | ! DESCRIPTION |
|---|
| 123 | ! Open a netCDF file for reading. |
|---|
| 124 | |
|---|
| 125 | ! AUTHORS & DATE |
|---|
| 126 | ! JB Clement, 12/2025 |
|---|
| 127 | |
|---|
| 128 | ! NOTES |
|---|
| 129 | ! |
|---|
| 130 | !----------------------------------------------------------------------- |
|---|
| 131 | |
|---|
| 132 | ! DECLARATION |
|---|
| 133 | ! ----------- |
|---|
| 134 | implicit none |
|---|
| 135 | |
|---|
| 136 | ! ARGUMENTS |
|---|
| 137 | ! --------- |
|---|
| 138 | character(*), intent(in) :: filename, mode |
|---|
| 139 | integer(di), intent(out), optional :: itime |
|---|
| 140 | |
|---|
| 141 | ! CODE |
|---|
| 142 | ! ---- |
|---|
| 143 | ! Diagevol logic |
|---|
| 144 | if (adjustl(trim(filename)) == diagevol_name) then |
|---|
| 145 | if (adjustl(trim(mode)) == 'read') then |
|---|
| 146 | call stop_clean(__FILE__,__LINE__,'opening mode "read" cannot be used with '//trim(filename)//'"!',1) |
|---|
| 147 | else if (adjustl(trim(mode)) == 'write') then |
|---|
| 148 | if (.not. open_diagevol) then |
|---|
| 149 | call check_nc(nf90_open(trim(filename),nf90_write,ncid_diagevol),'opening file '//trim(filename)//' to write') |
|---|
| 150 | if (present(itime)) call get_next_itime_nc('Time',itime) ! Next time index |
|---|
| 151 | open_diagevol = .true. |
|---|
| 152 | end if |
|---|
| 153 | return |
|---|
| 154 | else |
|---|
| 155 | call stop_clean(__FILE__,__LINE__,'opening mode "'//mode//'" unknown!',1) |
|---|
| 156 | end if |
|---|
| 157 | end if |
|---|
| 158 | |
|---|
| 159 | ! Standard logic |
|---|
| 160 | if (open_ncfile) then ! A file is already opened |
|---|
| 161 | call stop_clean(__FILE__,__LINE__,'a NetCDF file is already opened!',1) |
|---|
| 162 | else |
|---|
| 163 | if (adjustl(trim(mode)) == 'read') then |
|---|
| 164 | call check_nc(nf90_open(trim(filename),nf90_nowrite,ncid_file),'opening file '//trim(filename)//' to read') |
|---|
| 165 | else if (adjustl(trim(mode)) == 'write') then |
|---|
| 166 | call check_nc(nf90_open(trim(filename),nf90_write,ncid_file),'opening file '//trim(filename)//' to write') |
|---|
| 167 | if (present(itime)) call get_next_itime_nc('Time',itime) ! Next time index |
|---|
| 168 | else |
|---|
| 169 | call stop_clean(__FILE__,__LINE__,'opening mode "'//mode//'" unknown!',1) |
|---|
| 170 | end if |
|---|
| 171 | open_ncfile = .true. |
|---|
| 172 | end if |
|---|
| 173 | |
|---|
| 174 | END SUBROUTINE open_nc |
|---|
| 175 | !======================================================================= |
|---|
| 176 | |
|---|
| 177 | !======================================================================= |
|---|
| 178 | SUBROUTINE close_nc(filename) |
|---|
| 179 | !----------------------------------------------------------------------- |
|---|
| 180 | ! NAME |
|---|
| 181 | ! close_nc |
|---|
| 182 | ! |
|---|
| 183 | ! DESCRIPTION |
|---|
| 184 | ! Open a netCDF file. |
|---|
| 185 | |
|---|
| 186 | ! AUTHORS & DATE |
|---|
| 187 | ! JB Clement, 12/2025 |
|---|
| 188 | |
|---|
| 189 | ! NOTES |
|---|
| 190 | ! |
|---|
| 191 | !----------------------------------------------------------------------- |
|---|
| 192 | |
|---|
| 193 | ! DECLARATION |
|---|
| 194 | ! ----------- |
|---|
| 195 | implicit none |
|---|
| 196 | |
|---|
| 197 | ! ARGUMENTS |
|---|
| 198 | ! --------- |
|---|
| 199 | character(*), intent(in) :: filename |
|---|
| 200 | |
|---|
| 201 | ! CODE |
|---|
| 202 | ! ---- |
|---|
| 203 | if (adjustl(trim(filename)) == diagevol_name) then ! Diagevol logic |
|---|
| 204 | call check_nc(nf90_close(ncid_diagevol),'closing file '//trim(filename)) |
|---|
| 205 | open_diagevol = .false. |
|---|
| 206 | else ! Standard logic |
|---|
| 207 | call check_nc(nf90_close(ncid_file),'closing file '//trim(filename)) |
|---|
| 208 | open_ncfile = .false. |
|---|
| 209 | end if |
|---|
| 210 | |
|---|
| 211 | END SUBROUTINE close_nc |
|---|
| 212 | !======================================================================= |
|---|
| 213 | |
|---|
| 214 | !======================================================================= |
|---|
| 215 | SUBROUTINE get_dim_nc(dim_name,d,found) |
|---|
| 216 | !----------------------------------------------------------------------- |
|---|
| 217 | ! NAME |
|---|
| 218 | ! get_dim_nc |
|---|
| 219 | ! |
|---|
| 220 | ! DESCRIPTION |
|---|
| 221 | ! Read a 0D variable from open NetCDF file. |
|---|
| 222 | ! |
|---|
| 223 | ! AUTHORS & DATE |
|---|
| 224 | ! JB Clement, 12/2025 |
|---|
| 225 | ! |
|---|
| 226 | ! NOTES |
|---|
| 227 | ! |
|---|
| 228 | !----------------------------------------------------------------------- |
|---|
| 229 | |
|---|
| 230 | ! DECLARATION |
|---|
| 231 | ! ----------- |
|---|
| 232 | implicit none |
|---|
| 233 | |
|---|
| 234 | ! ARGUMENTS |
|---|
| 235 | ! --------- |
|---|
| 236 | character(*), intent(in) :: dim_name ! Variable name |
|---|
| 237 | integer(di), intent(out) :: d ! Output dimension |
|---|
| 238 | logical(k4), intent(out), optional :: found |
|---|
| 239 | |
|---|
| 240 | ! LOCAL VARIABLES |
|---|
| 241 | ! --------------- |
|---|
| 242 | integer(di) :: dimid ! Dimension ID |
|---|
| 243 | |
|---|
| 244 | ! CODE |
|---|
| 245 | ! ---- |
|---|
| 246 | if (present(found)) then |
|---|
| 247 | call check_nc(nf90_inq_dimid(ncid_file,dim_name,dimid),'inquiring '//dim_name,found) |
|---|
| 248 | call check_nc(nf90_inquire_dimension(ncid_file,dimid,len = d),'getting '//dim_name,found) |
|---|
| 249 | else |
|---|
| 250 | call check_nc(nf90_inq_dimid(ncid_file,dim_name,dimid),'inquiring '//dim_name) |
|---|
| 251 | call check_nc(nf90_inquire_dimension(ncid_file,dimid,len = d),'getting '//dim_name) |
|---|
| 252 | end if |
|---|
| 253 | |
|---|
| 254 | END SUBROUTINE get_dim_nc |
|---|
| 255 | !======================================================================= |
|---|
| 256 | |
|---|
| 257 | !======================================================================= |
|---|
| 258 | SUBROUTINE def_var_nc(var_name,title,units,dimids) |
|---|
| 259 | !----------------------------------------------------------------------- |
|---|
| 260 | ! NAME |
|---|
| 261 | ! def_var_nc |
|---|
| 262 | ! |
|---|
| 263 | ! DESCRIPTION |
|---|
| 264 | ! Define a variable into open NetCDF file. |
|---|
| 265 | ! |
|---|
| 266 | ! AUTHORS & DATE |
|---|
| 267 | ! JB Clement, 01/2026 |
|---|
| 268 | ! |
|---|
| 269 | ! NOTES |
|---|
| 270 | ! |
|---|
| 271 | !----------------------------------------------------------------------- |
|---|
| 272 | |
|---|
| 273 | ! DECLARATION |
|---|
| 274 | ! ----------- |
|---|
| 275 | implicit none |
|---|
| 276 | |
|---|
| 277 | ! ARGUMENTS |
|---|
| 278 | ! --------- |
|---|
| 279 | character(*), intent(in) :: var_name, title, units ! Variable|title|units name |
|---|
| 280 | integer(di), dimension(:), intent(in) :: dimids ! Dimensions IDs |
|---|
| 281 | |
|---|
| 282 | ! LOCAL VARIABLES |
|---|
| 283 | ! --------------- |
|---|
| 284 | logical(k4) :: found |
|---|
| 285 | integer(di) :: ncid |
|---|
| 286 | |
|---|
| 287 | ! CODE |
|---|
| 288 | ! ---- |
|---|
| 289 | ! Diagevol logic priority over standard logic |
|---|
| 290 | if (open_diagevol) then |
|---|
| 291 | ncid = ncid_diagevol |
|---|
| 292 | else |
|---|
| 293 | ncid = ncid_file |
|---|
| 294 | end if |
|---|
| 295 | |
|---|
| 296 | ! Check if the variable exists |
|---|
| 297 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name,found) |
|---|
| 298 | if (found) return ! Variable is already defined |
|---|
| 299 | |
|---|
| 300 | ! Enter define mode |
|---|
| 301 | call check_nc(nf90_redef(ncid),'entering define mode') |
|---|
| 302 | |
|---|
| 303 | ! Define variable |
|---|
| 304 | call check_nc(nf90_def_var(ncid,var_name,NF90_DOUBLE,dimids,varid),'defining variable '//var_name) |
|---|
| 305 | call check_nc(nf90_put_att(ncid,varid,'title',title),'putting title attribute for '//var_name) |
|---|
| 306 | call check_nc(nf90_put_att(ncid,varid,'units',units),'putting units attribute for '//var_name) |
|---|
| 307 | |
|---|
| 308 | ! Leave define mode |
|---|
| 309 | call check_nc(nf90_enddef(ncid),'leaving define mode') |
|---|
| 310 | |
|---|
| 311 | END SUBROUTINE def_var_nc |
|---|
| 312 | !======================================================================= |
|---|
| 313 | |
|---|
| 314 | !======================================================================= |
|---|
| 315 | SUBROUTINE get_var0d_nc(var_name,var,found) |
|---|
| 316 | !----------------------------------------------------------------------- |
|---|
| 317 | ! NAME |
|---|
| 318 | ! get_var0d_nc |
|---|
| 319 | ! |
|---|
| 320 | ! DESCRIPTION |
|---|
| 321 | ! Read a 0D variable from open NetCDF file. |
|---|
| 322 | ! |
|---|
| 323 | ! AUTHORS & DATE |
|---|
| 324 | ! JB Clement, 12/2025 |
|---|
| 325 | ! |
|---|
| 326 | ! NOTES |
|---|
| 327 | ! |
|---|
| 328 | !----------------------------------------------------------------------- |
|---|
| 329 | |
|---|
| 330 | ! DECLARATION |
|---|
| 331 | ! ----------- |
|---|
| 332 | implicit none |
|---|
| 333 | |
|---|
| 334 | ! ARGUMENTS |
|---|
| 335 | ! --------- |
|---|
| 336 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 337 | real(dp), intent(out) :: var ! Output variable |
|---|
| 338 | logical(k4), intent(out), optional :: found |
|---|
| 339 | |
|---|
| 340 | ! CODE |
|---|
| 341 | ! ---- |
|---|
| 342 | if (present(found)) then |
|---|
| 343 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) |
|---|
| 344 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found) |
|---|
| 345 | else |
|---|
| 346 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name) |
|---|
| 347 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) |
|---|
| 348 | end if |
|---|
| 349 | call check_valid_var_nc(var_name,var) |
|---|
| 350 | |
|---|
| 351 | END SUBROUTINE get_var0d_nc |
|---|
| 352 | !======================================================================= |
|---|
| 353 | |
|---|
| 354 | !======================================================================= |
|---|
| 355 | SUBROUTINE get_var1d_nc(var_name,var,found) |
|---|
| 356 | !----------------------------------------------------------------------- |
|---|
| 357 | ! NAME |
|---|
| 358 | ! get_var1d_nc |
|---|
| 359 | ! |
|---|
| 360 | ! DESCRIPTION |
|---|
| 361 | ! Read a 1D variable from open NetCDF file. |
|---|
| 362 | ! |
|---|
| 363 | ! AUTHORS & DATE |
|---|
| 364 | ! JB Clement, 12/2025 |
|---|
| 365 | ! |
|---|
| 366 | ! NOTES |
|---|
| 367 | ! |
|---|
| 368 | !----------------------------------------------------------------------- |
|---|
| 369 | |
|---|
| 370 | ! DECLARATION |
|---|
| 371 | ! ----------- |
|---|
| 372 | implicit none |
|---|
| 373 | |
|---|
| 374 | ! ARGUMENTS |
|---|
| 375 | ! --------- |
|---|
| 376 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 377 | real(dp), dimension(:), intent(out) :: var ! Output variable |
|---|
| 378 | logical(k4), intent(out), optional :: found |
|---|
| 379 | |
|---|
| 380 | ! CODE |
|---|
| 381 | ! ---- |
|---|
| 382 | if (present(found)) then |
|---|
| 383 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) |
|---|
| 384 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found) |
|---|
| 385 | else |
|---|
| 386 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name) |
|---|
| 387 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) |
|---|
| 388 | end if |
|---|
| 389 | call check_valid_var_nc(var_name,var) |
|---|
| 390 | |
|---|
| 391 | END SUBROUTINE get_var1d_nc |
|---|
| 392 | !======================================================================= |
|---|
| 393 | |
|---|
| 394 | !======================================================================= |
|---|
| 395 | SUBROUTINE get_var2d_nc(var_name,var,found) |
|---|
| 396 | !----------------------------------------------------------------------- |
|---|
| 397 | ! NAME |
|---|
| 398 | ! get_var2d_nc |
|---|
| 399 | ! |
|---|
| 400 | ! DESCRIPTION |
|---|
| 401 | ! Read a 2D variable from open NetCDF file. |
|---|
| 402 | ! |
|---|
| 403 | ! AUTHORS & DATE |
|---|
| 404 | ! JB Clement, 12/2025 |
|---|
| 405 | ! |
|---|
| 406 | ! NOTES |
|---|
| 407 | ! |
|---|
| 408 | !----------------------------------------------------------------------- |
|---|
| 409 | |
|---|
| 410 | ! DECLARATION |
|---|
| 411 | ! ----------- |
|---|
| 412 | implicit none |
|---|
| 413 | |
|---|
| 414 | ! ARGUMENTS |
|---|
| 415 | ! --------- |
|---|
| 416 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 417 | real(dp), dimension(:,:), intent(out) :: var ! Output variable |
|---|
| 418 | logical(k4), intent(out), optional :: found |
|---|
| 419 | |
|---|
| 420 | ! CODE |
|---|
| 421 | ! ---- |
|---|
| 422 | if (present(found)) then |
|---|
| 423 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) |
|---|
| 424 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found) |
|---|
| 425 | else |
|---|
| 426 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name) |
|---|
| 427 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) |
|---|
| 428 | end if |
|---|
| 429 | call check_valid_var_nc(var_name,var) |
|---|
| 430 | |
|---|
| 431 | END SUBROUTINE get_var2d_nc |
|---|
| 432 | !======================================================================= |
|---|
| 433 | |
|---|
| 434 | !======================================================================= |
|---|
| 435 | SUBROUTINE get_var3d_nc(var_name,var,found) |
|---|
| 436 | !----------------------------------------------------------------------- |
|---|
| 437 | ! NAME |
|---|
| 438 | ! get_var3d_nc |
|---|
| 439 | ! |
|---|
| 440 | ! DESCRIPTION |
|---|
| 441 | ! Read a 3D variable from open NetCDF file. |
|---|
| 442 | ! |
|---|
| 443 | ! AUTHORS & DATE |
|---|
| 444 | ! JB Clement, 12/2025 |
|---|
| 445 | ! |
|---|
| 446 | ! NOTES |
|---|
| 447 | ! |
|---|
| 448 | !----------------------------------------------------------------------- |
|---|
| 449 | |
|---|
| 450 | ! DECLARATION |
|---|
| 451 | ! ----------- |
|---|
| 452 | implicit none |
|---|
| 453 | |
|---|
| 454 | ! ARGUMENTS |
|---|
| 455 | ! --------- |
|---|
| 456 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 457 | real(dp), dimension(:,:,:), intent(out) :: var ! Output variable |
|---|
| 458 | logical(k4), intent(out), optional :: found |
|---|
| 459 | |
|---|
| 460 | ! CODE |
|---|
| 461 | ! ---- |
|---|
| 462 | if (present(found)) then |
|---|
| 463 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) |
|---|
| 464 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found) |
|---|
| 465 | else |
|---|
| 466 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name) |
|---|
| 467 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) |
|---|
| 468 | end if |
|---|
| 469 | call check_valid_var_nc(var_name,var) |
|---|
| 470 | |
|---|
| 471 | END SUBROUTINE get_var3d_nc |
|---|
| 472 | !======================================================================= |
|---|
| 473 | |
|---|
| 474 | !======================================================================= |
|---|
| 475 | SUBROUTINE get_var4d_nc(var_name,var,found) |
|---|
| 476 | !----------------------------------------------------------------------- |
|---|
| 477 | ! NAME |
|---|
| 478 | ! get_var4d_nc |
|---|
| 479 | ! |
|---|
| 480 | ! DESCRIPTION |
|---|
| 481 | ! Read a 4D variable from open NetCDF file. |
|---|
| 482 | ! |
|---|
| 483 | ! AUTHORS & DATE |
|---|
| 484 | ! JB Clement, 12/2025 |
|---|
| 485 | ! |
|---|
| 486 | ! NOTES |
|---|
| 487 | ! |
|---|
| 488 | !----------------------------------------------------------------------- |
|---|
| 489 | |
|---|
| 490 | ! DECLARATION |
|---|
| 491 | ! ----------- |
|---|
| 492 | implicit none |
|---|
| 493 | |
|---|
| 494 | ! ARGUMENTS |
|---|
| 495 | ! --------- |
|---|
| 496 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 497 | real(dp), dimension(:,:,:,:), intent(out) :: var ! Output variable |
|---|
| 498 | logical(k4), intent(out), optional :: found |
|---|
| 499 | |
|---|
| 500 | ! CODE |
|---|
| 501 | ! ---- |
|---|
| 502 | if (present(found)) then |
|---|
| 503 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) |
|---|
| 504 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found) |
|---|
| 505 | else |
|---|
| 506 | call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name) |
|---|
| 507 | call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) |
|---|
| 508 | end if |
|---|
| 509 | call check_valid_var_nc(var_name,var) |
|---|
| 510 | |
|---|
| 511 | END SUBROUTINE get_var4d_nc |
|---|
| 512 | !======================================================================= |
|---|
| 513 | |
|---|
| 514 | !======================================================================= |
|---|
| 515 | SUBROUTINE put_var0d_nc(var_name,var,itime) |
|---|
| 516 | !----------------------------------------------------------------------- |
|---|
| 517 | ! NAME |
|---|
| 518 | ! put_var0d_nc |
|---|
| 519 | ! |
|---|
| 520 | ! DESCRIPTION |
|---|
| 521 | ! Write a 0D variable into open NetCDF file. |
|---|
| 522 | ! |
|---|
| 523 | ! AUTHORS & DATE |
|---|
| 524 | ! JB Clement, 01/2026 |
|---|
| 525 | ! |
|---|
| 526 | ! NOTES |
|---|
| 527 | ! |
|---|
| 528 | !----------------------------------------------------------------------- |
|---|
| 529 | |
|---|
| 530 | ! DECLARATION |
|---|
| 531 | ! ----------- |
|---|
| 532 | implicit none |
|---|
| 533 | |
|---|
| 534 | ! ARGUMENTS |
|---|
| 535 | ! --------- |
|---|
| 536 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 537 | real(dp), intent(in) :: var ! Input variable |
|---|
| 538 | integer(di), intent(in), optional :: itime ! Current time index |
|---|
| 539 | |
|---|
| 540 | ! LOCAL VARIABLES |
|---|
| 541 | ! --------------- |
|---|
| 542 | integer(di) :: ndims, unlimdimid, ncid |
|---|
| 543 | integer(di), dimension(nf90_max_var_dims) :: dimids |
|---|
| 544 | logical(k4) :: has_time |
|---|
| 545 | |
|---|
| 546 | ! CODE |
|---|
| 547 | ! ---- |
|---|
| 548 | ! Diagevol logic priority over standard logic |
|---|
| 549 | if (open_diagevol) then |
|---|
| 550 | ncid = ncid_diagevol |
|---|
| 551 | else |
|---|
| 552 | ncid = ncid_file |
|---|
| 553 | end if |
|---|
| 554 | |
|---|
| 555 | ! Check if the variable holds a Time dimension (unlimited dim) |
|---|
| 556 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name) |
|---|
| 557 | call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name) |
|---|
| 558 | if (ndims > 0) then |
|---|
| 559 | call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim') |
|---|
| 560 | has_time = (dimids(ndims) == unlimdimid) |
|---|
| 561 | else |
|---|
| 562 | has_time = .false. |
|---|
| 563 | end if |
|---|
| 564 | |
|---|
| 565 | if (present(itime)) then ! Time-dependent write |
|---|
| 566 | if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1) |
|---|
| 567 | ! For 0D variable with time, just write at the time index |
|---|
| 568 | call check_nc(nf90_put_var(ncid,varid,var,start = (/itime/)),'putting '//var_name) |
|---|
| 569 | else ! Static write |
|---|
| 570 | if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1) |
|---|
| 571 | call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name) |
|---|
| 572 | end if |
|---|
| 573 | |
|---|
| 574 | END SUBROUTINE put_var0d_nc |
|---|
| 575 | !======================================================================= |
|---|
| 576 | |
|---|
| 577 | !======================================================================= |
|---|
| 578 | SUBROUTINE put_var1d_nc(var_name,var,itime) |
|---|
| 579 | !----------------------------------------------------------------------- |
|---|
| 580 | ! NAME |
|---|
| 581 | ! put_var1d_nc |
|---|
| 582 | ! |
|---|
| 583 | ! DESCRIPTION |
|---|
| 584 | ! Write a 1D variable into open NetCDF file. |
|---|
| 585 | ! |
|---|
| 586 | ! AUTHORS & DATE |
|---|
| 587 | ! JB Clement, 01/2026 |
|---|
| 588 | ! |
|---|
| 589 | ! NOTES |
|---|
| 590 | ! |
|---|
| 591 | !----------------------------------------------------------------------- |
|---|
| 592 | |
|---|
| 593 | ! DECLARATION |
|---|
| 594 | ! ----------- |
|---|
| 595 | implicit none |
|---|
| 596 | |
|---|
| 597 | ! ARGUMENTS |
|---|
| 598 | ! --------- |
|---|
| 599 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 600 | real(dp), dimension(:), intent(in) :: var ! Input variable |
|---|
| 601 | integer(di), intent(in), optional :: itime ! Current time index |
|---|
| 602 | |
|---|
| 603 | ! LOCAL VARIABLES |
|---|
| 604 | ! --------------- |
|---|
| 605 | integer(di) :: i, ndims, unlimdimid, ncid |
|---|
| 606 | integer(di), dimension(nf90_max_var_dims) :: dimids |
|---|
| 607 | integer(di), dimension(:), allocatable :: strt, cnt |
|---|
| 608 | logical(k4) :: has_time |
|---|
| 609 | |
|---|
| 610 | ! CODE |
|---|
| 611 | ! ---- |
|---|
| 612 | ! Diagevol logic priority over standard logic |
|---|
| 613 | if (open_diagevol) then |
|---|
| 614 | ncid = ncid_diagevol |
|---|
| 615 | else |
|---|
| 616 | ncid = ncid_file |
|---|
| 617 | end if |
|---|
| 618 | |
|---|
| 619 | ! Check if the variable holds a Time dimension (unlimited dim) |
|---|
| 620 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name) |
|---|
| 621 | call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name) |
|---|
| 622 | if (ndims > 0) then |
|---|
| 623 | call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim') |
|---|
| 624 | has_time = (dimids(ndims) == unlimdimid) |
|---|
| 625 | else |
|---|
| 626 | has_time = .false. |
|---|
| 627 | end if |
|---|
| 628 | |
|---|
| 629 | if (present(itime)) then ! Time-dependent write |
|---|
| 630 | if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1) |
|---|
| 631 | allocate(strt(ndims),cnt(ndims)) |
|---|
| 632 | strt(:) = 1 |
|---|
| 633 | cnt(:) = 1 |
|---|
| 634 | strt(ndims) = itime |
|---|
| 635 | do i = 1,ndims - 1 |
|---|
| 636 | cnt(i) = size(var,i) |
|---|
| 637 | end do |
|---|
| 638 | call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name) |
|---|
| 639 | deallocate(strt,cnt) |
|---|
| 640 | else ! Static write |
|---|
| 641 | if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1) |
|---|
| 642 | call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name) |
|---|
| 643 | end if |
|---|
| 644 | |
|---|
| 645 | END SUBROUTINE put_var1d_nc |
|---|
| 646 | !======================================================================= |
|---|
| 647 | |
|---|
| 648 | !======================================================================= |
|---|
| 649 | SUBROUTINE put_var2d_nc(var_name,var,itime) |
|---|
| 650 | !----------------------------------------------------------------------- |
|---|
| 651 | ! NAME |
|---|
| 652 | ! put_var2d_nc |
|---|
| 653 | ! |
|---|
| 654 | ! DESCRIPTION |
|---|
| 655 | ! Write a 2D variable into open NetCDF file. |
|---|
| 656 | ! |
|---|
| 657 | ! AUTHORS & DATE |
|---|
| 658 | ! JB Clement, 01/2026 |
|---|
| 659 | ! |
|---|
| 660 | ! NOTES |
|---|
| 661 | ! |
|---|
| 662 | !----------------------------------------------------------------------- |
|---|
| 663 | |
|---|
| 664 | ! DECLARATION |
|---|
| 665 | ! ----------- |
|---|
| 666 | implicit none |
|---|
| 667 | |
|---|
| 668 | ! ARGUMENTS |
|---|
| 669 | ! --------- |
|---|
| 670 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 671 | real(dp), dimension(:,:), intent(in) :: var ! Input variable |
|---|
| 672 | integer(di), intent(in), optional :: itime ! Current time index |
|---|
| 673 | |
|---|
| 674 | ! LOCAL VARIABLES |
|---|
| 675 | ! --------------- |
|---|
| 676 | integer(di) :: i, ndims, unlimdimid, ncid |
|---|
| 677 | integer(di), dimension(nf90_max_var_dims) :: dimids |
|---|
| 678 | integer(di), dimension(:), allocatable :: strt, cnt |
|---|
| 679 | logical(k4) :: has_time |
|---|
| 680 | |
|---|
| 681 | ! CODE |
|---|
| 682 | ! ---- |
|---|
| 683 | ! Diagevol logic priority over standard logic |
|---|
| 684 | if (open_diagevol) then |
|---|
| 685 | ncid = ncid_diagevol |
|---|
| 686 | else |
|---|
| 687 | ncid = ncid_file |
|---|
| 688 | end if |
|---|
| 689 | |
|---|
| 690 | ! Check if the variable holds a Time dimension (unlimited dim) |
|---|
| 691 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name) |
|---|
| 692 | call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name) |
|---|
| 693 | if (ndims > 0) then |
|---|
| 694 | call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim') |
|---|
| 695 | has_time = (dimids(ndims) == unlimdimid) |
|---|
| 696 | else |
|---|
| 697 | has_time = .false. |
|---|
| 698 | end if |
|---|
| 699 | |
|---|
| 700 | if (present(itime)) then ! Time-dependent write |
|---|
| 701 | if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1) |
|---|
| 702 | allocate(strt(ndims),cnt(ndims)) |
|---|
| 703 | strt(:) = 1 |
|---|
| 704 | cnt(:) = 1 |
|---|
| 705 | strt(ndims) = itime |
|---|
| 706 | do i = 1,ndims - 1 |
|---|
| 707 | cnt(i) = size(var,i) |
|---|
| 708 | end do |
|---|
| 709 | call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name) |
|---|
| 710 | deallocate(strt,cnt) |
|---|
| 711 | else ! Static write |
|---|
| 712 | if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1) |
|---|
| 713 | call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name) |
|---|
| 714 | end if |
|---|
| 715 | |
|---|
| 716 | END SUBROUTINE put_var2d_nc |
|---|
| 717 | !======================================================================= |
|---|
| 718 | |
|---|
| 719 | !======================================================================= |
|---|
| 720 | SUBROUTINE put_var3d_nc(var_name,var,itime) |
|---|
| 721 | !----------------------------------------------------------------------- |
|---|
| 722 | ! NAME |
|---|
| 723 | ! put_var3d_nc |
|---|
| 724 | ! |
|---|
| 725 | ! DESCRIPTION |
|---|
| 726 | ! Write a 3D variable into open NetCDF file. |
|---|
| 727 | ! |
|---|
| 728 | ! AUTHORS & DATE |
|---|
| 729 | ! JB Clement, 01/2026 |
|---|
| 730 | ! |
|---|
| 731 | ! NOTES |
|---|
| 732 | ! |
|---|
| 733 | !----------------------------------------------------------------------- |
|---|
| 734 | |
|---|
| 735 | ! DECLARATION |
|---|
| 736 | ! ----------- |
|---|
| 737 | implicit none |
|---|
| 738 | |
|---|
| 739 | ! ARGUMENTS |
|---|
| 740 | ! --------- |
|---|
| 741 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 742 | real(dp), dimension(:,:,:), intent(in) :: var ! Input variable |
|---|
| 743 | integer(di), intent(in), optional :: itime ! Current time index |
|---|
| 744 | |
|---|
| 745 | ! LOCAL VARIABLES |
|---|
| 746 | ! --------------- |
|---|
| 747 | integer(di) :: i, ndims, unlimdimid, ncid |
|---|
| 748 | integer(di), dimension(nf90_max_var_dims) :: dimids |
|---|
| 749 | integer(di), dimension(:), allocatable :: strt, cnt |
|---|
| 750 | logical(k4) :: has_time |
|---|
| 751 | |
|---|
| 752 | ! CODE |
|---|
| 753 | ! ---- |
|---|
| 754 | ! Diagevol logic priority over standard logic |
|---|
| 755 | if (open_diagevol) then |
|---|
| 756 | ncid = ncid_diagevol |
|---|
| 757 | else |
|---|
| 758 | ncid = ncid_file |
|---|
| 759 | end if |
|---|
| 760 | |
|---|
| 761 | ! Check if the variable holds a Time dimension (unlimited dim) |
|---|
| 762 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name) |
|---|
| 763 | call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name) |
|---|
| 764 | if (ndims > 0) then |
|---|
| 765 | call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim') |
|---|
| 766 | has_time = (dimids(ndims) == unlimdimid) |
|---|
| 767 | else |
|---|
| 768 | has_time = .false. |
|---|
| 769 | end if |
|---|
| 770 | |
|---|
| 771 | if (present(itime)) then ! Time-dependent write |
|---|
| 772 | if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1) |
|---|
| 773 | allocate(strt(ndims),cnt(ndims)) |
|---|
| 774 | strt(:) = 1 |
|---|
| 775 | cnt(:) = 1 |
|---|
| 776 | strt(ndims) = itime |
|---|
| 777 | do i = 1,ndims - 1 |
|---|
| 778 | cnt(i) = size(var,i) |
|---|
| 779 | end do |
|---|
| 780 | call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name) |
|---|
| 781 | deallocate(strt,cnt) |
|---|
| 782 | else ! Static write |
|---|
| 783 | if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1) |
|---|
| 784 | call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name) |
|---|
| 785 | end if |
|---|
| 786 | |
|---|
| 787 | END SUBROUTINE put_var3d_nc |
|---|
| 788 | !======================================================================= |
|---|
| 789 | |
|---|
| 790 | !======================================================================= |
|---|
| 791 | SUBROUTINE put_var4d_nc(var_name,var,itime) |
|---|
| 792 | !----------------------------------------------------------------------- |
|---|
| 793 | ! NAME |
|---|
| 794 | ! put_var4d_nc |
|---|
| 795 | ! |
|---|
| 796 | ! DESCRIPTION |
|---|
| 797 | ! Write a 4D variable into open NetCDF file. |
|---|
| 798 | ! |
|---|
| 799 | ! AUTHORS & DATE |
|---|
| 800 | ! JB Clement, 01/2026 |
|---|
| 801 | ! |
|---|
| 802 | ! NOTES |
|---|
| 803 | ! |
|---|
| 804 | !----------------------------------------------------------------------- |
|---|
| 805 | |
|---|
| 806 | ! DECLARATION |
|---|
| 807 | ! ----------- |
|---|
| 808 | implicit none |
|---|
| 809 | |
|---|
| 810 | ! ARGUMENTS |
|---|
| 811 | ! --------- |
|---|
| 812 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 813 | real(dp), dimension(:,:,:,:), intent(in) :: var ! Input variable |
|---|
| 814 | integer(di), intent(in), optional :: itime ! Current time index |
|---|
| 815 | |
|---|
| 816 | ! LOCAL VARIABLES |
|---|
| 817 | ! --------------- |
|---|
| 818 | integer(di) :: i, ndims, unlimdimid, ncid |
|---|
| 819 | integer(di), dimension(nf90_max_var_dims) :: dimids |
|---|
| 820 | integer(di), dimension(:), allocatable :: strt, cnt |
|---|
| 821 | logical(k4) :: has_time |
|---|
| 822 | |
|---|
| 823 | ! CODE |
|---|
| 824 | ! ---- |
|---|
| 825 | ! Diagevol logic priority over standard logic |
|---|
| 826 | if (open_diagevol) then |
|---|
| 827 | ncid = ncid_diagevol |
|---|
| 828 | else |
|---|
| 829 | ncid = ncid_file |
|---|
| 830 | end if |
|---|
| 831 | |
|---|
| 832 | ! Check if the variable holds a Time dimension (unlimited dim) |
|---|
| 833 | call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name) |
|---|
| 834 | call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name) |
|---|
| 835 | if (ndims > 0) then |
|---|
| 836 | call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim') |
|---|
| 837 | has_time = (dimids(ndims) == unlimdimid) |
|---|
| 838 | else |
|---|
| 839 | has_time = .false. |
|---|
| 840 | end if |
|---|
| 841 | |
|---|
| 842 | if (present(itime)) then ! Time-dependent write |
|---|
| 843 | if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1) |
|---|
| 844 | allocate(strt(ndims),cnt(ndims)) |
|---|
| 845 | strt(:) = 1 |
|---|
| 846 | cnt(:) = 1 |
|---|
| 847 | strt(ndims) = itime |
|---|
| 848 | do i = 1,ndims - 1 |
|---|
| 849 | cnt(i) = size(var,i) |
|---|
| 850 | end do |
|---|
| 851 | call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name) |
|---|
| 852 | deallocate(strt,cnt) |
|---|
| 853 | else ! Static write |
|---|
| 854 | if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1) |
|---|
| 855 | call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name) |
|---|
| 856 | end if |
|---|
| 857 | |
|---|
| 858 | END SUBROUTINE put_var4d_nc |
|---|
| 859 | !======================================================================= |
|---|
| 860 | |
|---|
| 861 | !======================================================================= |
|---|
| 862 | SUBROUTINE get_next_itime_nc(dim_name,itime) |
|---|
| 863 | !----------------------------------------------------------------------- |
|---|
| 864 | ! NAME |
|---|
| 865 | ! get_next_itime_nc |
|---|
| 866 | ! |
|---|
| 867 | ! DESCRIPTION |
|---|
| 868 | ! Get the next time index in a NetCDF file to record variables. |
|---|
| 869 | ! |
|---|
| 870 | ! AUTHORS & DATE |
|---|
| 871 | ! JB Clement, 01/2026 |
|---|
| 872 | ! |
|---|
| 873 | ! NOTES |
|---|
| 874 | ! In most cases, dim_name = 'Time'. |
|---|
| 875 | !----------------------------------------------------------------------- |
|---|
| 876 | |
|---|
| 877 | ! DECLARATION |
|---|
| 878 | ! ----------- |
|---|
| 879 | implicit none |
|---|
| 880 | |
|---|
| 881 | ! ARGUMENTS |
|---|
| 882 | ! --------- |
|---|
| 883 | character(*), intent(in) :: dim_name ! Dimension name |
|---|
| 884 | integer(di), intent(out) :: itime ! Time index |
|---|
| 885 | |
|---|
| 886 | ! LOCAL VARIABLES |
|---|
| 887 | ! --------------- |
|---|
| 888 | integer(di) :: length |
|---|
| 889 | logical(k4) :: found |
|---|
| 890 | |
|---|
| 891 | ! CODE |
|---|
| 892 | ! ---- |
|---|
| 893 | ! Get dimension length |
|---|
| 894 | call get_dim_nc(dim_name,length,found) |
|---|
| 895 | if (.not. found) call stop_clean(__FILE__,__LINE__,'dimension '//dim_name//' not found',1) |
|---|
| 896 | |
|---|
| 897 | ! Next time index |
|---|
| 898 | itime = length + 1 |
|---|
| 899 | |
|---|
| 900 | END SUBROUTINE get_next_itime_nc |
|---|
| 901 | !======================================================================= |
|---|
| 902 | |
|---|
| 903 | !======================================================================= |
|---|
| 904 | SUBROUTINE check_valid_var0d_nc(var_name,var) |
|---|
| 905 | !----------------------------------------------------------------------- |
|---|
| 906 | ! NAME |
|---|
| 907 | ! check_valid_var0d_nc |
|---|
| 908 | ! |
|---|
| 909 | ! DESCRIPTION |
|---|
| 910 | ! Check the validity of a 0D variable read from a NetCDF file. |
|---|
| 911 | ! |
|---|
| 912 | ! AUTHORS & DATE |
|---|
| 913 | ! JB Clement, 02/2026 |
|---|
| 914 | ! |
|---|
| 915 | ! NOTES |
|---|
| 916 | ! |
|---|
| 917 | !----------------------------------------------------------------------- |
|---|
| 918 | |
|---|
| 919 | ! DEPENDENCIES |
|---|
| 920 | ! ------------ |
|---|
| 921 | use numerics, only: largest_nb |
|---|
| 922 | use stoppage, only: stop_clean |
|---|
| 923 | use utility, only: real2str |
|---|
| 924 | |
|---|
| 925 | ! DECLARATION |
|---|
| 926 | ! ----------- |
|---|
| 927 | implicit none |
|---|
| 928 | |
|---|
| 929 | ! ARGUMENTS |
|---|
| 930 | ! --------- |
|---|
| 931 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 932 | real(dp), intent(in) :: var ! Input variable |
|---|
| 933 | |
|---|
| 934 | ! LOCAL VARIABLES |
|---|
| 935 | ! --------------- |
|---|
| 936 | logical(k4) :: has_fill, has_range, has_valid_min, has_valid_max |
|---|
| 937 | real(dp) :: fill_value, valid_min, valid_max |
|---|
| 938 | real(dp), dimension(2) :: valid_range |
|---|
| 939 | integer(di) :: ncid |
|---|
| 940 | |
|---|
| 941 | ! CODE |
|---|
| 942 | ! ---- |
|---|
| 943 | ! Diagevol logic priority over standard logic |
|---|
| 944 | if (open_diagevol) then |
|---|
| 945 | ncid = ncid_diagevol |
|---|
| 946 | else |
|---|
| 947 | ncid = ncid_file |
|---|
| 948 | end if |
|---|
| 949 | |
|---|
| 950 | ! NaN |
|---|
| 951 | if (var /= var) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1) |
|---|
| 952 | |
|---|
| 953 | ! Infinite |
|---|
| 954 | if (abs(var) > largest_nb) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1) |
|---|
| 955 | |
|---|
| 956 | ! Fill value |
|---|
| 957 | has_fill = .false. |
|---|
| 958 | call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill) |
|---|
| 959 | if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill) |
|---|
| 960 | if (has_fill) then |
|---|
| 961 | if (abs(var - fill_value) < minieps) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1) |
|---|
| 962 | end if |
|---|
| 963 | |
|---|
| 964 | ! Valid range |
|---|
| 965 | has_range = .false. |
|---|
| 966 | call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range) |
|---|
| 967 | if (has_range) then |
|---|
| 968 | valid_min = valid_range(1) |
|---|
| 969 | valid_max = valid_range(2) |
|---|
| 970 | else! valid_min / valid_max (fallback) |
|---|
| 971 | has_valid_min = .false. |
|---|
| 972 | has_valid_max = .false. |
|---|
| 973 | call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min) |
|---|
| 974 | call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max) |
|---|
| 975 | if (has_valid_min .or. has_valid_max) then |
|---|
| 976 | has_range = .true. |
|---|
| 977 | if (.not. has_valid_min) valid_min = -largest_nb |
|---|
| 978 | if (.not. has_valid_max) valid_max = largest_nb |
|---|
| 979 | end if |
|---|
| 980 | end if |
|---|
| 981 | if (has_range) then |
|---|
| 982 | if (var < valid_min .or. var > valid_max) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1) |
|---|
| 983 | end if |
|---|
| 984 | |
|---|
| 985 | END SUBROUTINE check_valid_var0d_nc |
|---|
| 986 | !======================================================================= |
|---|
| 987 | |
|---|
| 988 | !======================================================================= |
|---|
| 989 | SUBROUTINE check_valid_var1d_nc(var_name,var) |
|---|
| 990 | !----------------------------------------------------------------------- |
|---|
| 991 | ! NAME |
|---|
| 992 | ! check_valid_var1d_nc |
|---|
| 993 | ! |
|---|
| 994 | ! DESCRIPTION |
|---|
| 995 | ! Check the validity of a 1D variable read from a NetCDF file. |
|---|
| 996 | ! |
|---|
| 997 | ! AUTHORS & DATE |
|---|
| 998 | ! JB Clement, 02/2026 |
|---|
| 999 | ! |
|---|
| 1000 | ! NOTES |
|---|
| 1001 | ! |
|---|
| 1002 | !----------------------------------------------------------------------- |
|---|
| 1003 | |
|---|
| 1004 | ! DEPENDENCIES |
|---|
| 1005 | ! ------------ |
|---|
| 1006 | use numerics, only: largest_nb |
|---|
| 1007 | use stoppage, only: stop_clean |
|---|
| 1008 | use utility, only: real2str |
|---|
| 1009 | |
|---|
| 1010 | ! DECLARATION |
|---|
| 1011 | ! ----------- |
|---|
| 1012 | implicit none |
|---|
| 1013 | |
|---|
| 1014 | ! ARGUMENTS |
|---|
| 1015 | ! --------- |
|---|
| 1016 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 1017 | real(dp), dimension(:), intent(in) :: var ! Input variable |
|---|
| 1018 | |
|---|
| 1019 | ! LOCAL VARIABLES |
|---|
| 1020 | ! --------------- |
|---|
| 1021 | logical(k4) :: has_fill, has_range, has_valid_min, has_valid_max |
|---|
| 1022 | real(dp) :: fill_value, valid_min, valid_max |
|---|
| 1023 | real(dp), dimension(2) :: valid_range |
|---|
| 1024 | integer(di) :: ncid |
|---|
| 1025 | |
|---|
| 1026 | ! CODE |
|---|
| 1027 | ! ---- |
|---|
| 1028 | ! Diagevol logic priority over standard logic |
|---|
| 1029 | if (open_diagevol) then |
|---|
| 1030 | ncid = ncid_diagevol |
|---|
| 1031 | else |
|---|
| 1032 | ncid = ncid_file |
|---|
| 1033 | end if |
|---|
| 1034 | |
|---|
| 1035 | ! NaN |
|---|
| 1036 | if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1) |
|---|
| 1037 | |
|---|
| 1038 | ! Infinite |
|---|
| 1039 | if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1) |
|---|
| 1040 | |
|---|
| 1041 | ! Fill value |
|---|
| 1042 | has_fill = .false. |
|---|
| 1043 | call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill) |
|---|
| 1044 | if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill) |
|---|
| 1045 | if (has_fill) then |
|---|
| 1046 | if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1) |
|---|
| 1047 | end if |
|---|
| 1048 | |
|---|
| 1049 | ! Valid range |
|---|
| 1050 | has_range = .false. |
|---|
| 1051 | call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range) |
|---|
| 1052 | if (has_range) then |
|---|
| 1053 | valid_min = valid_range(1) |
|---|
| 1054 | valid_max = valid_range(2) |
|---|
| 1055 | else! valid_min / valid_max (fallback) |
|---|
| 1056 | has_valid_min = .false. |
|---|
| 1057 | has_valid_max = .false. |
|---|
| 1058 | call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min) |
|---|
| 1059 | call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max) |
|---|
| 1060 | if (has_valid_min .or. has_valid_max) then |
|---|
| 1061 | has_range = .true. |
|---|
| 1062 | if (.not. has_valid_min) valid_min = -largest_nb |
|---|
| 1063 | if (.not. has_valid_max) valid_max = largest_nb |
|---|
| 1064 | end if |
|---|
| 1065 | end if |
|---|
| 1066 | if (has_range) then |
|---|
| 1067 | if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1) |
|---|
| 1068 | end if |
|---|
| 1069 | |
|---|
| 1070 | END SUBROUTINE check_valid_var1d_nc |
|---|
| 1071 | !======================================================================= |
|---|
| 1072 | |
|---|
| 1073 | !======================================================================= |
|---|
| 1074 | SUBROUTINE check_valid_var2d_nc(var_name,var) |
|---|
| 1075 | !----------------------------------------------------------------------- |
|---|
| 1076 | ! NAME |
|---|
| 1077 | ! check_valid_var2d_nc |
|---|
| 1078 | ! |
|---|
| 1079 | ! DESCRIPTION |
|---|
| 1080 | ! Check the validity of a 1D variable read from a NetCDF file. |
|---|
| 1081 | ! |
|---|
| 1082 | ! AUTHORS & DATE |
|---|
| 1083 | ! JB Clement, 02/2026 |
|---|
| 1084 | ! |
|---|
| 1085 | ! NOTES |
|---|
| 1086 | ! |
|---|
| 1087 | !----------------------------------------------------------------------- |
|---|
| 1088 | |
|---|
| 1089 | ! DEPENDENCIES |
|---|
| 1090 | ! ------------ |
|---|
| 1091 | use numerics, only: largest_nb |
|---|
| 1092 | use stoppage, only: stop_clean |
|---|
| 1093 | use utility, only: real2str |
|---|
| 1094 | |
|---|
| 1095 | ! DECLARATION |
|---|
| 1096 | ! ----------- |
|---|
| 1097 | implicit none |
|---|
| 1098 | |
|---|
| 1099 | ! ARGUMENTS |
|---|
| 1100 | ! --------- |
|---|
| 1101 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 1102 | real(dp), dimension(:,:), intent(in) :: var ! Input variable |
|---|
| 1103 | |
|---|
| 1104 | ! LOCAL VARIABLES |
|---|
| 1105 | ! --------------- |
|---|
| 1106 | logical(k4) :: has_fill, has_range, has_valid_min, has_valid_max |
|---|
| 1107 | real(dp) :: fill_value, valid_min, valid_max |
|---|
| 1108 | real(dp), dimension(2) :: valid_range |
|---|
| 1109 | integer(di) :: ncid |
|---|
| 1110 | |
|---|
| 1111 | ! CODE |
|---|
| 1112 | ! ---- |
|---|
| 1113 | ! Diagevol logic priority over standard logic |
|---|
| 1114 | if (open_diagevol) then |
|---|
| 1115 | ncid = ncid_diagevol |
|---|
| 1116 | else |
|---|
| 1117 | ncid = ncid_file |
|---|
| 1118 | end if |
|---|
| 1119 | |
|---|
| 1120 | ! NaN |
|---|
| 1121 | if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1) |
|---|
| 1122 | |
|---|
| 1123 | ! Infinite |
|---|
| 1124 | if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1) |
|---|
| 1125 | |
|---|
| 1126 | ! Fill value |
|---|
| 1127 | has_fill = .false. |
|---|
| 1128 | call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill) |
|---|
| 1129 | if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill) |
|---|
| 1130 | if (has_fill) then |
|---|
| 1131 | if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1) |
|---|
| 1132 | end if |
|---|
| 1133 | |
|---|
| 1134 | ! Valid range |
|---|
| 1135 | has_range = .false. |
|---|
| 1136 | call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range) |
|---|
| 1137 | if (has_range) then |
|---|
| 1138 | valid_min = valid_range(1) |
|---|
| 1139 | valid_max = valid_range(2) |
|---|
| 1140 | else! valid_min / valid_max (fallback) |
|---|
| 1141 | has_valid_min = .false. |
|---|
| 1142 | has_valid_max = .false. |
|---|
| 1143 | call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min) |
|---|
| 1144 | call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max) |
|---|
| 1145 | if (has_valid_min .or. has_valid_max) then |
|---|
| 1146 | has_range = .true. |
|---|
| 1147 | if (.not. has_valid_min) valid_min = -largest_nb |
|---|
| 1148 | if (.not. has_valid_max) valid_max = largest_nb |
|---|
| 1149 | end if |
|---|
| 1150 | end if |
|---|
| 1151 | if (has_range) then |
|---|
| 1152 | if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1) |
|---|
| 1153 | end if |
|---|
| 1154 | |
|---|
| 1155 | END SUBROUTINE check_valid_var2d_nc |
|---|
| 1156 | !======================================================================= |
|---|
| 1157 | |
|---|
| 1158 | !======================================================================= |
|---|
| 1159 | SUBROUTINE check_valid_var3d_nc(var_name,var) |
|---|
| 1160 | !----------------------------------------------------------------------- |
|---|
| 1161 | ! NAME |
|---|
| 1162 | ! check_valid_var3d_nc |
|---|
| 1163 | ! |
|---|
| 1164 | ! DESCRIPTION |
|---|
| 1165 | ! Check the validity of a 1D variable read from a NetCDF file. |
|---|
| 1166 | ! |
|---|
| 1167 | ! AUTHORS & DATE |
|---|
| 1168 | ! JB Clement, 02/2026 |
|---|
| 1169 | ! |
|---|
| 1170 | ! NOTES |
|---|
| 1171 | ! |
|---|
| 1172 | !----------------------------------------------------------------------- |
|---|
| 1173 | |
|---|
| 1174 | ! DEPENDENCIES |
|---|
| 1175 | ! ------------ |
|---|
| 1176 | use numerics, only: largest_nb |
|---|
| 1177 | use stoppage, only: stop_clean |
|---|
| 1178 | use utility, only: real2str |
|---|
| 1179 | |
|---|
| 1180 | ! DECLARATION |
|---|
| 1181 | ! ----------- |
|---|
| 1182 | implicit none |
|---|
| 1183 | |
|---|
| 1184 | ! ARGUMENTS |
|---|
| 1185 | ! --------- |
|---|
| 1186 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 1187 | real(dp), dimension(:,:,:), intent(in) :: var ! Input variable |
|---|
| 1188 | |
|---|
| 1189 | ! LOCAL VARIABLES |
|---|
| 1190 | ! --------------- |
|---|
| 1191 | logical(k4) :: has_fill, has_range, has_valid_min, has_valid_max |
|---|
| 1192 | real(dp) :: fill_value, valid_min, valid_max |
|---|
| 1193 | real(dp), dimension(2) :: valid_range |
|---|
| 1194 | integer(di) :: ncid |
|---|
| 1195 | |
|---|
| 1196 | ! CODE |
|---|
| 1197 | ! ---- |
|---|
| 1198 | ! Diagevol logic priority over standard logic |
|---|
| 1199 | if (open_diagevol) then |
|---|
| 1200 | ncid = ncid_diagevol |
|---|
| 1201 | else |
|---|
| 1202 | ncid = ncid_file |
|---|
| 1203 | end if |
|---|
| 1204 | |
|---|
| 1205 | ! NaN |
|---|
| 1206 | if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1) |
|---|
| 1207 | |
|---|
| 1208 | ! Infinite |
|---|
| 1209 | if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1) |
|---|
| 1210 | |
|---|
| 1211 | ! Fill value |
|---|
| 1212 | has_fill = .false. |
|---|
| 1213 | call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill) |
|---|
| 1214 | if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill) |
|---|
| 1215 | if (has_fill) then |
|---|
| 1216 | if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1) |
|---|
| 1217 | end if |
|---|
| 1218 | |
|---|
| 1219 | ! Valid range |
|---|
| 1220 | has_range = .false. |
|---|
| 1221 | call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range) |
|---|
| 1222 | if (has_range) then |
|---|
| 1223 | valid_min = valid_range(1) |
|---|
| 1224 | valid_max = valid_range(2) |
|---|
| 1225 | else! valid_min / valid_max (fallback) |
|---|
| 1226 | has_valid_min = .false. |
|---|
| 1227 | has_valid_max = .false. |
|---|
| 1228 | call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min) |
|---|
| 1229 | call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max) |
|---|
| 1230 | if (has_valid_min .or. has_valid_max) then |
|---|
| 1231 | has_range = .true. |
|---|
| 1232 | if (.not. has_valid_min) valid_min = -largest_nb |
|---|
| 1233 | if (.not. has_valid_max) valid_max = largest_nb |
|---|
| 1234 | end if |
|---|
| 1235 | end if |
|---|
| 1236 | if (has_range) then |
|---|
| 1237 | if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1) |
|---|
| 1238 | end if |
|---|
| 1239 | |
|---|
| 1240 | END SUBROUTINE check_valid_var3d_nc |
|---|
| 1241 | !======================================================================= |
|---|
| 1242 | |
|---|
| 1243 | !======================================================================= |
|---|
| 1244 | SUBROUTINE check_valid_var4d_nc(var_name,var) |
|---|
| 1245 | !----------------------------------------------------------------------- |
|---|
| 1246 | ! NAME |
|---|
| 1247 | ! check_valid_var4d_nc |
|---|
| 1248 | ! |
|---|
| 1249 | ! DESCRIPTION |
|---|
| 1250 | ! Check the validity of a 1D variable read from a NetCDF file. |
|---|
| 1251 | ! |
|---|
| 1252 | ! AUTHORS & DATE |
|---|
| 1253 | ! JB Clement, 02/2026 |
|---|
| 1254 | ! |
|---|
| 1255 | ! NOTES |
|---|
| 1256 | ! |
|---|
| 1257 | !----------------------------------------------------------------------- |
|---|
| 1258 | |
|---|
| 1259 | ! DEPENDENCIES |
|---|
| 1260 | ! ------------ |
|---|
| 1261 | use numerics, only: largest_nb |
|---|
| 1262 | use stoppage, only: stop_clean |
|---|
| 1263 | use utility, only: real2str |
|---|
| 1264 | |
|---|
| 1265 | ! DECLARATION |
|---|
| 1266 | ! ----------- |
|---|
| 1267 | implicit none |
|---|
| 1268 | |
|---|
| 1269 | ! ARGUMENTS |
|---|
| 1270 | ! --------- |
|---|
| 1271 | character(*), intent(in) :: var_name ! Variable name |
|---|
| 1272 | real(dp), dimension(:,:,:,:), intent(in) :: var ! Input variable |
|---|
| 1273 | |
|---|
| 1274 | ! LOCAL VARIABLES |
|---|
| 1275 | ! --------------- |
|---|
| 1276 | logical(k4) :: has_fill, has_range, has_valid_min, has_valid_max |
|---|
| 1277 | real(dp) :: fill_value, valid_min, valid_max |
|---|
| 1278 | real(dp), dimension(2) :: valid_range |
|---|
| 1279 | integer(di) :: ncid |
|---|
| 1280 | |
|---|
| 1281 | ! CODE |
|---|
| 1282 | ! ---- |
|---|
| 1283 | ! Diagevol logic priority over standard logic |
|---|
| 1284 | if (open_diagevol) then |
|---|
| 1285 | ncid = ncid_diagevol |
|---|
| 1286 | else |
|---|
| 1287 | ncid = ncid_file |
|---|
| 1288 | end if |
|---|
| 1289 | |
|---|
| 1290 | ! NaN |
|---|
| 1291 | if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1) |
|---|
| 1292 | |
|---|
| 1293 | ! Infinite |
|---|
| 1294 | if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1) |
|---|
| 1295 | |
|---|
| 1296 | ! Fill value |
|---|
| 1297 | has_fill = .false. |
|---|
| 1298 | call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill) |
|---|
| 1299 | if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill) |
|---|
| 1300 | if (has_fill) then |
|---|
| 1301 | if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1) |
|---|
| 1302 | end if |
|---|
| 1303 | |
|---|
| 1304 | ! Valid range |
|---|
| 1305 | has_range = .false. |
|---|
| 1306 | call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range) |
|---|
| 1307 | if (has_range) then |
|---|
| 1308 | valid_min = valid_range(1) |
|---|
| 1309 | valid_max = valid_range(2) |
|---|
| 1310 | else! valid_min / valid_max (fallback) |
|---|
| 1311 | has_valid_min = .false. |
|---|
| 1312 | has_valid_max = .false. |
|---|
| 1313 | call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min) |
|---|
| 1314 | call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max) |
|---|
| 1315 | if (has_valid_min .or. has_valid_max) then |
|---|
| 1316 | has_range = .true. |
|---|
| 1317 | if (.not. has_valid_min) valid_min = -largest_nb |
|---|
| 1318 | if (.not. has_valid_max) valid_max = largest_nb |
|---|
| 1319 | end if |
|---|
| 1320 | end if |
|---|
| 1321 | if (has_range) then |
|---|
| 1322 | if(any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1) |
|---|
| 1323 | end if |
|---|
| 1324 | |
|---|
| 1325 | END SUBROUTINE check_valid_var4d_nc |
|---|
| 1326 | !======================================================================= |
|---|
| 1327 | |
|---|
| 1328 | END MODULE io_netcdf |
|---|