eccodes/html/set__bitmap_8f90-example.html

103 lines
5.8 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: set_bitmap.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&nbsp;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&nbsp;Pages</span></a></li>
<li><a href="examples.html"><span>Examples</span></a></li>
</ul>
</div>
<h1>set_bitmap.f90</h1>How to set and use a bitmap.<p>
<div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 ! Copyright 2005-2013 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">set</span> a bitmap in a grib message
<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 set_bitmap
<a name="l00016"></a>00016 use grib_api
<a name="l00017"></a>00017 implicit none
<a name="l00018"></a>00018 integer :: infile,outfile
<a name="l00019"></a>00019 integer :: igrib, iret
<a name="l00020"></a>00020 integer :: numberOfValues
<a name="l00021"></a>00021 real, dimension(:), allocatable :: values
<a name="l00022"></a>00022 real :: missingValue
<a name="l00023"></a>00023 logical :: grib1Example
<a name="l00024"></a>00024
<a name="l00025"></a>00025 grib1Example=.true.
<a name="l00026"></a>00026
<a name="l00027"></a>00027 <span class="keywordflow">if</span> (grib1Example) then
<a name="l00028"></a>00028 ! GRIB 1 example
<a name="l00029"></a>00029 call grib_open_file(infile,'../../data/regular_latlon_surface.grib1<span class="charliteral">','</span>r')
<a name="l00030"></a>00030 else
<a name="l00031"></a>00031 ! GRIB 2 example
<a name="l00032"></a>00032 call grib_open_file(infile,'../../data/regular_latlon_surface.grib2','r')
<a name="l00033"></a>00033 end if
<a name="l00034"></a>00034
<a name="l00035"></a>00035 call grib_open_file(outfile,'out.grib','w')
<a name="l00036"></a>00036
<a name="l00037"></a>00037 ! a new 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(infile,igrib)
<a name="l00040"></a>00040
<a name="l00041"></a>00041 ! The missingValue is not coded in the message.
<a name="l00042"></a>00042 ! It is a value we define as a placeholder for a missing value
<a name="l00043"></a>00043 ! in a point of the grid.
<a name="l00044"></a>00044 ! It should be choosen in a way that it cannot be confused
<a name="l00045"></a>00045 ! with a valid field value
<a name="l00046"></a>00046 missingValue=9999
<a name="l00047"></a>00047 call grib_set(igrib, 'missingValue',missingValue)
<a name="l00048"></a>00048 write(*,*) 'missingValue=',missingValue
<a name="l00049"></a>00049
<a name="l00050"></a>00050 ! get the size of the values array
<a name="l00051"></a>00051 call <a name="a0"></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="l00052"></a>00052 write(*,*) 'numberOfValues=',numberOfValues
<a name="l00053"></a>00053
<a name="l00054"></a>00054 allocate(values(numberOfValues), stat=iret)
<a name="l00055"></a>00055
<a name="l00056"></a>00056 ! get data values
<a name="l00057"></a>00057 call grib_get(igrib,'values',values)
<a name="l00058"></a>00058
<a name="l00059"></a>00059 ! enable bitmap
<a name="l00060"></a>00060 call grib_set(igrib,"bitmapPresent",1)
<a name="l00061"></a>00061
<a name="l00062"></a>00062 ! some values are missing
<a name="l00063"></a>00063 values(1:10) = missingValue
<a name="l00064"></a>00064
<a name="l00065"></a>00065 ! set the values (the bitmap will be automatically built)
<a name="l00066"></a>00066 call grib_set(igrib,'values', values)
<a name="l00067"></a>00067
<a name="l00068"></a>00068 ! write modified message to a file
<a name="l00069"></a>00069 call grib_write(igrib,outfile)
<a name="l00070"></a>00070
<a name="l00071"></a>00071 ! FREE MEMORY
<a name="l00072"></a>00072 call grib_release(igrib)
<a name="l00073"></a>00073
<a name="l00074"></a>00074 call grib_close_file(infile)
<a name="l00075"></a>00075
<a name="l00076"></a>00076 call grib_close_file(outfile)
<a name="l00077"></a>00077
<a name="l00078"></a>00078 deallocate(values)
<a name="l00079"></a>00079
<a name="l00080"></a>00080 end program set_bitmap
</pre></div> <hr size="1"><address style="text-align: right;"><small>Generated on Tue Sep 22 15:18:22 2009 for grib_api by&nbsp;
<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>