source: flexpart.git/flexpart_code/detectformat.f90 @ b398fb6

FPv9.3.1FPv9.3.1b_testingFPv9.3.2fp9.3.1-20161214-nc4grib2nc4_repair
Last change on this file since b398fb6 was 496c607, checked in by Don Morton <Don.Morton@…>, 8 years ago

Initial commit of FPv9.3.1

Currently, this is a clone of snapshot FPv9.3.0

  • Property mode set to 100644
File size: 3.9 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine detectformat(metdata_format)
23
24  !*****************************************************************************
25  !                                                                            *
26  !   This routine reads the 1st file with windfields to determine             *
27  !   the format.                                                              *
28  !                                                                            *
29  !     Authors: M. Harustak                                                   *
30  !                                                                            *
31  !     6 May 2015                                                             *
32  !                                                                            *
33  !*****************************************************************************
34  !                                                                            *
35  ! Variables:                                                                 *
36  ! fname                file name of file to check                            *
37  !                                                                            *
38  !*****************************************************************************
39
40  use par_mod
41  use com_mod
42  use grib_api
43
44
45  implicit none
46
47  character(len=255) :: filename
48  character(len=255) :: wfname1(maxwf)
49  integer :: metdata_format
50  integer :: gfileid, status, gribid, centrenum
51
52  if ( maxwf.le.0 ) then
53    print*,'No wind file available'
54    metdata_format = unknown_metdata
55    return
56  endif
57
58  filename = path(3)(1:length(3)) // trim(wfname(1))
59 
60  ! Try to open the grib file - if successful,
61  ! gfileid will be the file handle
62  call grib_open_file(gfileid, TRIM(filename), 'r', status)
63  if (status.ne. GRIB_SUCCESS) then
64    print *, 'Unable to open: ', filename
65    metdata_format = unknown_metdata
66    return
67  endif
68
69  ! Get the GRIB id for the first message in the GRIB file.  This message
70  ! (and all the others) should have the "centre" field in it.  "gribid"
71  ! will be the handle for this message
72  call grib_new_from_file(gfileid, gribid, status)
73
74  ! Now, from the "gribid" message, get the value of the "centre" field
75  call grib_get(gribid, 'centre', centrenum)
76
77  if (centrenum .EQ. 7) then
78    metdata_format = gfs_metdata
79  elseif (centrenum == 98) then
80    metdata_format = ecmwf_metdata
81  else
82    metdata_format = unknown_metdata
83  endif
84
85  call grib_close_file(gfileid)
86
87
88end subroutine detectformat
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG