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 pv values in a GRIB message
|
2013-03-25 12:04:10 +00:00
|
|
|
!
|
|
|
|
!
|
2015-03-17 17:39:55 +00:00
|
|
|
program grib_set_pv
|
2021-02-14 18:14:39 +00:00
|
|
|
use eccodes
|
|
|
|
implicit none
|
|
|
|
integer :: numberOfLevels
|
|
|
|
integer :: numberOfCoefficients
|
|
|
|
integer :: outfile, igrib
|
|
|
|
integer :: i, ios
|
|
|
|
real, dimension(:), allocatable :: pv
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
numberOfLevels = 60
|
|
|
|
numberOfCoefficients = 2*(numberOfLevels + 1)
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
allocate (pv(numberOfCoefficients))
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! read the model level coefficients from file
|
|
|
|
open (unit=1, file="../../data/60_model_levels", &
|
|
|
|
form="formatted", action="read")
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
do i = 1, numberOfCoefficients, 2
|
|
|
|
read (unit=1, fmt=*, iostat=ios) pv(i), pv(i + 1)
|
|
|
|
if (ios /= 0) then
|
|
|
|
print *, "I/O error: ", ios
|
|
|
|
exit
|
|
|
|
end if
|
|
|
|
end do
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! print coefficients
|
|
|
|
!do i=1,numberOfCoefficients,2
|
|
|
|
! print *," a=",pv(i)," b=",pv(i+1)
|
|
|
|
!end do
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
close (unit=1)
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_open_file(outfile, 'out.pv.grib1', 'w')
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-02-14 18:14:39 +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_samples(igrib, "reduced_gg_sfc_grib1")
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! set levtype to ml (model level)
|
|
|
|
call codes_set(igrib, 'typeOfLevel', 'hybrid')
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! set level
|
|
|
|
call codes_set(igrib, 'level', 2)
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! set PVPresent as an integer
|
|
|
|
call codes_set(igrib, 'PVPresent', 1)
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_set(igrib, 'pv', pv)
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
! write modified message to a file
|
|
|
|
call codes_write(igrib, outfile)
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2021-11-30 13:01:07 +00:00
|
|
|
! free memory
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_release(igrib)
|
|
|
|
deallocate (pv)
|
2013-03-25 12:04:10 +00:00
|
|
|
|
2021-02-14 18:14:39 +00:00
|
|
|
call codes_close_file(outfile)
|
2017-02-06 16:45:30 +00:00
|
|
|
|
2015-03-17 17:39:55 +00:00
|
|
|
end program grib_set_pv
|