mirror of https://github.com/ecmwf/eccodes.git
135 lines
8.7 KiB
HTML
135 lines
8.7 KiB
HTML
<!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.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.f90</h1>How to get values through the key names.<p>
|
|
<div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 ! Copyright 2005-2014 ECMWF
|
|
<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> values <span class="keyword">using</span> keys.
|
|
<a name="l00010"></a>00010 !
|
|
<a name="l00011"></a>00011 ! Author: Enrico Fucile
|
|
<a name="l00012"></a>00012 !
|
|
<a name="l00013"></a>00013 !
|
|
<a name="l00014"></a>00014 program <span class="keyword">get</span>
|
|
<a name="l00015"></a>00015 use grib_api
|
|
<a name="l00016"></a>00016 implicit none
|
|
<a name="l00017"></a>00017
|
|
<a name="l00018"></a>00018 integer :: ifile
|
|
<a name="l00019"></a>00019 integer :: iret
|
|
<a name="l00020"></a>00020 integer :: igrib
|
|
<a name="l00021"></a>00021 real :: latitudeOfFirstPointInDegrees
|
|
<a name="l00022"></a>00022 real :: longitudeOfFirstPointInDegrees
|
|
<a name="l00023"></a>00023 real :: latitudeOfLastPointInDegrees
|
|
<a name="l00024"></a>00024 real :: longitudeOfLastPointInDegrees
|
|
<a name="l00025"></a>00025 integer :: numberOfPointsAlongAParallel
|
|
<a name="l00026"></a>00026 integer :: numberOfPointsAlongAMeridian
|
|
<a name="l00027"></a>00027 real, dimension(:), allocatable :: values
|
|
<a name="l00028"></a>00028 integer :: numberOfValues
|
|
<a name="l00029"></a>00029 real :: average,min_val, max_val
|
|
<a name="l00030"></a>00030 integer :: is_missing
|
|
<a name="l00031"></a>00031
|
|
<a name="l00032"></a>00032 call grib_open_file(ifile, &
|
|
<a name="l00033"></a>00033 '../../data/reduced_latlon_surface.grib1<span class="charliteral">','</span>r')
|
|
<a name="l00034"></a>00034
|
|
<a name="l00035"></a>00035 ! Loop on all the messages in a file.
|
|
<a name="l00036"></a>00036
|
|
<a name="l00037"></a>00037 ! a <span class="keyword">new</span> grib message is loaded from file
|
|
<a name="l00038"></a>00038 ! igrib is the grib <span class="keywordtype">id</span> to be used in subsequent calls
|
|
<a name="l00039"></a>00039 call grib_new_from_file(ifile,igrib, iret)
|
|
<a name="l00040"></a>00040
|
|
<a name="l00041"></a>00041 LOOP: 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="l00042"></a>00042
|
|
<a name="l00043"></a>00043 !check if the value of the key is MISSING
|
|
<a name="l00044"></a>00044 is_missing=0;
|
|
<a name="l00045"></a>00045 call grib_is_missing(igrib,'numberOfPointsAlongAParallel', &
|
|
<a name="l00046"></a>00046 is_missing);
|
|
<a name="l00047"></a>00047 if ( is_missing /= 1 ) then
|
|
<a name="l00048"></a>00048 ! get as a integer
|
|
<a name="l00049"></a>00049 call grib_get(igrib,'numberOfPointsAlongAParallel', &
|
|
<a name="l00050"></a>00050 numberOfPointsAlongAParallel)
|
|
<a name="l00051"></a>00051 write(*,*) 'numberOfPointsAlongAParallel=', &
|
|
<a name="l00052"></a>00052 numberOfPointsAlongAParallel
|
|
<a name="l00053"></a>00053 else
|
|
<a name="l00054"></a>00054 write(*,*) 'numberOfPointsAlongAParallel is missing'
|
|
<a name="l00055"></a>00055 endif
|
|
<a name="l00056"></a>00056 ! get as a integer
|
|
<a name="l00057"></a>00057 call grib_get(igrib,'numberOfPointsAlongAMeridian', &
|
|
<a name="l00058"></a>00058 numberOfPointsAlongAMeridian)
|
|
<a name="l00059"></a>00059 write(*,*) 'numberOfPointsAlongAMeridian=', &
|
|
<a name="l00060"></a>00060 numberOfPointsAlongAMeridian
|
|
<a name="l00061"></a>00061
|
|
<a name="l00062"></a>00062 ! get as a real
|
|
<a name="l00063"></a>00063 call grib_get(igrib, 'latitudeOfFirstGridPointInDegrees', &
|
|
<a name="l00064"></a>00064 latitudeOfFirstPointInDegrees)
|
|
<a name="l00065"></a>00065 write(*,*) 'latitudeOfFirstGridPointInDegrees=', &
|
|
<a name="l00066"></a>00066 latitudeOfFirstPointInDegrees
|
|
<a name="l00067"></a>00067
|
|
<a name="l00068"></a>00068 ! get as a real
|
|
<a name="l00069"></a>00069 call grib_get(igrib, 'longitudeOfFirstGridPointInDegrees', &
|
|
<a name="l00070"></a>00070 longitudeOfFirstPointInDegrees)
|
|
<a name="l00071"></a>00071 write(*,*) 'longitudeOfFirstGridPointInDegrees=', &
|
|
<a name="l00072"></a>00072 longitudeOfFirstPointInDegrees
|
|
<a name="l00073"></a>00073
|
|
<a name="l00074"></a>00074 ! get as a real
|
|
<a name="l00075"></a>00075 call grib_get(igrib, 'latitudeOfLastGridPointInDegrees', &
|
|
<a name="l00076"></a>00076 latitudeOfLastPointInDegrees)
|
|
<a name="l00077"></a>00077 write(*,*) 'latitudeOfLastGridPointInDegrees=', &
|
|
<a name="l00078"></a>00078 latitudeOfLastPointInDegrees
|
|
<a name="l00079"></a>00079
|
|
<a name="l00080"></a>00080 ! get as a real
|
|
<a name="l00081"></a>00081 call grib_get(igrib, 'longitudeOfLastGridPointInDegrees', &
|
|
<a name="l00082"></a>00082 longitudeOfLastPointInDegrees)
|
|
<a name="l00083"></a>00083 write(*,*) 'longitudeOfLastGridPointInDegrees=', &
|
|
<a name="l00084"></a>00084 longitudeOfLastPointInDegrees
|
|
<a name="l00085"></a>00085
|
|
<a name="l00086"></a>00086
|
|
<a name="l00087"></a>00087 ! get the size of the values array
|
|
<a name="l00088"></a>00088 call <a name="a1"></a><a class="code" href="group__get__set.html#g18b622ed86b24d5e5fcab70c309fc245" title="Get the number of coded value from a key, if several keys of the same name are present...">grib_get_size</a>(igrib,'values',numberOfValues)
|
|
<a name="l00089"></a>00089 write(*,*) 'numberOfValues=',numberOfValues
|
|
<a name="l00090"></a>00090
|
|
<a name="l00091"></a>00091 allocate(values(numberOfValues), stat=iret)
|
|
<a name="l00092"></a>00092 ! get data values
|
|
<a name="l00093"></a>00093 call grib_get(igrib,'values',values)
|
|
<a name="l00094"></a>00094 call grib_get(igrib,'min',min_val) ! can also be obtained through minval(values)
|
|
<a name="l00095"></a>00095 call grib_get(igrib,'max',max_val) ! can also be obtained through maxval(values)
|
|
<a name="l00096"></a>00096 call grib_get(igrib,'average',average) ! can also be obtained through maxval(values)
|
|
<a name="l00097"></a>00097
|
|
<a name="l00098"></a>00098 write(*,*)'There are ',numberOfValues, &
|
|
<a name="l00099"></a>00099 ' average is ',average, &
|
|
<a name="l00100"></a>00100 ' min is ', min_val, &
|
|
<a name="l00101"></a>00101 ' max is ', max_val
|
|
<a name="l00102"></a>00102
|
|
<a name="l00103"></a>00103 call grib_release(igrib)
|
|
<a name="l00104"></a>00104
|
|
<a name="l00105"></a>00105 call grib_new_from_file(ifile,igrib, iret)
|
|
<a name="l00106"></a>00106
|
|
<a name="l00107"></a>00107 end do LOOP
|
|
<a name="l00108"></a>00108
|
|
<a name="l00109"></a>00109 call grib_close_file(ifile)
|
|
<a name="l00110"></a>00110
|
|
<a name="l00111"></a>00111 deallocate(values)
|
|
<a name="l00112"></a>00112 end program get
|
|
</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>
|