eccodes/examples/F90/grib_set_bitmap.f90

79 lines
2.4 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 set a bitmap in a GRIB message
2016-11-21 13:10:21 +00:00
! to encode missing values in the data
2013-03-25 12:04:10 +00:00
!
!
program set_bitmap
use eccodes
implicit none
integer :: infile, outfile
integer :: igrib, iret
integer :: numberOfValues
real, dimension(:), allocatable :: values
real :: missingValue
logical :: grib1Example
2013-03-25 12:04:10 +00:00
grib1Example = .true.
2013-03-25 12:04:10 +00:00
if (grib1Example) then
! GRIB 1 example
call codes_open_file(infile, '../../data/regular_latlon_surface.grib1', 'r')
else
! GRIB 2 example
call codes_open_file(infile, '../../data/regular_latlon_surface.grib2', 'r')
end if
2017-02-06 16:45:30 +00:00
call codes_open_file(outfile, 'out.set_bitmap_f.grib', 'w')
2017-02-06 16:45:30 +00:00
2021-11-30 13:01:07 +00:00
! 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(infile, igrib)
2017-02-06 16:45:30 +00:00
2021-11-30 13:01:07 +00:00
! the missingValue is not coded in the message.
! It is a value we define as a placeholder for a missing value
! at a point in the grid.
! It should be chosen so that it cannot be confused
! with a valid field value
missingValue = 1.0e36
call codes_set(igrib, 'missingValue', missingValue)
write (*, *) 'missingValue=', missingValue
2013-03-25 12:04:10 +00:00
! get the size of the values array
call codes_get_size(igrib, 'values', numberOfValues)
write (*, *) 'numberOfValues=', numberOfValues
2017-02-06 16:45:30 +00:00
allocate (values(numberOfValues), stat=iret)
2013-03-25 12:04:10 +00:00
! get data values
call codes_get(igrib, 'values', values)
2017-02-06 16:45:30 +00:00
! enable bitmap
call codes_set(igrib, 'bitmapPresent', 1)
2013-03-25 12:04:10 +00:00
! set some values to be missing
values(1:10) = missingValue
2013-03-25 12:04:10 +00:00
! set the values (the bitmap will be automatically built)
call codes_set(igrib, 'values', values)
2013-03-25 12:04:10 +00:00
! write modified message to a file
call codes_write(igrib, outfile)
2017-02-06 16:45:30 +00:00
! free memory
call codes_release(igrib)
2017-02-06 16:45:30 +00:00
call codes_close_file(infile)
call codes_close_file(outfile)
2013-03-25 12:04:10 +00:00
deallocate (values)
2013-03-25 12:04:10 +00:00
end program set_bitmap