File: grib_get_keys.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 (110 lines) | stat: -rw-r--r-- 4,030 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
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
! (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.
!
!
!  Description: how to get values using keys from GRIB messages
!
!
program grib_get_keys
  use eccodes
  implicit none

  integer                            ::  ifile
  integer                            ::  iret
  integer                            ::  igrib
  real                               ::  latitudeOfFirstPointInDegrees
  real                               ::  longitudeOfFirstPointInDegrees
  real                               ::  latitudeOfLastPointInDegrees
  real                               ::  longitudeOfLastPointInDegrees
  integer                            ::  numberOfPointsAlongAParallel
  integer                            ::  numberOfPointsAlongAMeridian
  real, dimension(:), allocatable    ::  values
  integer                            ::  numberOfValues
  real                               ::  average,min_val, max_val
  integer                            ::  is_missing
  character(len=10)                  ::  open_mode='r'

  call codes_open_file(ifile, &
       '../../data/reduced_latlon_surface.grib1', open_mode)

  ! Loop on all the messages in a file.

  ! 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(ifile,igrib, iret)

  LOOP: DO WHILE (iret /= CODES_END_OF_FILE)

    ! Check if the value of the key is MISSING
    is_missing=0;
    call codes_is_missing(igrib,'Ni',is_missing);
    if ( is_missing /= 1 ) then
        ! Key value is not missing so get as an integer
        call codes_get(igrib,'Ni',numberOfPointsAlongAParallel)
        write(*,*) 'numberOfPointsAlongAParallel=', &
             numberOfPointsAlongAParallel
    else
        write(*,*) 'numberOfPointsAlongAParallel is missing'
    endif

    ! Get as an integer
    call codes_get(igrib,'Nj',numberOfPointsAlongAMeridian)
    write(*,*) 'numberOfPointsAlongAMeridian=', &
         numberOfPointsAlongAMeridian

    ! Get as a real
    call codes_get(igrib, 'latitudeOfFirstGridPointInDegrees', &
          latitudeOfFirstPointInDegrees)
    write(*,*) 'latitudeOfFirstGridPointInDegrees=', &
          latitudeOfFirstPointInDegrees

    ! Get as a real
    call codes_get(igrib, 'longitudeOfFirstGridPointInDegrees', &
          longitudeOfFirstPointInDegrees)
    write(*,*) 'longitudeOfFirstGridPointInDegrees=', &
          longitudeOfFirstPointInDegrees

    ! Get as a real
    call codes_get(igrib, 'latitudeOfLastGridPointInDegrees', &
          latitudeOfLastPointInDegrees)
    write(*,*) 'latitudeOfLastGridPointInDegrees=', &
          latitudeOfLastPointInDegrees

    ! Get as a real
    call codes_get(igrib, 'longitudeOfLastGridPointInDegrees', &
          longitudeOfLastPointInDegrees)
    write(*,*) 'longitudeOfLastGridPointInDegrees=', &
          longitudeOfLastPointInDegrees

    ! Get the size of the values array
    call codes_get_size(igrib,'values',numberOfValues)
    write(*,*) 'numberOfValues=',numberOfValues

    allocate(values(numberOfValues), stat=iret)
    ! Get data values
    call codes_get(igrib,'values',values)
    call codes_get(igrib,'min',min_val) ! can also be obtained through minval(values)
    call codes_get(igrib,'max',max_val) ! can also be obtained through maxval(values)
    call codes_get(igrib,'average',average) ! can also be obtained through maxval(values)

    deallocate(values)

    write(*,*)'There are ',numberOfValues, &
          ' average is ',average, &
          ' min is ',  min_val, &
          ' max is ',  max_val

    call codes_release(igrib)

    call codes_grib_new_from_file(ifile,igrib, iret)

  end do LOOP

  call codes_close_file(ifile)

end program grib_get_keys