File: get_fortran.f90

package info (click to toggle)
eccodes 2.45.4-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 154,456 kB
  • sloc: cpp: 162,953; ansic: 26,308; sh: 21,742; f90: 6,854; perl: 6,361; python: 5,172; java: 2,226; javascript: 1,427; yacc: 854; fortran: 543; lex: 359; makefile: 278; xml: 183; awk: 66
file content (117 lines) | stat: -rw-r--r-- 4,317 bytes parent folder | download | duplicates (4)
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
111
112
113
114
115
116
117
! (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.
!
!
program get
   use eccodes
   implicit none

   integer                                         ::  ifile
   integer                                         ::  iret
   integer                                         ::  igrib
   integer                                         ::  i
   real(kind=8)                                    ::  latitudeOfFirstPointInDegrees
   real(kind=8)                                    ::  longitudeOfFirstPointInDegrees
   real(kind=8)                                    ::  latitudeOfLastPointInDegrees
   real(kind=8)                                    ::  longitudeOfLastPointInDegrees
   real(kind=8)                                    ::  jDirectionIncrementInDegrees
   real(kind=8)                                    ::  iDirectionIncrementInDegrees
   integer(kind=4)                               ::  numberOfPointsAlongAParallel
   integer(kind=4)                               ::  numberOfPointsAlongAMeridian
   real(kind=8), dimension(:), allocatable         ::  values
   integer(kind=4)                               ::  numberOfValues
   real(kind=8)                                    ::  average

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

   ! 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)

   ! get as a integer
   call codes_get(igrib, 'numberOfPointsAlongAParallel', &
                  numberOfPointsAlongAParallel)
   write (*, *) 'numberOfPointsAlongAParallel=', &
      numberOfPointsAlongAParallel

   ! get as a integer
   call codes_get(igrib, 'numberOfPointsAlongAMeridian', &
                  numberOfPointsAlongAMeridian)
   write (*, *) 'numberOfPointsAlongAMeridian=', &
      numberOfPointsAlongAMeridian

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

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

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

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

   ! get as a real8
   call codes_get(igrib, &
                  'jDirectionIncrementInDegrees', &
                  jDirectionIncrementInDegrees)
   write (*, *) 'jDirectionIncrementInDegrees=', &
      jDirectionIncrementInDegrees

   ! get as a real8
   call codes_get(igrib, &
                  'iDirectionIncrementInDegrees', &
                  iDirectionIncrementInDegrees)
   write (*, *) 'iDirectionIncrementInDegrees=', &
      iDirectionIncrementInDegrees

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

   allocate (values(2*numberOfValues), stat=iret)
   ! get data values
   print *, size(values)
   call codes_get(igrib, 'values', values)

   average = 0
   do i = 1, numberOfValues
      average = average + values(i);
   end do

   average = average/numberOfValues

   write (*, *) 'There are ', numberOfValues, &
      ' average is ', average

   call codes_release(igrib)

   call codes_close_file(ifile)

   deallocate (values)
end program get