eccodes/examples/F90/grib_get_data.f90

82 lines
2.3 KiB
Fortran
Raw Normal View History

2017-01-03 11:03:48 +00:00
! Copyright 2005-2017 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.
!
!
! Description: how to get lat/lon/values.
!
!
program get_data
2015-01-28 22:54:42 +00:00
use eccodes
2013-03-25 12:04:10 +00:00
implicit none
integer :: ifile
integer :: iret,i
real(kind=8),dimension(:),allocatable :: lats,lons,values
2015-12-29 16:37:03 +00:00
integer,dimension(:),allocatable :: bitmap
2013-03-25 12:04:10 +00:00
integer(4) :: numberOfPoints
2015-12-29 17:29:00 +00:00
logical :: is_missing_value
integer :: count1=0, count2=0, bitmapPresent=0, bmp_len=0
2013-03-25 12:04:10 +00:00
! Message identifier.
2013-03-25 12:04:10 +00:00
integer :: igrib
ifile=5
call codes_open_file(ifile, &
2013-03-25 12:04:10 +00:00
'../../data/reduced_latlon_surface.grib1','R')
2015-12-29 16:37:03 +00:00
! Loop on all the messages in a file.
2015-03-04 17:11:18 +00:00
call codes_grib_new_from_file(ifile,igrib,iret)
2013-03-25 12:04:10 +00:00
do while (iret/=CODES_END_OF_FILE)
2013-03-25 12:04:10 +00:00
count1=count1+1
print *, "===== Message #",count1
call codes_get(igrib,'numberOfPoints',numberOfPoints)
2015-12-29 16:37:03 +00:00
call codes_get(igrib,'bitmapPresent',bitmapPresent)
2013-03-25 12:04:10 +00:00
allocate(lats(numberOfPoints))
allocate(lons(numberOfPoints))
allocate(values(numberOfPoints))
2015-12-29 16:37:03 +00:00
if (bitmapPresent == 1) then
! get the bitmap
call codes_get_size(igrib, 'bitmap', bmp_len)
allocate(bitmap(bmp_len))
call codes_get(igrib,'bitmap', bitmap)
end if
2013-03-25 12:04:10 +00:00
call codes_grib_get_data(igrib,lats,lons,values)
2013-03-25 12:04:10 +00:00
do i=1,numberOfPoints
2015-12-29 16:37:03 +00:00
! Consult bitmap to see if the i'th value is missing
2015-12-29 17:29:00 +00:00
is_missing_value=.false.
2015-12-29 16:37:03 +00:00
if (bitmapPresent == 1 .and. bitmap(i) == 0) then
2015-12-29 17:29:00 +00:00
is_missing_value=.true.
2015-12-29 16:37:03 +00:00
end if
! Only print non-missing values
2015-12-29 17:29:00 +00:00
if (.not. is_missing_value) then
2013-03-25 12:04:10 +00:00
print *, lats(i),lons(i),values(i)
2015-12-29 17:29:00 +00:00
count2=count2+1
2013-03-25 12:04:10 +00:00
end if
enddo
2015-12-29 17:29:00 +00:00
print *, 'count of non-missing values=',count2
if (count2 /= 214661) then
call codes_check(-2, 'incorrect number of missing', '')
end if
2013-03-25 12:04:10 +00:00
deallocate(lats)
deallocate(lons)
deallocate(values)
call codes_release(igrib)
2015-03-04 17:11:18 +00:00
call codes_grib_new_from_file(ifile,igrib, iret)
2013-03-25 12:04:10 +00:00
2017-02-06 16:45:30 +00:00
end do
2013-03-25 12:04:10 +00:00
call codes_close_file(ifile)
2013-03-25 12:04:10 +00:00
2017-02-06 16:45:30 +00:00
end program