precision.f90

How to control precision when coding a grib field.

00001 ! Copyright 2005-2016 ECMWF
00002 ! This software is licensed under the terms of the Apache Licence Version 2.0
00003 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
00004 ! 
00005 ! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
00006 ! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
00007 !
00008 !
00009 !
00010 !  Description: how to control decimal precision when packing fields.
00011 !
00012 !
00013 !  Author: Enrico Fucile 
00014 !
00015 !
00016 !
00017 program precision
00018   use grib_api
00019   implicit none
00020   integer(kind = 4)                             :: size
00021   integer                                       :: infile,outfile
00022   integer                                       :: igrib
00023   real(kind = 8), dimension(:), allocatable     :: values1
00024   real(kind = 8), dimension(:), allocatable     :: values2
00025   real(kind = 8)                                ::  maxa,a,maxv,minv,maxr,r
00026   integer( kind = 4)                            :: decimalPrecision,bitsPerValue1,bitsPerValue2
00027   integer                                       :: i, iret
00028 
00029   call grib_open_file(infile, &
00030        '../../data/regular_latlon_surface_constant.grib1','r')
00031 
00032   call grib_open_file(outfile, &
00033        '../../data/regular_latlon_surface_prec.grib1','w')
00034 
00035   !     a new grib message is loaded from file
00036   !     igrib is the grib id to be used in subsequent calls
00037   call grib_new_from_file(infile,igrib)
00038 
00039   !     bitsPerValue before changing the packing parameters
00040   call grib_get(igrib,'bitsPerValue',bitsPerValue1)
00041 
00042   !     get the size of the values array
00043   call grib_get_size(igrib,"values",size)
00044 
00045   allocate(values1(size), stat=iret)
00046   allocate(values2(size), stat=iret)
00047   !     get data values before changing the packing parameters*/
00048   call grib_get(igrib,"values",values1)
00049 
00050   !     setting decimal precision=2 means that 2 decimal digits
00051   !     are preserved when packing.
00052   decimalPrecision=2
00053   call grib_set(igrib,"changeDecimalPrecision", &
00054        decimalPrecision)
00055 
00056   !     bitsPerValue after changing the packing parameters
00057   call grib_get(igrib,"bitsPerValue",bitsPerValue2)
00058 
00059   !     get data values after changing the packing parameters
00060   call grib_get(igrib,"values",values2)
00061 
00062   !     computing error
00063   maxa=0
00064   maxr=0
00065   maxv=values2(1)
00066   minv=maxv
00067   do i=1,size
00068      a=abs(values2(i)-values1(i))
00069      if ( values2(i) .gt. maxv ) maxv=values2(i)
00070      if ( values2(i) .lt. maxv ) minv=values2(i)
00071      if ( values2(i) .ne. 0 ) then
00072         r=abs((values2(i)-values1(i))/values2(i))
00073      endif
00074      if ( a .gt. maxa ) maxa=a
00075      if ( r .gt. maxr ) maxr=r
00076   enddo
00077   write(*,*) "max absolute error = ",maxa
00078   write(*,*) "max relative error = ",maxr
00079   write(*,*) "min value = ",minv
00080   write(*,*) "max value = ",maxv
00081 
00082   write(*,*) "old number of bits per value=",bitsPerValue1
00083   write(*,*) "new number of bits per value=",bitsPerValue2
00084 
00085   !     write modified message to a file
00086   call grib_write(igrib,outfile)
00087 
00088   call grib_release(igrib)
00089 
00090   call grib_close_file(infile)
00091 
00092   call grib_close_file(outfile)
00093 
00094   deallocate(values1)
00095   deallocate(values2)
00096 end program precision
00097 

Generated on Tue Sep 22 15:18:21 2009 for grib_api by  doxygen 1.5.3