Changeset c0884a8 in flexpart.git for src/getfields.f90


Ignore:
Timestamp:
Jul 20, 2018, 2:40:39 PM (6 years ago)
Author:
pesei <petra seibert at univie ac at>
Branches:
univie
Children:
7ca2ef4
Parents:
f251e57
Message:

replace CTBTO code for checking type of GRIB

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/getfields.f90

    r6ecb30a rc0884a8  
    2020!**********************************************************************
    2121
    22 subroutine getfields(itime,nstop,metdata_format)
     22subroutine getfields(itime,nstop,id_centre)
    2323!                       i     o
    2424!*****************************************************************************
     
    4444!   Unified ECMWF and GFS builds                                             *
    4545!   Marian Harustak, 12.5.2017                                               *
    46 !     - Added passing of metdata_format as it was needed by called routines  *
     46!     - Added passing of id_centre as it was needed by called routines       *
     47!                                                                            *
     48!  Petra Seibert, 2018-06-26: simplified version met data format detection   *
     49!                                                                            *
    4750!*****************************************************************************
    4851!                                                                            *
     
    6265! tt(0:nxmax,0:nymax,nuvzmax,2)   temperature [K]                            *
    6366! ps(0:nxmax,0:nymax,2)           surface pressure [Pa]                      *
    64 ! metdata_format     format of metdata (ecmwf/gfs)                           *
     67! id_centre            format of metdata (ecmwf/gfs)                         *
    6568!                                                                            *
    6669! Constants:                                                                 *
     
    7275  use par_mod
    7376  use com_mod
    74   use class_gribfile
     77  use check_gribfile_mod
    7578
    7679  implicit none
    7780
    7881  integer :: indj,itime,nstop,memaux
    79   integer :: metdata_format
     82  integer :: id_centre
    8083
    8184  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
     
    132135    do indj=indmin,numbwf-1
    133136      if (ldirect*wftime(indj+1).gt.ldirect*itime) then
    134         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     137        if (id_centre.eq.icg_id_ecmwf) then
    135138          call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
    136139        else
     
    138141        end if
    139142        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
    140         call calcpar(memind(2),uuh,vvh,pvh,metdata_format)
    141         call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
    142         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     143        call calcpar(memind(2),uuh,vvh,pvh,id_centre)
     144        call calcpar_nests(memind(2),uuhn,vvhn,pvhn,id_centre)
     145        if (id_centre.eq.icg_id_ecmwf) then
    143146          call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
    144147        else
     
    167170           (ldirect*wftime(indj+1).gt.ldirect*itime)) then
    168171        memind(1)=1
    169         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     172        if (id_centre.eq.icg_id_ecmwf) then
    170173          call readwind_ecmwf(indj,memind(1),uuh,vvh,wwh)
    171174        else
     
    173176        end if
    174177        call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
    175         call calcpar(memind(1),uuh,vvh,pvh,metdata_format)
    176         call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
    177         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     178        call calcpar(memind(1),uuh,vvh,pvh,id_centre)
     179        call calcpar_nests(memind(1),uuhn,vvhn,pvhn,id_centre)
     180        if (id_centre.eq.icg_id_ecmwf) then
    178181          call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh)
    179182        else
     
    183186        memtime(1)=wftime(indj)
    184187        memind(2)=2
    185         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     188        if (id_centre.eq.icg_id_ecmwf) then
    186189          call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
    187190        else
     
    189192        end if
    190193        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
    191         call calcpar(memind(2),uuh,vvh,pvh,metdata_format)
    192         call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
    193         if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     194        call calcpar(memind(2),uuh,vvh,pvh,id_centre)
     195        call calcpar_nests(memind(2),uuhn,vvhn,pvhn,id_centre)
     196        if (id_centre.eq.icg_id_ecmwf) then
    194197          call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
    195198        else
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG