File: bufr_read_temp.f90

package info (click to toggle)
eccodes 2.20.0-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 400,332 kB
  • sloc: ansic: 167,977; makefile: 21,348; sh: 10,719; f90: 5,927; python: 4,831; perl: 3,031; javascript: 1,427; yacc: 818; lex: 356; awk: 66
file content (79 lines) | stat: -rw-r--r-- 3,586 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
! (C) Copyright 2005- ECMWF.
!
! 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.
!
! 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.
!
!
! FORTRAN 90 Implementation: bufr_read_temp
!
! Description: how to read levels from TEMP BUFR messages.
!
! Please note that TEMP reports can be encoded in various ways in BUFR. Therefore the code
! below might not work directly for other types of TEMP messages than the one used in the
! example. It is advised to use bufr_dump first to understand the structure of these messages.
!
program bufr_read_temp
  use eccodes
  implicit none
  integer            :: ifile
  integer            :: iret
  integer            :: ibufr
  integer            :: i, count=0
  integer(kind=4),dimension(:), allocatable  :: timePeriod,extendedVerticalSoundingSignificance
  integer(kind=4)  :: blockNumber,stationNumber
  real(kind=8),dimension(:), allocatable :: pressure,airTemperature,dewpointTemperature
  real(kind=8),dimension(:), allocatable :: geopotentialHeight,latitudeDisplacement,longitudeDisplacement
  real(kind=8),dimension(:), allocatable :: windDirection,windSpeed

  call codes_open_file(ifile,'../../data/bufr/PraticaTemp.bufr','r')

  ! The first bufr message is loaded from file,
  ! ibufr is the bufr id to be used in subsequent calls
  call codes_bufr_new_from_file(ifile,ibufr,iret)
  do while (iret/=CODES_END_OF_FILE)
    write(*,*) 'message: ',count
    call codes_set(ibufr,'unpack',1)
    call codes_get(ibufr,'timePeriod',timePeriod)
    call codes_get(ibufr,'pressure',pressure)
    call codes_get(ibufr,'extendedVerticalSoundingSignificance',extendedVerticalSoundingSignificance)
    call codes_get(ibufr,'nonCoordinateGeopotentialHeight',geopotentialHeight)
    call codes_get(ibufr,'latitudeDisplacement',latitudeDisplacement)
    call codes_get(ibufr,'longitudeDisplacement',longitudeDisplacement)
    call codes_get(ibufr,'airTemperature',airTemperature)
    call codes_get(ibufr,'dewpointTemperature',dewpointTemperature)
    call codes_get(ibufr,'windDirection',windDirection)
    call codes_get(ibufr,'windSpeed',windSpeed)
    call codes_get(ibufr,'blockNumber',blockNumber)
    call codes_get(ibufr,'stationNumber',stationNumber)
    print *,'station',blockNumber,stationNumber
    print *,'timePeriod pressure geopotentialHeight latitudeDisplacement &
          &longitudeDisplacement airTemperature windDirection windSpeed significance'
    do i=1,size(windSpeed)
      write(*,'(I5,6X,F9.1,2X,F9.2,10X,F8.2,14X,F8.2,16X,F8.2,6X,F8.2,4X,F8.2,4X,I0)') timePeriod(i),pressure(i),&
          &geopotentialHeight(i),latitudeDisplacement(i),&
          &longitudeDisplacement(i),airTemperature(i),windDirection(i),windSpeed(i),extendedVerticalSoundingSignificance(i)
    enddo
    ! Free arrays
    deallocate(timePeriod)
    deallocate(pressure)
    deallocate(geopotentialHeight)
    deallocate(latitudeDisplacement)
    deallocate(longitudeDisplacement)
    deallocate(airTemperature)
    deallocate(dewpointTemperature)
    deallocate(windDirection)
    deallocate(windSpeed)
    deallocate(extendedVerticalSoundingSignificance)
    ! Release the bufr message
    call codes_release(ibufr)
    ! Load the next bufr message
    call codes_bufr_new_from_file(ifile,ibufr,iret)
    count=count+1
  end do
  ! Close file
  call codes_close_file(ifile)

end program bufr_read_temp