2013-04-02 14:02:10 +00:00
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN">
< html > < head > < meta http-equiv = "Content-Type" content = "text/html;charset=UTF-8" >
< title > grib_api: get_data.f90< / title >
< link href = "doxygen.css" rel = "stylesheet" type = "text/css" >
< link href = "tabs.css" rel = "stylesheet" type = "text/css" >
< / head > < body >
<!-- Generated by Doxygen 1.5.3 -->
< div class = "tabs" >
< ul >
< li > < a href = "index.html" > < span > Main Page< / span > < / a > < / li >
< li > < a href = "modules.html" > < span > Modules< / span > < / a > < / li >
< li > < a href = "files.html" > < span > Files< / span > < / a > < / li >
< li > < a href = "pages.html" > < span > Related Pages< / span > < / a > < / li >
< li > < a href = "examples.html" > < span > Examples< / span > < / a > < / li >
< / ul >
< / div >
< h1 > get_data.f90< / h1 > How to get latitude/longitude/values.< p >
2015-01-05 15:45:46 +00:00
< div class = "fragment" > < pre class = "fragment" > < a name = "l00001" > < / a > 00001 ! Copyright 2005-2015 ECMWF
2013-04-02 14:02:10 +00:00
< a name = "l00002" > < / a > 00002 ! This software is licensed under the terms of the Apache Licence Version 2.0
< a name = "l00003" > < / a > 00003 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
< a name = "l00004" > < / a > 00004 !
< a name = "l00005" > < / a > 00005 ! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
< a name = "l00006" > < / a > 00006 ! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
< a name = "l00007" > < / a > 00007 !
< a name = "l00008" > < / a > 00008 !
< a name = "l00009" > < / a > 00009 ! Description: how to < span class = "keyword" > get< / span > lat/lon/values.
< a name = "l00010" > < / a > 00010 !
< a name = "l00011" > < / a > 00011 !
< a name = "l00012" > < / a > 00012 ! Author: Enrico Fucile
< a name = "l00013" > < / a > 00013 !
< a name = "l00014" > < / a > 00014 !
< a name = "l00015" > < / a > 00015 program get_data
< a name = "l00016" > < / a > 00016 use grib_api
< a name = "l00017" > < / a > 00017 implicit none
< a name = "l00018" > < / a > 00018 integer :: ifile
< a name = "l00019" > < / a > 00019 integer :: iret,i
< a name = "l00020" > < / a > 00020 real(kind=8),dimension(:),allocatable :: lats,lons,values
< a name = "l00021" > < / a > 00021 integer(4) :: numberOfPoints
< a name = "l00022" > < / a > 00022 real(8) :: missingValue=9999
< a name = "l00023" > < / a > 00023 integer :: count=0
< a name = "l00024" > < / a > 00024 character(len=256) :: filename
< a name = "l00025" > < / a > 00025
< a name = "l00026" > < / a > 00026 ! Message identifier.
< a name = "l00027" > < / a > 00027 integer :: igrib
< a name = "l00028" > < / a > 00028
< a name = "l00029" > < / a > 00029 ifile=5
< a name = "l00030" > < / a > 00030
< a name = "l00031" > < / a > 00031 call grib_open_file(ifile, &
< a name = "l00032" > < / a > 00032 '../../data/reduced_latlon_surface.grib1','R')
< a name = "l00033" > < / a > 00033
< a name = "l00034" > < / a > 00034 ! Loop on all the messages in a file.
< a name = "l00035" > < / a > 00035
< a name = "l00036" > < / a > 00036 call grib_new_from_file(ifile,igrib,iret)
< a name = "l00037" > < / a > 00037
< a name = "l00038" > < / a > 00038 do while (iret/=< a name = "a0" > < / a > < a class = "code" href = "grib__api_8h.html#3bd3d72fe8bc116ca08c2d4b99203768" title = "End of ressource reached." > GRIB_END_OF_FILE< / a > )
< a name = "l00039" > < / a > 00039 count=count+1
< a name = "l00040" > < / a > 00040 print *, < span class = "stringliteral" > "===== Message #"< / span > ,count
< a name = "l00041" > < / a > 00041 call grib_get(igrib,'numberOfPoints',numberOfPoints)
< a name = "l00042" > < / a > 00042 call grib_set(igrib,'missingValue',missingValue)
< a name = "l00043" > < / a > 00043
< a name = "l00044" > < / a > 00044 allocate(lats(numberOfPoints))
< a name = "l00045" > < / a > 00045 allocate(lons(numberOfPoints))
< a name = "l00046" > < / a > 00046 allocate(values(numberOfPoints))
< a name = "l00047" > < / a > 00047
< a name = "l00048" > < / a > 00048 call grib_get_data(igrib,lats,lons,values)
< a name = "l00049" > < / a > 00049
< a name = "l00050" > < / a > 00050 do i=1,numberOfPoints
< a name = "l00051" > < / a > 00051 if (values(i) /= missingValue) then
< a name = "l00052" > < / a > 00052 print *, lats(i),lons(i),values(i)
< a name = "l00053" > < / a > 00053 end if
< a name = "l00054" > < / a > 00054 enddo
< a name = "l00055" > < / a > 00055
< a name = "l00056" > < / a > 00056 deallocate(lats)
< a name = "l00057" > < / a > 00057 deallocate(lons)
< a name = "l00058" > < / a > 00058 deallocate(values)
< a name = "l00059" > < / a > 00059
< a name = "l00060" > < / a > 00060 call grib_release(igrib)
< a name = "l00061" > < / a > 00061 call grib_new_from_file(ifile,igrib, iret)
< a name = "l00062" > < / a > 00062
< a name = "l00063" > < / a > 00063 end do
< a name = "l00064" > < / a > 00064
< a name = "l00065" > < / a > 00065
< a name = "l00066" > < / a > 00066 call grib_close_file(ifile)
< a name = "l00067" > < / a > 00067
< a name = "l00068" > < / a > 00068 end program
< / pre > < / div > < hr size = "1" > < address style = "text-align: right;" > < small > Generated on Tue Sep 22 15:18:21 2009 for grib_api by
< a href = "http://www.doxygen.org/index.html" >
< img src = "doxygen.png" alt = "doxygen" align = "middle" border = "0" > < / a > 1.5.3 < / small > < / address >
< / body >
< / html >