mirror of https://github.com/ecmwf/eccodes.git
104 lines
8.4 KiB
HTML
104 lines
8.4 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: iterator.c</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>iterator.c</h1>iterator.c How to use an iterator on latitude, longitude, values.<p>
|
|
<div class="fragment"><pre class="fragment"><a name="l00001"></a>00001
|
|
<a name="l00010"></a>00010 <span class="comment">/*</span>
|
|
<a name="l00011"></a>00011 <span class="comment"> * C Implementation: iterator</span>
|
|
<a name="l00012"></a>00012 <span class="comment"> *</span>
|
|
<a name="l00013"></a>00013 <span class="comment"> * Description: how to use an iterator on lat/lon/values.</span>
|
|
<a name="l00014"></a>00014 <span class="comment"> *</span>
|
|
<a name="l00015"></a>00015 <span class="comment"> *</span>
|
|
<a name="l00016"></a>00016 <span class="comment"> * Author: Enrico Fucile</span>
|
|
<a name="l00017"></a>00017 <span class="comment"> *</span>
|
|
<a name="l00018"></a>00018 <span class="comment"> *</span>
|
|
<a name="l00019"></a>00019 <span class="comment"> */</span>
|
|
<a name="l00020"></a>00020
|
|
<a name="l00021"></a>00021 <span class="preprocessor">#include <stdio.h></span>
|
|
<a name="l00022"></a>00022 <span class="preprocessor">#include <stdlib.h></span>
|
|
<a name="l00023"></a>00023 <span class="preprocessor">#include <string.h></span>
|
|
<a name="l00024"></a>00024
|
|
<a name="l00025"></a>00025 <span class="preprocessor">#include "<a class="code" href="grib__api_8h.html" title="Copyright 2005-2016 ECMWF.">grib_api.h</a>"</span>
|
|
<a name="l00026"></a>00026
|
|
<a name="l00027"></a>00027 <span class="keywordtype">void</span> usage(<span class="keywordtype">char</span>* prog) {
|
|
<a name="l00028"></a>00028 printf(<span class="stringliteral">"Usage: %s grib_file\n"</span>,prog);
|
|
<a name="l00029"></a>00029 exit(1);
|
|
<a name="l00030"></a>00030 }
|
|
<a name="l00031"></a>00031
|
|
<a name="l00032"></a>00032 <span class="keywordtype">int</span> main(<span class="keywordtype">int</span> argc, <span class="keywordtype">char</span>** argv) {
|
|
<a name="l00033"></a>00033 FILE* in = NULL;
|
|
<a name="l00034"></a>00034 <span class="keywordtype">int</span> err = 0;
|
|
<a name="l00035"></a>00035 <span class="keywordtype">double</span> lat,lon,value,missingValue=0;
|
|
<a name="l00036"></a>00036 <span class="keywordtype">int</span> n=0;
|
|
<a name="l00037"></a>00037 <span class="keywordtype">char</span>* filename = NULL;
|
|
<a name="l00038"></a>00038
|
|
<a name="l00039"></a>00039 <span class="comment">/* Message handle. Required in all the grib_api calls acting on a message.*/</span>
|
|
<a name="l00040"></a>00040 <a name="a0"></a><a class="code" href="group__grib__handle.html#g309a5ee24f4c730646d3f80ad0ef5f1b">grib_handle</a> *h = NULL;
|
|
<a name="l00041"></a>00041 <span class="comment">/* Iterator on lat/lon/values.*/</span>
|
|
<a name="l00042"></a>00042 <a name="a1"></a><a class="code" href="grib__api_8h.html#8f20a42a04122a55dec71774c70a51c5">grib_iterator</a>* iter=NULL;
|
|
<a name="l00043"></a>00043
|
|
<a name="l00044"></a>00044 <span class="keywordflow">if</span> (argc != 2) usage(argv[0]);
|
|
<a name="l00045"></a>00045
|
|
<a name="l00046"></a>00046 filename=strdup(argv[1]);
|
|
<a name="l00047"></a>00047
|
|
<a name="l00048"></a>00048 in = fopen(filename,<span class="stringliteral">"r"</span>);
|
|
<a name="l00049"></a>00049 <span class="keywordflow">if</span>(!in) {
|
|
<a name="l00050"></a>00050 printf(<span class="stringliteral">"ERROR: unable to open file %s\n"</span>,filename);
|
|
<a name="l00051"></a>00051 <span class="keywordflow">return</span> 1;
|
|
<a name="l00052"></a>00052 }
|
|
<a name="l00053"></a>00053
|
|
<a name="l00054"></a>00054 <span class="comment">/* Loop on all the messages in a file.*/</span>
|
|
<a name="l00055"></a>00055 <span class="keywordflow">while</span> ((h = <a name="a2"></a><a class="code" href="group__grib__handle.html#g5e24f8499aa7e4178ccc25a5de3145c5" title="Create a handle from a file resource.">grib_handle_new_from_file</a>(0,in,&err)) != NULL ) {
|
|
<a name="l00056"></a>00056 <span class="comment">/* Check of errors after reading a message. */</span>
|
|
<a name="l00057"></a>00057 <span class="keywordflow">if</span> (err != <a name="a3"></a><a class="code" href="grib__api_8h.html#5ec59f24fc07a0e9d05768e908b9eb41" title="No error.">GRIB_SUCCESS</a>) GRIB_CHECK(err,0);
|
|
<a name="l00058"></a>00058
|
|
<a name="l00059"></a>00059 <span class="comment">/* Get the double representing the missing value in the field. */</span>
|
|
<a name="l00060"></a>00060 GRIB_CHECK(<a name="a4"></a><a class="code" href="group__get__set.html#g5d9eeda38bf59ee3fa9ce3f92e65009e" title="Get a double value from a key, if several keys of the same name are present, the...">grib_get_double</a>(h,<span class="stringliteral">"missingValue"</span>,&missingValue),0);
|
|
<a name="l00061"></a>00061
|
|
<a name="l00062"></a>00062 <span class="comment">/* A new iterator on lat/lon/values is created from the message handle h. */</span>
|
|
<a name="l00063"></a>00063 iter=<a name="a5"></a><a class="code" href="group__iterators.html#gefb1f87110bdce732edc9154cf0e7d58" title="Create a new iterator from a handle, using current geometry and values.">grib_iterator_new</a>(h,0,&err);
|
|
<a name="l00064"></a>00064 <span class="keywordflow">if</span> (err != <a class="code" href="grib__api_8h.html#5ec59f24fc07a0e9d05768e908b9eb41" title="No error.">GRIB_SUCCESS</a>) GRIB_CHECK(err,0);
|
|
<a name="l00065"></a>00065
|
|
<a name="l00066"></a>00066 n = 0;
|
|
<a name="l00067"></a>00067 <span class="comment">/* Loop on all the lat/lon/values. */</span>
|
|
<a name="l00068"></a>00068 <span class="keywordflow">while</span>(<a name="a6"></a><a class="code" href="group__iterators.html#g4f73056dbfdda3de0060559b9b39ea34" title="Get the next value from an iterator.">grib_iterator_next</a>(iter,&lat,&lon,&value)) {
|
|
<a name="l00069"></a>00069 <span class="comment">/* You can now print lat and lon, */</span>
|
|
<a name="l00070"></a>00070 printf(<span class="stringliteral">"- %d - lat=%f lon=%f value="</span>,n,lat,lon);
|
|
<a name="l00071"></a>00071 <span class="comment">/* decide what to print if a missing value is found. */</span>
|
|
<a name="l00072"></a>00072 <span class="keywordflow">if</span> (value == missingValue ) printf(<span class="stringliteral">"missing\n"</span>);
|
|
<a name="l00073"></a>00073 <span class="comment">/* and print the value if is not missing. */</span>
|
|
<a name="l00074"></a>00074 <span class="keywordflow">else</span> printf(<span class="stringliteral">"%f\n"</span>,value);
|
|
<a name="l00075"></a>00075 n++;
|
|
<a name="l00076"></a>00076 }
|
|
<a name="l00077"></a>00077
|
|
<a name="l00078"></a>00078 <span class="comment">/* At the end the iterator is deleted to free memory. */</span>
|
|
<a name="l00079"></a>00079 <a name="a7"></a><a class="code" href="group__iterators.html#gd46ed73a16af56e6f3b46fe86ee8a759" title="Frees an iterator from memory.">grib_iterator_delete</a>(iter);
|
|
<a name="l00080"></a>00080
|
|
<a name="l00081"></a>00081 <span class="comment">/* At the end the grib_handle is deleted to free memory. */</span>
|
|
<a name="l00082"></a>00082 <a name="a8"></a><a class="code" href="group__grib__handle.html#g0e4b2585f22247c49b930c1579257677" title="Frees a handle, also frees the message if it is not a user message.">grib_handle_delete</a>(h);
|
|
<a name="l00083"></a>00083 }
|
|
<a name="l00084"></a>00084
|
|
<a name="l00085"></a>00085
|
|
<a name="l00086"></a>00086 fclose(in);
|
|
<a name="l00087"></a>00087
|
|
<a name="l00088"></a>00088 <span class="keywordflow">return</span> 0;
|
|
<a name="l00089"></a>00089 }
|
|
</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>
|