eccodes/examples/F90/grib_count_messages.f90

109 lines
3.9 KiB
Fortran
Raw Permalink 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: Count messages before processing
2013-03-25 12:04:10 +00:00
!
!
2024-01-13 17:14:08 +00:00
program grib_count_messages
use eccodes
implicit none
integer :: ifile
integer :: iret
integer :: n
integer :: i
2024-01-13 17:14:08 +00:00
integer, dimension(:), allocatable :: igrib
real :: latitudeOfFirstPointInDegrees
real :: longitudeOfFirstPointInDegrees
real :: latitudeOfLastPointInDegrees
real :: longitudeOfLastPointInDegrees
integer :: numberOfPointsAlongAParallel
integer :: numberOfPointsAlongAMeridian
real, dimension(:), allocatable :: values
2024-01-13 17:14:08 +00:00
integer(8) :: numberOfValues
real :: average, min_val, max_val
2024-01-13 17:14:08 +00:00
call codes_open_file(ifile, '../../data/tigge_pf_ecmwf.grib2', 'r')
! count the messages in the file
call codes_count_in_file(ifile, n)
allocate (igrib(n))
igrib = -1
2021-11-30 13:01:07 +00:00
! load the messages from the file.
DO i = 1, n
call codes_grib_new_from_file(ifile, igrib(i), iret)
END DO
! we can close the file
call codes_close_file(ifile)
2021-11-30 13:01:07 +00:00
! loop on all the messages in memory
DO i = 1, n
write (*, *) 'processing message number ', i
2021-11-30 13:01:07 +00:00
! get as a integer
call codes_get(igrib(i), 'Ni', numberOfPointsAlongAParallel)
write (*, *) 'numberOfPointsAlongAParallel=', &
numberOfPointsAlongAParallel
2021-11-30 13:01:07 +00:00
! get as a integer
call codes_get(igrib(i), 'Nj', numberOfPointsAlongAMeridian)
write (*, *) 'numberOfPointsAlongAMeridian=', &
numberOfPointsAlongAMeridian
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib(i), 'latitudeOfFirstGridPointInDegrees', &
latitudeOfFirstPointInDegrees)
write (*, *) 'latitudeOfFirstGridPointInDegrees=', &
latitudeOfFirstPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib(i), 'longitudeOfFirstGridPointInDegrees', &
longitudeOfFirstPointInDegrees)
write (*, *) 'longitudeOfFirstGridPointInDegrees=', &
longitudeOfFirstPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib(i), 'latitudeOfLastGridPointInDegrees', &
latitudeOfLastPointInDegrees)
write (*, *) 'latitudeOfLastGridPointInDegrees=', &
latitudeOfLastPointInDegrees
2021-11-30 13:01:07 +00:00
! get as a real
call codes_get(igrib(i), 'longitudeOfLastGridPointInDegrees', &
longitudeOfLastPointInDegrees)
write (*, *) 'longitudeOfLastGridPointInDegrees=', &
longitudeOfLastPointInDegrees
2021-11-30 13:01:07 +00:00
! get the size of the values array
call codes_get_size(igrib(i), '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(i), 'values', values)
call codes_get(igrib(i), 'min', min_val) ! can also be obtained through minval(values)
call codes_get(igrib(i), 'max', max_val) ! can also be obtained through maxval(values)
call codes_get(igrib(i), 'average', average) ! can also be obtained through maxval(values)
write (*, *) 'There are ', numberOfValues, &
' average is ', average, &
' min is ', min_val, &
' max is ', max_val
write (*, *) '---------------------'
deallocate (values)
END DO
DO i = 1, n
call codes_release(igrib(i))
END DO
deallocate (igrib)
2013-03-25 12:04:10 +00:00
2024-01-13 17:14:08 +00:00
end program