eccodes/html/precision_8f90-example.html

120 lines
7.2 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: precision.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>precision.f90</h1>How to control precision when coding a grib field.<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 !
<a name="l00010"></a>00010 ! Description: how to control decimal precision when packing fields.
<a name="l00011"></a>00011 !
<a name="l00012"></a>00012 !
<a name="l00013"></a>00013 ! Author: Enrico Fucile
<a name="l00014"></a>00014 !
<a name="l00015"></a>00015 !
<a name="l00016"></a>00016 !
<a name="l00017"></a>00017 program precision
<a name="l00018"></a>00018 use grib_api
<a name="l00019"></a>00019 implicit none
<a name="l00020"></a>00020 integer(kind = 4) :: size
<a name="l00021"></a>00021 integer :: infile,outfile
<a name="l00022"></a>00022 integer :: igrib
<a name="l00023"></a>00023 real(kind = 8), dimension(:), allocatable :: values1
<a name="l00024"></a>00024 real(kind = 8), dimension(:), allocatable :: values2
<a name="l00025"></a>00025 real(kind = 8) :: maxa,a,maxv,minv,maxr,r
<a name="l00026"></a>00026 integer( kind = 4) :: decimalPrecision,bitsPerValue1,bitsPerValue2
<a name="l00027"></a>00027 integer :: i, iret
<a name="l00028"></a>00028
<a name="l00029"></a>00029 call grib_open_file(infile, &amp;
<a name="l00030"></a>00030 '../../data/regular_latlon_surface_constant.grib1','r')
<a name="l00031"></a>00031
<a name="l00032"></a>00032 call grib_open_file(outfile, &amp;
<a name="l00033"></a>00033 '../../data/regular_latlon_surface_prec.grib1','w')
<a name="l00034"></a>00034
<a name="l00035"></a>00035 ! a new grib message is loaded from file
<a name="l00036"></a>00036 ! igrib is the grib id to be used in subsequent calls
<a name="l00037"></a>00037 call grib_new_from_file(infile,igrib)
<a name="l00038"></a>00038
<a name="l00039"></a>00039 ! bitsPerValue before changing the packing parameters
<a name="l00040"></a>00040 call grib_get(igrib,'bitsPerValue',bitsPerValue1)
<a name="l00041"></a>00041
<a name="l00042"></a>00042 ! get the size of the values array
<a name="l00043"></a>00043 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,<span class="stringliteral">"values"</span>,size)
<a name="l00044"></a>00044
<a name="l00045"></a>00045 allocate(values1(size), stat=iret)
<a name="l00046"></a>00046 allocate(values2(size), stat=iret)
<a name="l00047"></a>00047 ! get data values before changing the packing parameters*/
<a name="l00048"></a>00048 call grib_get(igrib,<span class="stringliteral">"values"</span>,values1)
<a name="l00049"></a>00049
<a name="l00050"></a>00050 ! setting decimal precision=2 means that 2 decimal digits
<a name="l00051"></a>00051 ! are preserved when packing.
<a name="l00052"></a>00052 decimalPrecision=2
<a name="l00053"></a>00053 call grib_set(igrib,<span class="stringliteral">"changeDecimalPrecision"</span>, &amp;
<a name="l00054"></a>00054 decimalPrecision)
<a name="l00055"></a>00055
<a name="l00056"></a>00056 ! bitsPerValue after changing the packing parameters
<a name="l00057"></a>00057 call grib_get(igrib,<span class="stringliteral">"bitsPerValue"</span>,bitsPerValue2)
<a name="l00058"></a>00058
<a name="l00059"></a>00059 ! get data values after changing the packing parameters
<a name="l00060"></a>00060 call grib_get(igrib,<span class="stringliteral">"values"</span>,values2)
<a name="l00061"></a>00061
<a name="l00062"></a>00062 ! computing error
<a name="l00063"></a>00063 maxa=0
<a name="l00064"></a>00064 maxr=0
<a name="l00065"></a>00065 maxv=values2(1)
<a name="l00066"></a>00066 minv=maxv
<a name="l00067"></a>00067 do i=1,size
<a name="l00068"></a>00068 a=abs(values2(i)-values1(i))
<a name="l00069"></a>00069 if ( values2(i) .gt. maxv ) maxv=values2(i)
<a name="l00070"></a>00070 if ( values2(i) .lt. maxv ) minv=values2(i)
<a name="l00071"></a>00071 if ( values2(i) .ne. 0 ) then
<a name="l00072"></a>00072 r=abs((values2(i)-values1(i))/values2(i))
<a name="l00073"></a>00073 endif
<a name="l00074"></a>00074 if ( a .gt. maxa ) maxa=a
<a name="l00075"></a>00075 if ( r .gt. maxr ) maxr=r
<a name="l00076"></a>00076 enddo
<a name="l00077"></a>00077 write(*,*) <span class="stringliteral">"max absolute error = "</span>,maxa
<a name="l00078"></a>00078 write(*,*) <span class="stringliteral">"max relative error = "</span>,maxr
<a name="l00079"></a>00079 write(*,*) <span class="stringliteral">"min value = "</span>,minv
<a name="l00080"></a>00080 write(*,*) <span class="stringliteral">"max value = "</span>,maxv
<a name="l00081"></a>00081
<a name="l00082"></a>00082 write(*,*) <span class="stringliteral">"old number of bits per value="</span>,bitsPerValue1
<a name="l00083"></a>00083 write(*,*) <span class="stringliteral">"new number of bits per value="</span>,bitsPerValue2
<a name="l00084"></a>00084
<a name="l00085"></a>00085 ! write modified message to a file
<a name="l00086"></a>00086 call grib_write(igrib,outfile)
<a name="l00087"></a>00087
<a name="l00088"></a>00088 call grib_release(igrib)
<a name="l00089"></a>00089
<a name="l00090"></a>00090 call grib_close_file(infile)
<a name="l00091"></a>00091
<a name="l00092"></a>00092 call grib_close_file(outfile)
<a name="l00093"></a>00093
<a name="l00094"></a>00094 deallocate(values1)
<a name="l00095"></a>00095 deallocate(values2)
<a name="l00096"></a>00096 end program precision
<a name="l00097"></a>00097
</pre></div> <hr size="1"><address style="text-align: right;"><small>Generated on Tue Sep 22 15:18:21 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>