eccodes/html/clone_8f90-example.html

102 lines
5.6 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: clone.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>clone.f90</h1>How to clone a message.<p>
<div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 ! Copyright 2005-2015 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 create a <span class="keyword">new</span> GRIB message by cloning
<a name="l00010"></a>00010 ! an existing message.
<a name="l00011"></a>00011 !
<a name="l00012"></a>00012 !
<a name="l00013"></a>00013 ! Author: Anne Fouilloux
<a name="l00014"></a>00014 !
<a name="l00015"></a>00015 !
<a name="l00016"></a>00016 program clone
<a name="l00017"></a>00017 use grib_api
<a name="l00018"></a>00018 implicit none
<a name="l00019"></a>00019 integer :: err,i,iret
<a name="l00020"></a>00020 integer :: nx, ny
<a name="l00021"></a>00021 integer :: infile,outfile
<a name="l00022"></a>00022 integer :: igrib_in
<a name="l00023"></a>00023 integer :: igrib_out
<a name="l00024"></a>00024 character(len=2) :: step
<a name="l00025"></a>00025 double precision, dimension(:,:), allocatable :: field2D
<a name="l00026"></a>00026
<a name="l00027"></a>00027
<a name="l00028"></a>00028 call grib_open_file(infile,'../../data/constant_field.grib1','r')
<a name="l00029"></a>00029 call grib_open_file(outfile,'out.grib1','w')
<a name="l00030"></a>00030
<a name="l00031"></a>00031 ! a new grib message is loaded from file
<a name="l00032"></a>00032 ! igrib is the grib id to be used in subsequent calls
<a name="l00033"></a>00033 call grib_new_from_file(infile,igrib_in)
<a name="l00034"></a>00034
<a name="l00035"></a>00035 call grib_get(igrib_in,<span class="stringliteral">"numberOfPointsAlongAParallel"</span>, nx)
<a name="l00036"></a>00036
<a name="l00037"></a>00037 call grib_get(igrib_in,<span class="stringliteral">"numberOfPointsAlongAMeridian"</span>,ny)
<a name="l00038"></a>00038
<a name="l00039"></a>00039 allocate(field2D(nx,ny),stat=err)
<a name="l00040"></a>00040
<a name="l00041"></a>00041 if (err .ne. 0) then
<a name="l00042"></a>00042 print*, 'Failed to allocate ', nx*ny, ' values'
<a name="l00043"></a>00043 STOP
<a name="l00044"></a>00044 end if
<a name="l00045"></a>00045 ! clone the constant field to create 4 new GRIB messages
<a name="l00046"></a>00046 do i=0,18,6
<a name="l00047"></a>00047 call grib_clone(igrib_in, igrib_out)
<a name="l00048"></a>00048 write(step,'(i2)') i
<a name="l00049"></a>00049 ! Careful: stepRange is a string (could be 0-6, 12-24, etc.)
<a name="l00050"></a>00050 ! use adjustl to remove blank from the left.
<a name="l00051"></a>00051 call grib_set(igrib_out,'stepRange',adjustl(step))
<a name="l00052"></a>00052
<a name="l00053"></a>00053 call generate_field(field2D)
<a name="l00054"></a>00054
<a name="l00055"></a>00055 ! use pack to create 1D values
<a name="l00056"></a>00056 call grib_set(igrib_out,'values',pack(field2D, mask=.true.))
<a name="l00057"></a>00057
<a name="l00058"></a>00058 ! write cloned messages to a file
<a name="l00059"></a>00059 call grib_write(igrib_out,outfile)
<a name="l00060"></a>00060 call grib_release(igrib_out)
<a name="l00061"></a>00061 end do
<a name="l00062"></a>00062
<a name="l00063"></a>00063 call grib_release(igrib_in)
<a name="l00064"></a>00064
<a name="l00065"></a>00065 call grib_close_file(infile)
<a name="l00066"></a>00066
<a name="l00067"></a>00067 call grib_close_file(outfile)
<a name="l00068"></a>00068 deallocate(field2D)
<a name="l00069"></a>00069
<a name="l00070"></a>00070 contains
<a name="l00071"></a>00071 !======================================
<a name="l00072"></a>00072 subroutine generate_field(gfield2D)
<a name="l00073"></a>00073 double precision, dimension(:,:) :: gfield2D
<a name="l00074"></a>00074
<a name="l00075"></a>00075 call random_number(gfield2D)
<a name="l00076"></a>00076 end subroutine generate_field
<a name="l00077"></a>00077 !======================================
<a name="l00078"></a>00078
<a name="l00079"></a>00079 end program clone
</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>