mirror of https://github.com/ecmwf/eccodes.git
243 lines
17 KiB
HTML
243 lines
17 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: grib_dump</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><a class="anchor" name="grib_dump">grib_dump</a></h1><h2><a class="anchor" name="DESCRIPTION">
|
|
DESCRIPTION</a></h2>
|
|
Dump the content of a grib file in different formats.<h2><a class="anchor" name="USAGE">
|
|
USAGE</a></h2>
|
|
grib_dump [options] grib_file grib_file ...<h2><a class="anchor" name="OPTIONS">
|
|
OPTIONS</a></h2>
|
|
-O <br>
|
|
Octet mode. WMO documentation style dump. <br>
|
|
<br>
|
|
-D <br>
|
|
Debug mode. <br>
|
|
<br>
|
|
-P key[:{s/d/l}],key[:{s/d/l}],... <br>
|
|
As -p adding the declared keys to the default list. <br>
|
|
<br>
|
|
-d <br>
|
|
Print all data values. Available only in C mode <br>
|
|
<br>
|
|
-C <br>
|
|
C code mode. A C code program generating the grib message is dumped. <br>
|
|
<br>
|
|
-t <br>
|
|
Print type information. <br>
|
|
<br>
|
|
-H <br>
|
|
Print octet content in hexadecimal format. <br>
|
|
<br>
|
|
-a <br>
|
|
Dump aliases. <br>
|
|
<br>
|
|
-w key[:{s/d/l}]{=/!=}value,key[:{s/d/l}]{=/!=}value,... <br>
|
|
Where clause. Grib messages are processed only if they match all the key/value constraints. A valid constraint is of type key=value or key!=value. For each key a string (key:s) or a double (key:d) or a long (key:l) type can be specified. Default type is string. <br>
|
|
<br>
|
|
-M <br>
|
|
Multi-grib support off. Turn off support for multiple fields in single grib message <br>
|
|
<br>
|
|
-7 <br>
|
|
Does not fail when the message has wrong length <br>
|
|
<br>
|
|
-V <br>
|
|
Version. <br>
|
|
<br>
|
|
-G <br>
|
|
GRIBEX compatibility mode. <br>
|
|
<br>
|
|
<h2><a class="anchor" name="grib_dump_examples">
|
|
grib_dump examples</a></h2>
|
|
<ol type=1>
|
|
<li>To dump in a WMO documentation style with hexadecimal octet values (-H)<br>
|
|
<div class="fragment"><pre class="fragment">
|
|
>grib_dump -H ../data/reduced_gaussian_model_level.grib1
|
|
</pre></div><br>
|
|
</li><li>To obtain all the key names available in a grib file.<br>
|
|
<div class="fragment"><pre class="fragment">
|
|
> grib_dump -D ../data/regular_latlon_surface.grib1
|
|
</pre></div><br>
|
|
</li><li>To obtain a C code example from a grib file.<br>
|
|
<div class="fragment"><pre class="fragment">>grib_dump -C ../data/regular_latlon_surface.grib1
|
|
<span class="preprocessor">#include <<a class="code" href="grib__api_8h.html" title="Copyright 2005-2015 ECMWF.">grib_api.h</a>></span>
|
|
|
|
<span class="comment">/* This code was generated automatically */</span>
|
|
|
|
|
|
<span class="keywordtype">int</span> main(<span class="keywordtype">int</span> argc,<span class="keyword">const</span> <span class="keywordtype">char</span>** argv)
|
|
{
|
|
<a class="code" href="group__grib__handle.html#g309a5ee24f4c730646d3f80ad0ef5f1b">grib_handle</a> *h = NULL;
|
|
<span class="keywordtype">size_t</span> size = 0;
|
|
<span class="keywordtype">double</span>* vdouble = NULL;
|
|
<span class="keywordtype">long</span>* vlong = NULL;
|
|
FILE* f = NULL;
|
|
<span class="keyword">const</span> <span class="keywordtype">char</span>* p = NULL;
|
|
<span class="keyword">const</span> <span class="keywordtype">void</span>* buffer = NULL;
|
|
|
|
<span class="keywordflow">if</span>(argc != 2) {
|
|
fprintf(stderr,<span class="stringliteral">"usage: %s out\n"</span>,argv[0]);
|
|
exit(1);
|
|
}
|
|
|
|
h = <a class="code" href="group__grib__handle.html#gadefac64c19fb5ff06cf805ad4af06ff" title="Create a handle from a message contained in a samples directory.">grib_handle_new_from_samples</a>(NULL,<span class="stringliteral">"GRIB1"</span>);
|
|
<span class="keywordflow">if</span>(!h) {
|
|
fprintf(stderr,<span class="stringliteral">"Cannot create grib handle\n"</span>);
|
|
exit(1);
|
|
}
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"editionNumber"</span>,1),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"table2Version"</span>,128),0);
|
|
|
|
<span class="comment">/* 98 = European Center for Medium-Range Weather Forecasts (grib1/0.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"centre"</span>,98),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"generatingProcessIdentifier"</span>,130),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"gridDefinition"</span>,255),0);
|
|
|
|
<span class="comment">/* 128 = 10000000</span>
|
|
<span class="comment"> (1=1) Section 2 included</span>
|
|
<span class="comment"> (2=0) Section 3 omited</span>
|
|
<span class="comment"> See grib1/1.table */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"section1Flags"</span>,128),0);
|
|
|
|
|
|
<span class="comment">/* 167 = 2 metre temperature (K) (grib1/2.98.128.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"indicatorOfParameter"</span>,167),0);
|
|
|
|
|
|
<span class="comment">/* 1 = Surface (of the Earth, which includes sea surface) (grib1/3.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"indicatorOfTypeOfLevel"</span>,1),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"level"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"yearOfCentury"</span>,8),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"month"</span>,2),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"day"</span>,6),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"hour"</span>,12),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"minute"</span>,0),0);
|
|
|
|
<span class="comment">/* 1 = Hour (grib1/4.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"unitOfTimeRange"</span>,1),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"P1"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"P2"</span>,0),0);
|
|
|
|
<span class="comment">/* 0 = Forecast product valid at reference time + P1 (P1>0) (grib1/5.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"timeRangeIndicator"</span>,0),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"numberIncludedInAverage"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"numberMissingFromAveragesOrAccumulations"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"centuryOfReferenceTimeOfData"</span>,21),0);
|
|
|
|
<span class="comment">/* 0 = Unknown code table entry (grib1/0.ecmf.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"subCentre"</span>,0),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"decimalScaleFactor"</span>,0),0);
|
|
|
|
<span class="comment">/* 1 = MARS labelling or ensemble forecast data (grib1/localDefinitionNumber.98.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"localDefinitionNumber"</span>,1),0);
|
|
|
|
|
|
<span class="comment">/* 1 = Operational archive (mars/class.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"marsClass"</span>,1),0);
|
|
|
|
|
|
<span class="comment">/* 2 = Analysis (mars/type.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"marsType"</span>,2),0);
|
|
|
|
|
|
<span class="comment">/* 1025 = Atmospheric model (mars/stream.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"marsStream"</span>,1025),0);
|
|
|
|
p = <span class="stringliteral">"0001"</span>;
|
|
size = strlen(p)+1;
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g78143cc64571b454b0aba14246e9a53a" title="Set a string value from a key.">grib_set_string</a>(h,<span class="stringliteral">"experimentVersionNumber"</span>,p,&size),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"perturbationNumber"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"numberOfForecastsInEnsemble"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"numberOfVerticalCoordinateValues"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"pvlLocation"</span>,255),0);
|
|
|
|
<span class="comment">/* 0 = Latitude/Longitude Grid (grib1/6.table) */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"dataRepresentationType"</span>,0),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"Ni"</span>,16),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"Nj"</span>,31),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"latitudeOfFirstGridPoint"</span>,60000),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"longitudeOfFirstGridPoint"</span>,0),0);
|
|
|
|
<span class="comment">/* 128 = 10000000</span>
|
|
<span class="comment"> (1=1) Direction increments given</span>
|
|
<span class="comment"> (2=0) Earth assumed spherical with radius = 6367.47 km</span>
|
|
<span class="comment"> (5=0) u and v components resolved relative to easterly and northerly directions</span>
|
|
<span class="comment"> See grib1/7.table */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"resolutionAndComponentFlags"</span>,128),0);
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"latitudeOfLastGridPoint"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"longitudeOfLastGridPoint"</span>,30000),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"iDirectionIncrement"</span>,2000),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"jDirectionIncrement"</span>,2000),0);
|
|
|
|
<span class="comment">/* 0 = 00000000</span>
|
|
<span class="comment"> (1=0) Points scan in +i direction</span>
|
|
<span class="comment"> (2=0) Points scan in -j direction</span>
|
|
<span class="comment"> (3=0) Adjacent points in i direction are consecutive </span>
|
|
<span class="comment"> See grib1/8.table */</span>
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"scanningMode"</span>,0),0);
|
|
|
|
|
|
<span class="comment">/* ITERATOR */</span>
|
|
|
|
|
|
<span class="comment">/* NEAREST */</span>
|
|
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"bitsPerValue"</span>,16),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"sphericalHarmonics"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"complexPacking"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"integerPointValues"</span>,0),0);
|
|
GRIB_CHECK(<a class="code" href="group__get__set.html#g94c33cfe90c3aa887fb8e14f0bd87fe2" title="Set a long value from a key.">grib_set_long</a>(h,<span class="stringliteral">"additionalFlagPresent"</span>,0),0);
|
|
|
|
<span class="comment">/* gribSection5 */</span>
|
|
|
|
<span class="comment">/* Save the message */</span>
|
|
|
|
f = fopen(argv[1],<span class="stringliteral">"w"</span>);
|
|
<span class="keywordflow">if</span>(!f) {
|
|
perror(argv[1]);
|
|
exit(1);
|
|
}
|
|
|
|
GRIB_CHECK(<a class="code" href="group__handling__coded__messages.html#g9d654bd4fc5f422c161edd0a140ea185" title="getting the message attached to a handle">grib_get_message</a>(h,&buffer,&size),0);
|
|
|
|
<span class="keywordflow">if</span>(fwrite(buffer,1,size,f) != size) {
|
|
perror(argv[1]);
|
|
exit(1);
|
|
}
|
|
|
|
<span class="keywordflow">if</span>(fclose(f)) {
|
|
perror(argv[1]);
|
|
exit(1);
|
|
}
|
|
|
|
<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);
|
|
<span class="keywordflow">return</span> 0;
|
|
}
|
|
</pre></div><br>
|
|
</li></ol>
|
|
<hr size="1"><address style="text-align: right;"><small>Generated on Tue Sep 22 15:18:22 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>
|