Changeset 6ecb30a in flexpart.git for src/getfields.f90


Ignore:
Timestamp:
Aug 17, 2017, 4:39:17 PM (7 years ago)
Author:
Espen Sollum ATMOS <eso@…>
Branches:
master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug, univie
Children:
5b34509
Parents:
61e07ba
Message:

Merged changes from CTBTO project

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/getfields.f90

    rd1a8707 r6ecb30a  
    2020!**********************************************************************
    2121
    22 subroutine getfields(itime,nstop)
     22subroutine getfields(itime,nstop,metdata_format)
    2323!                       i     o
    2424!*****************************************************************************
     
    3838!                                                                            *
    3939!*****************************************************************************
    40 !  Changes, Bernd C. Krueger, Feb. 2001:
    41 !   Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block.
    42 !   Function of nstop extended.
     40!  Changes, Bernd C. Krueger, Feb. 2001:                                     *
     41!   Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block.        *
     42!   Function of nstop extended.                                              *
     43!                                                                            *
     44!   Unified ECMWF and GFS builds                                             *
     45!   Marian Harustak, 12.5.2017                                               *
     46!     - Added passing of metdata_format as it was needed by called routines  *
    4347!*****************************************************************************
    4448!                                                                            *
     
    5862! tt(0:nxmax,0:nymax,nuvzmax,2)   temperature [K]                            *
    5963! ps(0:nxmax,0:nymax,2)           surface pressure [Pa]                      *
     64! metdata_format     format of metdata (ecmwf/gfs)                           *
    6065!                                                                            *
    6166! Constants:                                                                 *
     
    6772  use par_mod
    6873  use com_mod
     74  use class_gribfile
    6975
    7076  implicit none
    7177
    7278  integer :: indj,itime,nstop,memaux
     79  integer :: metdata_format
    7380
    7481  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
     
    125132    do indj=indmin,numbwf-1
    126133      if (ldirect*wftime(indj+1).gt.ldirect*itime) then
    127         call readwind(indj+1,memind(2),uuh,vvh,wwh)
     134        if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     135          call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
     136        else
     137          call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh)
     138        end if
    128139        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
    129         call calcpar(memind(2),uuh,vvh,pvh)
    130         call calcpar_nests(memind(2),uuhn,vvhn,pvhn)
    131         call verttransform(memind(2),uuh,vvh,wwh,pvh)
     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 verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
     144        else
     145          call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh)
     146        end if
    132147        call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
    133148        memtime(2)=wftime(indj+1)
     
    152167           (ldirect*wftime(indj+1).gt.ldirect*itime)) then
    153168        memind(1)=1
    154         call readwind(indj,memind(1),uuh,vvh,wwh)
     169        if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     170          call readwind_ecmwf(indj,memind(1),uuh,vvh,wwh)
     171        else
     172          call readwind_gfs(indj,memind(1),uuh,vvh,wwh)
     173        end if
    155174        call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
    156         call calcpar(memind(1),uuh,vvh,pvh)
    157         call calcpar_nests(memind(1),uuhn,vvhn,pvhn)
    158         call verttransform(memind(1),uuh,vvh,wwh,pvh)
     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 verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh)
     179        else
     180          call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh)
     181        end if
    159182        call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
    160183        memtime(1)=wftime(indj)
    161184        memind(2)=2
    162         call readwind(indj+1,memind(2),uuh,vvh,wwh)
     185        if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
     186          call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
     187        else
     188          call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh)
     189        end if
    163190        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
    164         call calcpar(memind(2),uuh,vvh,pvh)
    165         call calcpar_nests(memind(2),uuhn,vvhn,pvhn)
    166         call verttransform(memind(2),uuh,vvh,wwh,pvh)
     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 verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
     195        else
     196          call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh)
     197        end if
    167198        call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
    168199        memtime(2)=wftime(indj+1)
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG