eccodes/html/nearest_8f90-example.html

105 lines
5.9 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: nearest.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>nearest.f90</h1>How to find the nearest grid points.<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 use grib_find_nearest and grib_get_element
<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 !
<a name="l00016"></a>00016 program find
<a name="l00017"></a>00017 use grib_api
<a name="l00018"></a>00018 implicit none
<a name="l00019"></a>00019 integer :: npoints
<a name="l00020"></a>00020 integer :: infile
<a name="l00021"></a>00021 integer :: igrib, ios, i
<a name="l00022"></a>00022 real(8), dimension(:), allocatable :: lats, lons
<a name="l00023"></a>00023 real(8), dimension(:), allocatable :: nearest_lats, nearest_lons
<a name="l00024"></a>00024 real(8), dimension(:), allocatable :: distances, values, lsm_values
<a name="l00025"></a>00025 integer(kind=kindOfInt), dimension(:), allocatable :: indexes
<a name="l00026"></a>00026 real(kind=8) :: value
<a name="l00027"></a>00027
<a name="l00028"></a>00028 ! initialization
<a name="l00029"></a>00029 open( unit=1, file=<span class="stringliteral">"../../data/list_points"</span>,form=<span class="stringliteral">"formatted"</span>,action=<span class="stringliteral">"read"</span>)
<a name="l00030"></a>00030 read(unit=1,fmt=*) npoints
<a name="l00031"></a>00031 allocate(lats(npoints))
<a name="l00032"></a>00032 allocate(lons(npoints))
<a name="l00033"></a>00033 allocate(nearest_lats(npoints))
<a name="l00034"></a>00034 allocate(nearest_lons(npoints))
<a name="l00035"></a>00035 allocate(distances(npoints))
<a name="l00036"></a>00036 allocate(lsm_values(npoints))
<a name="l00037"></a>00037 allocate(values(npoints))
<a name="l00038"></a>00038 allocate(indexes(npoints))
<a name="l00039"></a>00039 do i=1,npoints
<a name="l00040"></a>00040 read(unit=1,fmt=*, iostat=ios) lats(i), lons(i)
<a name="l00041"></a>00041 if (ios /= 0) then
<a name="l00042"></a>00042 npoints = i - 1
<a name="l00043"></a>00043 exit
<a name="l00044"></a>00044 end if
<a name="l00045"></a>00045 end do
<a name="l00046"></a>00046 close(unit=1)
<a name="l00047"></a>00047 call grib_open_file(infile, &amp;
<a name="l00048"></a>00048 '../../data/reduced_gaussian_lsm.grib1','r')
<a name="l00049"></a>00049
<a name="l00050"></a>00050 ! a new grib message is loaded from file
<a name="l00051"></a>00051 ! igrib is the grib id to be used in subsequent calls
<a name="l00052"></a>00052 call grib_new_from_file(infile,igrib)
<a name="l00053"></a>00053
<a name="l00054"></a>00054
<a name="l00055"></a>00055 call grib_find_nearest(igrib, .true., lats, lons, nearest_lats, nearest_lons,lsm_values, distances, indexes)
<a name="l00056"></a>00056 call grib_release(igrib)
<a name="l00057"></a>00057
<a name="l00058"></a>00058 call grib_close_file(infile)
<a name="l00059"></a>00059
<a name="l00060"></a>00060 ! will apply it to another GRIB
<a name="l00061"></a>00061 call grib_open_file(infile, &amp;
<a name="l00062"></a>00062 '../../data/reduced_gaussian_pressure_level.grib1','r')
<a name="l00063"></a>00063 call grib_new_from_file(infile,igrib)
<a name="l00064"></a>00064
<a name="l00065"></a>00065 call grib_get_element(igrib,<span class="stringliteral">"values"</span>, indexes, values)
<a name="l00066"></a>00066 call grib_release(igrib)
<a name="l00067"></a>00067 call grib_close_file(infile)
<a name="l00068"></a>00068
<a name="l00069"></a>00069 do i=1, npoints
<a name="l00070"></a>00070 print*,lats(i), lons(i), nearest_lats(i), nearest_lons(i), distances(i), lsm_values(i), values(i)
<a name="l00071"></a>00071 end do
<a name="l00072"></a>00072
<a name="l00073"></a>00073 deallocate(lats)
<a name="l00074"></a>00074 deallocate(lons)
<a name="l00075"></a>00075 deallocate(nearest_lats)
<a name="l00076"></a>00076 deallocate(nearest_lons)
<a name="l00077"></a>00077 deallocate(distances)
<a name="l00078"></a>00078 deallocate(lsm_values)
<a name="l00079"></a>00079 deallocate(values)
<a name="l00080"></a>00080 deallocate(indexes)
<a name="l00081"></a>00081
<a name="l00082"></a>00082 end program find
</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>