eccodes/examples/F90/grib_get_keys.f90

112 lines
4.2 KiB
Fortran
Raw Normal View History

2020-01-28 14:32:34 +00:00
! (C) Copyright 2005- ECMWF.
2013-03-25 12:04:10 +00:00
!
! This software is licensed under the terms of the Apache Licence Version 2.0
! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
2017-02-06 16:45:30 +00:00
!
2013-03-25 12:04:10 +00:00
! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
!
!
2021-02-23 17:14:11 +00:00
! Description: How to get values using keys from GRIB messages
2013-03-25 12:04:10 +00:00
!
!
program grib_get_keys
use eccodes
implicit none
integer :: ifile
integer :: iret
integer :: igrib
real :: latitudeOfFirstPointInDegrees
real :: longitudeOfFirstPointInDegrees
real :: latitudeOfLastPointInDegrees
real :: longitudeOfLastPointInDegrees
2024-01-13 17:14:08 +00:00
integer(4) :: numberOfPointsAlongAParallel
integer(8) :: numberOfPointsAlongAMeridian
real, dimension(:), allocatable :: values
integer :: numberOfValues
real :: average, min_val, max_val
2023-08-23 19:07:46 +00:00
integer :: is_missing, is_defined
character(len=10) :: open_mode = 'r'
2021-11-30 13:01:07 +00:00
call codes_open_file(ifile, '../../data/reduced_latlon_surface.grib1', open_mode)
! loop on all the messages in a file
LOOP: DO WHILE (.true.)
! a new GRIB message is loaded from file.
! igrib is the grib id to be used in subsequent calls
call codes_grib_new_from_file(ifile, igrib, iret)
if (iret == CODES_END_OF_FILE) exit LOOP
2023-08-23 19:07:46 +00:00
! check key is defined
is_defined = 0
call codes_is_defined(igrib, 'Ni', is_defined)
write (*, *) 'Key Ni is defined? ', is_defined
2021-11-30 13:01:07 +00:00
! check if the value of the key is MISSING
2023-08-23 19:07:46 +00:00
is_missing = 0
call codes_is_missing(igrib, 'Ni', is_missing)
if (is_missing /= 1) then
2021-11-30 13:01:07 +00:00
! key value is not missing so get as an integer
call codes_get(igrib, 'Ni', numberOfPointsAlongAParallel)
write (*, *) 'numberOfPointsAlongAParallel=', &
numberOfPointsAlongAParallel
else
write (*, *) 'numberOfPointsAlongAParallel is missing'
end if
2021-11-30 13:01:07 +00:00
! get as an integer
call codes_get(igrib, 'Nj', numberOfPointsAlongAMeridian)
write (*, *) 'numberOfPointsAlongAMeridian=', &
numberOfPointsAlongAMeridian
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib, 'latitudeOfFirstGridPointInDegrees', &
latitudeOfFirstPointInDegrees)
write (*, *) 'latitudeOfFirstGridPointInDegrees=', &
latitudeOfFirstPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib, 'longitudeOfFirstGridPointInDegrees', &
longitudeOfFirstPointInDegrees)
write (*, *) 'longitudeOfFirstGridPointInDegrees=', &
longitudeOfFirstPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib, 'latitudeOfLastGridPointInDegrees', &
latitudeOfLastPointInDegrees)
write (*, *) 'latitudeOfLastGridPointInDegrees=', &
latitudeOfLastPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib, 'longitudeOfLastGridPointInDegrees', &
longitudeOfLastPointInDegrees)
write (*, *) 'longitudeOfLastGridPointInDegrees=', &
longitudeOfLastPointInDegrees
2013-03-25 12:04:10 +00:00
2021-11-30 13:01:07 +00:00
! get the size of the values array
call codes_get_size(igrib, 'values', numberOfValues)
write (*, *) 'numberOfValues=', numberOfValues
allocate (values(numberOfValues), stat=iret)
2021-11-30 13:01:07 +00:00
! get data values
call codes_get(igrib, 'values', values)
call codes_get(igrib, 'min', min_val) ! can also be obtained through minval(values)
call codes_get(igrib, 'max', max_val) ! can also be obtained through maxval(values)
call codes_get(igrib, 'average', average) ! can also be obtained through maxval(values)
deallocate (values)
2016-07-13 16:10:34 +00:00
write (*, *) 'There are ', numberOfValues, &
' average is ', average, &
' min is ', min_val, &
' max is ', max_val
2016-07-13 16:10:34 +00:00
call codes_release(igrib)
2016-07-13 16:10:34 +00:00
end do LOOP
call codes_close_file(ifile)
end program grib_get_keys