File: grib_get_data.f90

package info (click to toggle)
eccodes 2.44.2-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 150,256 kB
  • sloc: cpp: 163,056; ansic: 26,308; sh: 21,602; f90: 6,854; perl: 6,363; python: 5,087; java: 2,226; javascript: 1,427; yacc: 854; fortran: 543; lex: 359; makefile: 274; xml: 183; awk: 66
file content (78 lines) | stat: -rw-r--r-- 2,484 bytes parent folder | download | duplicates (2)
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
! (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 lat/lon/values.
!
!
program get_data
   use eccodes
   implicit none
   integer            :: ifile
   integer            :: iret, i
   real(kind=8), dimension(:), allocatable     :: lats, lons, values
   integer, dimension(:), allocatable          :: bitmap
   integer(4)        :: numberOfPoints
   logical            :: is_missing_value
   integer            :: count1 = 0, count2 = 0, bitmapPresent = 0, bmp_len = 0
   integer            :: igrib  ! message identifier

   ifile = 5

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

   ! loop on all the messages in the file
   do while (.true.)
      call codes_grib_new_from_file(ifile, igrib, iret)
      if (iret == CODES_END_OF_FILE) exit

      count1 = count1 + 1
      print *, "===== Message #", count1
      call codes_get(igrib, 'numberOfPoints', numberOfPoints)
      call codes_get(igrib, 'bitmapPresent', bitmapPresent)

      allocate (lats(numberOfPoints))
      allocate (lons(numberOfPoints))
      allocate (values(numberOfPoints))
      if (bitmapPresent == 1) then
         ! get the bitmap
         call codes_get_size(igrib, 'bitmap', bmp_len)
         allocate (bitmap(bmp_len))
         call codes_get(igrib, 'bitmap', bitmap)
      end if

      call codes_grib_get_data(igrib, lats, lons, values)

      do i = 1, numberOfPoints
         ! consult bitmap to see if the i'th value is missing
         is_missing_value = .false.
         if (bitmapPresent == 1 .and. bitmap(i) == 0) then
            is_missing_value = .true.
         end if
         ! only print non-missing values
         if (.not. is_missing_value) then
            print *, lats(i), lons(i), values(i)
            count2 = count2 + 1
         end if
      end do
      print *, 'count of non-missing values=', count2
      if (count2 /= 214661) then
         call codes_check(-2, 'incorrect number of missing', '')
      end if

      deallocate (lats)
      deallocate (lons)
      deallocate (values)

      call codes_release(igrib)

   end do

   call codes_close_file(ifile)

end program