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
|
|
|
!
|
|
|
|
!
|
2015-03-17 13:19:09 +00:00
|
|
|
program grib_get_keys
|
2021-02-14 18:14:39 +00:00
|
|
|
use eccodes
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer :: ifile
|
|
|
|
integer :: iret
|
|
|
|
integer :: igrib
|
|
|
|
real :: latitudeOfFirstPointInDegrees
|
|
|
|
real :: longitudeOfFirstPointInDegrees
|
|
|
|
real :: latitudeOfLastPointInDegrees
|
|
|
|
real :: longitudeOfLastPointInDegrees
|
|
|
|
integer :: numberOfPointsAlongAParallel
|
|
|
|
integer :: numberOfPointsAlongAMeridian
|
|
|
|
real, dimension(:), allocatable :: values
|
|
|
|
integer :: numberOfValues
|
|
|
|
real :: average, min_val, max_val
|
|
|
|
integer :: is_missing
|
|
|
|
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)
|
2021-02-14 18:14:39 +00:00
|
|
|
|
2023-08-17 10:37:42 +00:00
|
|
|
! 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
|
2021-02-14 18:14:39 +00:00
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! check if the value of the key is MISSING
|
2022-12-22 16:26:01 +00:00
|
|
|
is_missing = 0;
|
|
|
|
call codes_is_missing(igrib, 'Ni', is_missing);
|
2021-02-14 18:14:39 +00:00
|
|
|
if (is_missing /= 1) then
|
2021-11-30 13:01:07 +00:00
|
|
|
! key value is not missing so get as an integer
|
2021-02-14 18:14:39 +00:00
|
|
|
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
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_get(igrib, 'Nj', numberOfPointsAlongAMeridian)
|
|
|
|
write (*, *) 'numberOfPointsAlongAMeridian=', &
|
2015-05-13 12:00:31 +00:00
|
|
|
numberOfPointsAlongAMeridian
|
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! get as a real
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_get(igrib, 'latitudeOfFirstGridPointInDegrees', &
|
|
|
|
latitudeOfFirstPointInDegrees)
|
|
|
|
write (*, *) 'latitudeOfFirstGridPointInDegrees=', &
|
|
|
|
latitudeOfFirstPointInDegrees
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! get as a real
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_get(igrib, 'longitudeOfFirstGridPointInDegrees', &
|
|
|
|
longitudeOfFirstPointInDegrees)
|
|
|
|
write (*, *) 'longitudeOfFirstGridPointInDegrees=', &
|
|
|
|
longitudeOfFirstPointInDegrees
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! get as a real
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_get(igrib, 'latitudeOfLastGridPointInDegrees', &
|
|
|
|
latitudeOfLastPointInDegrees)
|
|
|
|
write (*, *) 'latitudeOfLastGridPointInDegrees=', &
|
|
|
|
latitudeOfLastPointInDegrees
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! get as a real
|
2021-02-14 18:14:39 +00:00
|
|
|
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
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_get_size(igrib, 'values', numberOfValues)
|
|
|
|
write (*, *) 'numberOfValues=', numberOfValues
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
allocate (values(numberOfValues), stat=iret)
|
2021-11-30 13:01:07 +00:00
|
|
|
! get data values
|
2021-02-14 18:14:39 +00:00
|
|
|
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)
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
deallocate (values)
|
2016-07-13 16:10:34 +00:00
|
|
|
|
2021-02-14 18:14:39 +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
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_release(igrib)
|
2016-07-13 16:10:34 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
end do LOOP
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_close_file(ifile)
|
2015-05-13 12:00:31 +00:00
|
|
|
|
2015-03-17 13:19:09 +00:00
|
|
|
end program grib_get_keys
|