00001 ! Copyright 2005-2016 ECMWF 00002 ! This software is licensed under the terms of the Apache Licence Version 2.0 00003 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. 00004 ! 00005 ! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by 00006 ! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction. 00007 ! 00008 ! 00009 ! Description: how to get values using keys. 00010 ! 00011 ! Author: Enrico Fucile 00012 ! 00013 ! 00014 program get 00015 use grib_api 00016 implicit none 00017 00018 integer :: ifile 00019 integer :: iret 00020 integer :: igrib 00021 real :: latitudeOfFirstPointInDegrees 00022 real :: longitudeOfFirstPointInDegrees 00023 real :: latitudeOfLastPointInDegrees 00024 real :: longitudeOfLastPointInDegrees 00025 integer :: numberOfPointsAlongAParallel 00026 integer :: numberOfPointsAlongAMeridian 00027 real, dimension(:), allocatable :: values 00028 integer :: numberOfValues 00029 real :: average,min_val, max_val 00030 integer :: is_missing 00031 00032 call grib_open_file(ifile, & 00033 '../../data/reduced_latlon_surface.grib1','r') 00034 00035 ! Loop on all the messages in a file. 00036 00037 ! a new grib message is loaded from file 00038 ! igrib is the grib id to be used in subsequent calls 00039 call grib_new_from_file(ifile,igrib, iret) 00040 00041 LOOP: DO WHILE (iret /= GRIB_END_OF_FILE) 00042 00043 !check if the value of the key is MISSING 00044 is_missing=0; 00045 call grib_is_missing(igrib,'numberOfPointsAlongAParallel', & 00046 is_missing); 00047 if ( is_missing /= 1 ) then 00048 ! get as a integer 00049 call grib_get(igrib,'numberOfPointsAlongAParallel', & 00050 numberOfPointsAlongAParallel) 00051 write(*,*) 'numberOfPointsAlongAParallel=', & 00052 numberOfPointsAlongAParallel 00053 else 00054 write(*,*) 'numberOfPointsAlongAParallel is missing' 00055 endif 00056 ! get as a integer 00057 call grib_get(igrib,'numberOfPointsAlongAMeridian', & 00058 numberOfPointsAlongAMeridian) 00059 write(*,*) 'numberOfPointsAlongAMeridian=', & 00060 numberOfPointsAlongAMeridian 00061 00062 ! get as a real 00063 call grib_get(igrib, 'latitudeOfFirstGridPointInDegrees', & 00064 latitudeOfFirstPointInDegrees) 00065 write(*,*) 'latitudeOfFirstGridPointInDegrees=', & 00066 latitudeOfFirstPointInDegrees 00067 00068 ! get as a real 00069 call grib_get(igrib, 'longitudeOfFirstGridPointInDegrees', & 00070 longitudeOfFirstPointInDegrees) 00071 write(*,*) 'longitudeOfFirstGridPointInDegrees=', & 00072 longitudeOfFirstPointInDegrees 00073 00074 ! get as a real 00075 call grib_get(igrib, 'latitudeOfLastGridPointInDegrees', & 00076 latitudeOfLastPointInDegrees) 00077 write(*,*) 'latitudeOfLastGridPointInDegrees=', & 00078 latitudeOfLastPointInDegrees 00079 00080 ! get as a real 00081 call grib_get(igrib, 'longitudeOfLastGridPointInDegrees', & 00082 longitudeOfLastPointInDegrees) 00083 write(*,*) 'longitudeOfLastGridPointInDegrees=', & 00084 longitudeOfLastPointInDegrees 00085 00086 00087 ! get the size of the values array 00088 call grib_get_size(igrib,'values',numberOfValues) 00089 write(*,*) 'numberOfValues=',numberOfValues 00090 00091 allocate(values(numberOfValues), stat=iret) 00092 ! get data values 00093 call grib_get(igrib,'values',values) 00094 call grib_get(igrib,'min',min_val) ! can also be obtained through minval(values) 00095 call grib_get(igrib,'max',max_val) ! can also be obtained through maxval(values) 00096 call grib_get(igrib,'average',average) ! can also be obtained through maxval(values) 00097 00098 write(*,*)'There are ',numberOfValues, & 00099 ' average is ',average, & 00100 ' min is ', min_val, & 00101 ' max is ', max_val 00102 00103 call grib_release(igrib) 00104 00105 call grib_new_from_file(ifile,igrib, iret) 00106 00107 end do LOOP 00108 00109 call grib_close_file(ifile) 00110 00111 deallocate(values) 00112 end program get