source: flexpart.git/src/readwind_gfs.f90 @ e2132cf

dev
Last change on this file since e2132cf was e2132cf, checked in by Ignacio Pisso <ip@…>, 3 years ago

fix permissions

  • Property mode set to 100644
File size: 25.1 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine readwind_gfs(indj,n,uuh,vvh,wwh)
5
6!***********************************************************************
7!*                                                                     *
8!*             TRAJECTORY MODEL SUBROUTINE READWIND                    *
9!*                                                                     *
10!***********************************************************************
11!*                                                                     *
12!*             AUTHOR:      G. WOTAWA                                  *
13!*             DATE:        1997-08-05                                 *
14!*             LAST UPDATE: 2000-10-17, Andreas Stohl                  *
15!*             CHANGE: 01/02/2001, Bernd C. Krueger, Variables tth and *
16!*                     qvh (on eta coordinates) in common block        *
17!*             CHANGE: 16/11/2005, Caroline Forster, GFS data          *
18!*             CHANGE: 11/01/2008, Harald Sodemann, Input of GRIB1/2   *
19!*                     data with the ECMWF grib_api library            *
20!*             CHANGE: 03/12/2008, Harald Sodemann, update to f90 with *
21!*                                 ECMWF grib_api                      *
22!                                                                      *
23!   Unified ECMWF and GFS builds                                       *
24!   Marian Harustak, 12.5.2017                                         *
25!     - Renamed routine from readwind to readwind_gfs                  *
26!*                                                                     *
27!***********************************************************************
28!*                                                                     *
29!* DESCRIPTION:                                                        *
30!*                                                                     *
31!* READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE   *
32!* INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE          *
33!*                                                                     *
34!* INPUT:                                                              *
35!* indj               indicates number of the wind field to be read in *
36!* n                  temporal index for meteorological fields (1 to 3)*
37!*                                                                     *
38!* IMPORTANT VARIABLES FROM COMMON BLOCK:                              *
39!*                                                                     *
40!* wfname             File name of data to be read in                  *
41!* nx,ny,nuvz,nwz     expected field dimensions                        *
42!* nlev_ec            number of vertical levels ecmwf model            *
43!* uu,vv,ww           wind fields                                      *
44!* tt,qv              temperature and specific humidity                *
45!* ps                 surface pressure                                 *
46!*                                                                     *
47!***********************************************************************
48
49  use eccodes
50  use par_mod
51  use com_mod
52
53  implicit none
54
55!HSO  new parameters for grib_api
56  integer :: ifile
57  integer :: iret
58  integer :: igrib
59  integer :: gribVer,parCat,parNum,typSurf,discipl,valSurf
60!HSO end edits
61  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
62  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
63  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
64  integer :: ii,indj,i,j,k,n,levdiff2,ifield,iumax,iwmax
65
66! NCEP
67  integer :: numpt,numpu,numpv,numpw,numprh,numpclwch
68  real :: help, temp, ew
69  real :: elev
70  real :: ulev1(0:nxmax-1,0:nymax-1),vlev1(0:nxmax-1,0:nymax-1)
71  real :: tlev1(0:nxmax-1,0:nymax-1)
72  real :: qvh2(0:nxmax-1,0:nymax-1)
73
74  integer :: i180
75
76! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
77!HSO kept isec1, isec2 and zsec4 for consistency with gribex GRIB input
78
79  integer :: isec1(8),isec2(3)
80  real           :: xsec18
81  real(kind=4) :: zsec4(jpunp)
82  real(kind=4) :: xaux,yaux,xaux0,yaux0
83  real,parameter :: eps=spacing(2.0_4*360.0_4)
84  real(kind=8) :: xauxin,yauxin
85  real(kind=4) :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1)
86  real :: plev1,hlev1,ff10m,fflev1
87
88  logical :: hflswitch,strswitch
89
90!HSO  for grib api error messages
91  character(len=24) :: gribErrorMsg = 'Error reading grib file'
92  character(len=20) :: gribFunction = 'readwind_gfs'
93  character(len=20) :: shortname
94
95
96  hflswitch=.false.
97  strswitch=.false.
98  levdiff2=nlev_ec-nwz+1
99  iumax=0
100  iwmax=0
101
102
103! OPENING OF DATA FILE (GRIB CODE)
104
105!HSO
106  call grib_open_file(ifile,path(3)(1:length(3)) &
107       //trim(wfname(indj)),'r',iret)
108  if (iret.ne.GRIB_SUCCESS) then
109    goto 888   ! ERROR DETECTED
110  endif
111!turn on support for multi fields messages
112  call grib_multi_support_on
113
114  numpt=0
115  numpu=0
116  numpv=0
117  numpw=0
118  numprh=0
119  numpclwch=0
120  ifield=0
12110 ifield=ifield+1
122!
123! GET NEXT FIELDS
124!
125  call grib_new_from_file(ifile,igrib,iret)
126  if (iret.eq.GRIB_END_OF_FILE)  then
127    goto 50    ! EOF DETECTED
128  elseif (iret.ne.GRIB_SUCCESS) then
129    goto 888   ! ERROR DETECTED
130  endif
131
132!first see if we read GRIB1 or GRIB2
133  call grib_get_int(igrib,'editionNumber',gribVer,iret)
134!  call grib_check(iret,gribFunction,gribErrorMsg)
135
136  if (gribVer.eq.1) then ! GRIB Edition 1
137
138!read the grib1 identifiers
139
140    call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
141    call grib_check(iret,gribFunction,gribErrorMsg)
142    call grib_get_int(igrib,'indicatorOfTypeOfLevel',isec1(7),iret)
143    call grib_check(iret,gribFunction,gribErrorMsg)
144    call grib_get_int(igrib,'level',isec1(8),iret)
145    call grib_check(iret,gribFunction,gribErrorMsg)
146!JMA / SH: isec1(8) not evaluated any more below
147!b/c with GRIB 2 this may be a real variable
148    xsec18 = real(isec1(8))
149
150  else ! GRIB Edition 2
151
152!read the grib2 identifiers
153   
154    call grib_get_string(igrib,'shortName',shortname,iret)
155
156    call grib_get_int(igrib,'discipline',discipl,iret)
157    call grib_check(iret,gribFunction,gribErrorMsg)
158    call grib_get_int(igrib,'parameterCategory',parCat,iret)
159    call grib_check(iret,gribFunction,gribErrorMsg)
160    call grib_get_int(igrib,'parameterNumber',parNum,iret)
161    call grib_check(iret,gribFunction,gribErrorMsg)
162    call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
163    call grib_check(iret,gribFunction,gribErrorMsg)
164!
165    call grib_get_int(igrib,'scaledValueOfFirstFixedSurface', &
166         valSurf,iret)
167    call grib_check(iret,gribFunction,gribErrorMsg)
168
169!    write(*,*) 'Field: ',ifield,parCat,parNum,typSurf,shortname
170!convert to grib1 identifiers
171    isec1(6:8)=-1
172    xsec18  =-1.0
173!JMA / SH: isec1(8) not evaluated any more below
174!b/c with GRIB 2 this may be a real variable
175    if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.100)) then ! T
176      isec1(6)=11          ! indicatorOfParameter
177      isec1(7)=100         ! indicatorOfTypeOfLevel
178      xsec18=valSurf/100.0 ! level, convert to hPa
179    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.100)) then ! U
180      isec1(6)=33          ! indicatorOfParameter
181      isec1(7)=100         ! indicatorOfTypeOfLevel
182      xsec18=valSurf/100.0 ! level, convert to hPa
183    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.100)) then ! V
184      isec1(6)=34          ! indicatorOfParameter
185      isec1(7)=100         ! indicatorOfTypeOfLevel
186      xsec18=valSurf/100.0 ! level, convert to hPa
187    elseif ((parCat.eq.2).and.(parNum.eq.8).and.(typSurf.eq.100)) then ! W
188      isec1(6)=39          ! indicatorOfParameter
189      isec1(7)=100         ! indicatorOfTypeOfLevel
190      xsec18=valSurf/100.0 ! level, convert to hPa
191    elseif ((parCat.eq.1).and.(parNum.eq.1).and.(typSurf.eq.100)) then ! RH
192      isec1(6)=52          ! indicatorOfParameter
193      isec1(7)=100         ! indicatorOfTypeOfLevel
194      xsec18=valSurf/100.0 ! level, convert to hPa
195    elseif ((parCat.eq.1).and.(parNum.eq.1).and.(typSurf.eq.103)) then ! RH2
196      isec1(6)=52          ! indicatorOfParameter
197      isec1(7)=105         ! indicatorOfTypeOfLevel
198      xsec18=real(2)
199    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! T2
200      isec1(6)=11          ! indicatorOfParameter
201      isec1(7)=105         ! indicatorOfTypeOfLevel
202      xsec18=real(2)
203    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! U10
204      isec1(6)=33          ! indicatorOfParameter
205      isec1(7)=105         ! indicatorOfTypeOfLevel
206      xsec18=real(10)
207    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! V10
208      isec1(6)=34          ! indicatorOfParameter
209      isec1(7)=105         ! indicatorOfTypeOfLevel
210      xsec18=real(10)
211    elseif ((parCat.eq.1).and.(parNum.eq.22).and.(typSurf.eq.100)) then ! CLWMR Cloud Mixing Ratio [kg/kg]:
212      isec1(6)=153         ! indicatorOfParameter
213      isec1(7)=100         ! indicatorOfTypeOfLevel
214      xsec18=valSurf/100.0 ! level, convert to hPa
215    elseif ((parCat.eq.3).and.(parNum.eq.1).and.(typSurf.eq.101)) then ! SLP
216      isec1(6)=2           ! indicatorOfParameter
217      isec1(7)=102         ! indicatorOfTypeOfLevel
218      xsec18=real(0)
219   elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1).and. &
220        (discipl.eq.0)) then ! SP
221      isec1(6)=1           ! indicatorOfParameter
222      isec1(7)=1           ! indicatorOfTypeOfLevel
223      xsec18=real(0)
224    elseif ((parCat.eq.1).and.(parNum.eq.13).and.(typSurf.eq.1)) then ! SNOW
225      isec1(6)=66          ! indicatorOfParameter
226      isec1(7)=1           ! indicatorOfTypeOfLevel
227      xsec18=real(0)
228    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.104)) then ! T sigma 0
229      isec1(6)=11          ! indicatorOfParameter
230      isec1(7)=107         ! indicatorOfTypeOfLevel
231      xsec18=0.995         ! lowest sigma level
232    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.104)) then ! U sigma 0
233      isec1(6)=33          ! indicatorOfParameter
234      isec1(7)=107         ! indicatorOfTypeOfLevel
235      xsec18=0.995         ! lowest sigma level
236    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.104)) then ! V sigma 0
237      isec1(6)=34          ! indicatorOfParameter
238      isec1(7)=107         ! indicatorOfTypeOfLevel
239      xsec18=0.995         ! lowest sigma level
240    elseif ((parCat.eq.3).and.(parNum.eq.5).and.(typSurf.eq.1)) then ! TOPO
241      isec1(6)=7           ! indicatorOfParameter
242      isec1(7)=1           ! indicatorOfTypeOfLevel
243      xsec18=real(0)
244    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.1) &
245         .and.(discipl.eq.2)) then ! LSM
246      isec1(6)=81          ! indicatorOfParameter
247      isec1(7)=1           ! indicatorOfTypeOfLevel
248      xsec18=real(0)
249    elseif ((parCat.eq.3).and.(parNum.eq.196).and.(typSurf.eq.1)) then ! BLH
250      isec1(6)=221         ! indicatorOfParameter
251      isec1(7)=1           ! indicatorOfTypeOfLevel
252      xsec18=real(0)
253    elseif ((parCat.eq.1).and.(parNum.eq.7).and.(typSurf.eq.1)) then ! LSP/TP
254      isec1(6)=62          ! indicatorOfParameter
255      isec1(7)=1           ! indicatorOfTypeOfLevel
256      xsec18=real(0)
257    elseif ((parCat.eq.1).and.(parNum.eq.196).and.(typSurf.eq.1)) then ! CP
258      isec1(6)=63          ! indicatorOfParameter
259      isec1(7)=1           ! indicatorOfTypeOfLevel
260      xsec18=real(0)
261    endif
262
263  endif ! gribVer
264
265  if (isec1(6).ne.-1) then
266!  get the size and data of the values array
267    call grib_get_real4_array(igrib,'values',zsec4,iret)
268    call grib_check(iret,gribFunction,gribErrorMsg)
269  endif
270
271  if(ifield.eq.1) then
272
273!get the required fields from section 2
274!store compatible to gribex input
275    call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
276         isec2(2),iret)
277    call grib_check(iret,gribFunction,gribErrorMsg)
278    call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
279         isec2(3),iret)
280    call grib_check(iret,gribFunction,gribErrorMsg)
281    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
282         xauxin,iret)
283    call grib_check(iret,gribFunction,gribErrorMsg)
284    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
285         yauxin,iret)
286    call grib_check(iret,gribFunction,gribErrorMsg)
287    xaux=xauxin+real(nxshift)*dx
288    yaux=yauxin
289
290! CHECK GRID SPECIFICATIONS
291
292    if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
293    if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
294    if(xaux.eq.0.) xaux=-180.0     ! NCEP DATA
295    xaux0=xlon0
296    yaux0=ylat0
297    if(xaux.lt.0.) xaux=xaux+360.
298    if(yaux.lt.0.) yaux=yaux+360.
299    if(xaux0.lt.0.) xaux0=xaux0+360.
300    if(yaux0.lt.0.) yaux0=yaux0+360.
301    if(abs(xaux-xaux0).gt.eps) then
302      write (*, *) xaux, xaux0
303      stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
304    endif
305    if(abs(yaux-yaux0).gt.eps) then
306      write (*, *) yaux, yaux0
307      stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
308    end if
309  endif
310!HSO end of edits
311
312  i180=nint(180./dx)
313
314  if (isec1(6).ne.-1) then
315
316! write (*, *) 'nxfield: ', nxfield, i180
317
318    do j=0,nymin1
319      do i=0,nxfield-1
320        if((isec1(6).eq.011).and.(isec1(7).eq.100)) then
321! TEMPERATURE
322          if((i.eq.0).and.(j.eq.0)) then
323            numpt=minloc(abs(xsec18*100.0-akz),dim=1)
324          endif
325          help=zsec4(nxfield*(ny-j-1)+i+1)
326          if (help.eq.0) then
327            write (*, *) 'i, j: ', i, j
328            stop 'help == 0.0'
329          endif
330          if(i.lt.i180) then
331            tth(i180+i,j,numpt,n)=help
332          else
333            tth(i-i180,j,numpt,n)=help
334          endif
335        endif
336        if((isec1(6).eq.033).and.(isec1(7).eq.100)) then
337! U VELOCITY
338          if((i.eq.0).and.(j.eq.0)) then
339            numpu=minloc(abs(xsec18*100.0-akz),dim=1)
340          endif
341          help=zsec4(nxfield*(ny-j-1)+i+1)
342          if(i.lt.i180) then
343            uuh(i180+i,j,numpu)=help
344          else
345            uuh(i-i180,j,numpu)=help
346          endif
347        endif
348        if((isec1(6).eq.034).and.(isec1(7).eq.100)) then
349! V VELOCITY
350          if((i.eq.0).and.(j.eq.0)) then
351            numpv=minloc(abs(xsec18*100.0-akz),dim=1)
352          endif
353          help=zsec4(nxfield*(ny-j-1)+i+1)
354          if(i.lt.i180) then
355            vvh(i180+i,j,numpv)=help
356          else
357            vvh(i-i180,j,numpv)=help
358          endif
359        endif
360        if((isec1(6).eq.052).and.(isec1(7).eq.100)) then
361! RELATIVE HUMIDITY -> CONVERT TO SPECIFIC HUMIDITY LATER
362          if((i.eq.0).and.(j.eq.0)) then
363            numprh=minloc(abs(xsec18*100.0-akz),dim=1)
364          endif
365          help=zsec4(nxfield*(ny-j-1)+i+1)
366          if(i.lt.i180) then
367            qvh(i180+i,j,numprh,n)=help
368          else
369            qvh(i-i180,j,numprh,n)=help
370          endif
371        endif
372        if((isec1(6).eq.001).and.(isec1(7).eq.001)) then
373! SURFACE PRESSURE
374          help=zsec4(nxfield*(ny-j-1)+i+1)
375          if(i.lt.i180) then
376            ps(i180+i,j,1,n)=help
377          else
378            ps(i-i180,j,1,n)=help
379          endif
380        endif
381        if((isec1(6).eq.039).and.(isec1(7).eq.100)) then
382! W VELOCITY
383          if((i.eq.0).and.(j.eq.0)) numpw=minloc(abs(xsec18*100.0-akz),dim=1)
384          help=zsec4(nxfield*(ny-j-1)+i+1)
385          if(i.lt.i180) then
386            wwh(i180+i,j,numpw)=help
387          else
388            wwh(i-i180,j,numpw)=help
389          endif
390        endif
391        if((isec1(6).eq.066).and.(isec1(7).eq.001)) then
392! SNOW DEPTH
393          help=zsec4(nxfield*(ny-j-1)+i+1)
394          if(i.lt.i180) then
395            sd(i180+i,j,1,n)=help
396          else
397            sd(i-i180,j,1,n)=help
398          endif
399        endif
400        if((isec1(6).eq.002).and.(isec1(7).eq.102)) then
401! MEAN SEA LEVEL PRESSURE
402          help=zsec4(nxfield*(ny-j-1)+i+1)
403          if(i.lt.i180) then
404            msl(i180+i,j,1,n)=help
405          else
406            msl(i-i180,j,1,n)=help
407          endif
408        endif
409        if((isec1(6).eq.071).and.(isec1(7).eq.244)) then
410! TOTAL CLOUD COVER
411          help=zsec4(nxfield*(ny-j-1)+i+1)
412          if(i.lt.i180) then
413            tcc(i180+i,j,1,n)=help
414          else
415            tcc(i-i180,j,1,n)=help
416          endif
417        endif
418        if((isec1(6).eq.033).and.(isec1(7).eq.105).and. &
419             (nint(xsec18).eq.10)) then
420! 10 M U VELOCITY
421          help=zsec4(nxfield*(ny-j-1)+i+1)
422          if(i.lt.i180) then
423            u10(i180+i,j,1,n)=help
424          else
425            u10(i-i180,j,1,n)=help
426          endif
427        endif
428        if((isec1(6).eq.034).and.(isec1(7).eq.105).and. &
429             (nint(xsec18).eq.10)) then
430! 10 M V VELOCITY
431          help=zsec4(nxfield*(ny-j-1)+i+1)
432          if(i.lt.i180) then
433            v10(i180+i,j,1,n)=help
434          else
435            v10(i-i180,j,1,n)=help
436          endif
437        endif
438        if((isec1(6).eq.011).and.(isec1(7).eq.105).and. &
439             (nint(xsec18).eq.2)) then
440! 2 M TEMPERATURE
441          help=zsec4(nxfield*(ny-j-1)+i+1)
442          if(i.lt.i180) then
443            tt2(i180+i,j,1,n)=help
444          else
445            tt2(i-i180,j,1,n)=help
446          endif
447        endif
448        if((isec1(6).eq.017).and.(isec1(7).eq.105).and. &
449             (nint(xsec18).eq.2)) then
450! 2 M DEW POINT TEMPERATURE
451          help=zsec4(nxfield*(ny-j-1)+i+1)
452          if(i.lt.i180) then
453            td2(i180+i,j,1,n)=help
454          else
455            td2(i-i180,j,1,n)=help
456          endif
457        endif
458        if((isec1(6).eq.062).and.(isec1(7).eq.001)) then
459! LARGE SCALE PREC.
460          help=zsec4(nxfield*(ny-j-1)+i+1)
461          if(i.lt.i180) then
462            lsprec(i180+i,j,1,n)=help
463          else
464            lsprec(i-i180,j,1,n)=help
465          endif
466        endif
467        if((isec1(6).eq.063).and.(isec1(7).eq.001)) then
468! CONVECTIVE PREC.
469          help=zsec4(nxfield*(ny-j-1)+i+1)
470          if(i.lt.i180) then
471            convprec(i180+i,j,1,n)=help
472          else
473            convprec(i-i180,j,1,n)=help
474          endif
475        endif
476        if((isec1(6).eq.007).and.(isec1(7).eq.001)) then
477! TOPOGRAPHY
478          help=zsec4(nxfield*(ny-j-1)+i+1)
479          if(i.lt.i180) then
480            oro(i180+i,j)=help
481            excessoro(i180+i,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
482          else
483            oro(i-i180,j)=help
484            excessoro(i-i180,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
485          endif
486        endif
487        if((isec1(6).eq.081).and.(isec1(7).eq.001)) then
488! LAND SEA MASK
489          help=zsec4(nxfield*(ny-j-1)+i+1)
490          if(i.lt.i180) then
491            lsm(i180+i,j)=help
492          else
493            lsm(i-i180,j)=help
494          endif
495        endif
496        if((isec1(6).eq.221).and.(isec1(7).eq.001)) then
497! MIXING HEIGHT
498          help=zsec4(nxfield*(ny-j-1)+i+1)
499          if(i.lt.i180) then
500            hmix(i180+i,j,1,n)=help
501          else
502            hmix(i-i180,j,1,n)=help
503          endif
504        endif
505        if((isec1(6).eq.052).and.(isec1(7).eq.105).and. &
506             (nint(xsec18).eq.02)) then
507! 2 M RELATIVE HUMIDITY
508          help=zsec4(nxfield*(ny-j-1)+i+1)
509          if(i.lt.i180) then
510            qvh2(i180+i,j)=help
511          else
512            qvh2(i-i180,j)=help
513          endif
514        endif
515        if((isec1(6).eq.011).and.(isec1(7).eq.107)) then
516! TEMPERATURE LOWEST SIGMA LEVEL
517          help=zsec4(nxfield*(ny-j-1)+i+1)
518          if(i.lt.i180) then
519            tlev1(i180+i,j)=help
520          else
521            tlev1(i-i180,j)=help
522          endif
523        endif
524        if((isec1(6).eq.033).and.(isec1(7).eq.107)) then
525! U VELOCITY LOWEST SIGMA LEVEL
526          help=zsec4(nxfield*(ny-j-1)+i+1)
527          if(i.lt.i180) then
528            ulev1(i180+i,j)=help
529          else
530            ulev1(i-i180,j)=help
531          endif
532        endif
533        if((isec1(6).eq.034).and.(isec1(7).eq.107)) then
534! V VELOCITY LOWEST SIGMA LEVEL
535          help=zsec4(nxfield*(ny-j-1)+i+1)
536          if(i.lt.i180) then
537            vlev1(i180+i,j)=help
538          else
539            vlev1(i-i180,j)=help
540          endif
541        endif
542! SEC & IP 12/2018 read GFS clouds
543        if((isec1(6).eq.153).and.(isec1(7).eq.100)) then  !! CLWCR  Cloud liquid water content [kg/kg]
544          if((i.eq.0).and.(j.eq.0)) then
545            numpclwch=minloc(abs(xsec18*100.0-akz),dim=1)
546          endif
547          help=zsec4(nxfield*(ny-j-1)+i+1)
548          if(i.lt.i180) then
549            clwch(i180+i,j,numpclwch,n)=help
550          else
551            clwch(i-i180,j,numpclwch,n)=help
552          endif
553          readclouds=.true.
554          sumclouds=.true.
555        endif
556
557
558      end do
559    end do
560
561  endif
562
563  if((isec1(6).eq.33).and.(isec1(7).eq.100)) then
564! NCEP ISOBARIC LEVELS
565    iumax=iumax+1
566  endif
567
568  call grib_release(igrib)
569  goto 10                      !! READ NEXT LEVEL OR PARAMETER
570!
571! CLOSING OF INPUT DATA FILE
572!
573
574!HSO close grib file
57550 continue
576  call grib_close_file(ifile)
577
578! SENS. HEAT FLUX
579  sshf(:,:,1,n)=0.0     ! not available from gfs.tccz.pgrbfxx files
580  hflswitch=.false.    ! Heat flux not available
581! SOLAR RADIATIVE FLUXES
582  ssr(:,:,1,n)=0.0      ! not available from gfs.tccz.pgrbfxx files
583! EW SURFACE STRESS
584  ewss=0.0         ! not available from gfs.tccz.pgrbfxx files
585! NS SURFACE STRESS
586  nsss=0.0         ! not available from gfs.tccz.pgrbfxx files
587  strswitch=.false.    ! stress not available
588
589! CONVERT TP TO LSP (GRIB2 only)
590  if (gribVer.eq.2) then
591    lsprec(0:nxfield-1, 0:nymin1, 1, n) = max( 0.0 , lsprec(0:nxfield-1,0:nymin1,1,n)-convprec(0:nxfield-1,0:nymin1,1,n) )
592  endif
593!HSO end edits
594
595
596! TRANSFORM RH TO SPECIFIC HUMIDITY
597
598  do j=0,ny-1
599    do i=0,nxfield-1
600      do k=1,nuvz
601        help=qvh(i,j,k,n)
602        temp=tth(i,j,k,n)
603        if (temp .eq. 0.0) then
604          write (*, *) i, j, k, n
605!          temp = 273.0
606          stop
607        endif
608        plev1=akm(k)+bkm(k)*ps(i,j,1,n)
609        elev=ew(temp)*help/100.0
610        qvh(i,j,k,n)=xmwml*(elev/(plev1-((1.0-xmwml)*elev)))
611      end do
612    end do
613  end do
614
615! CALCULATE 2 M DEW POINT FROM 2 M RELATIVE HUMIDITY
616! USING BOLTON'S (1980) FORMULA
617! BECAUSE td2 IS NOT AVAILABLE FROM NCEP GFS DATA
618
619  do j=0,ny-1
620    do i=0,nxfield-1
621      help=qvh2(i,j)
622      temp=tt2(i,j,1,n)
623      if (temp .eq. 0.0) then
624        write (*, *) i, j, n
625        temp = 273.0
626      endif
627      elev=ew(temp)/100.*help/100.   !vapour pressure in hPa
628      td2(i,j,1,n)=243.5/(17.67/log(elev/6.112)-1)+273.
629      if (help.le.0.) td2(i,j,1,n)=tt2(i,j,1,n)
630    end do
631  end do
632
633  if(levdiff2.eq.0) then
634    iwmax=nlev_ec+1
635    do i=0,nxmin1
636      do j=0,nymin1
637        wwh(i,j,nlev_ec+1)=0.
638      end do
639    end do
640  endif
641
642
643! For global fields, assign the leftmost data column also to the rightmost
644! data column; if required, shift whole grid by nxshift grid points
645!*************************************************************************
646
647  if (xglobal) then
648    call shift_field_0(ewss,nxfield,ny)
649    call shift_field_0(nsss,nxfield,ny)
650    call shift_field_0(oro,nxfield,ny)
651    call shift_field_0(excessoro,nxfield,ny)
652    call shift_field_0(lsm,nxfield,ny)
653    call shift_field_0(ulev1,nxfield,ny)
654    call shift_field_0(vlev1,nxfield,ny)
655    call shift_field_0(tlev1,nxfield,ny)
656    call shift_field_0(qvh2,nxfield,ny)
657    call shift_field(ps,nxfield,ny,1,1,2,n)
658    call shift_field(sd,nxfield,ny,1,1,2,n)
659    call shift_field(msl,nxfield,ny,1,1,2,n)
660    call shift_field(tcc,nxfield,ny,1,1,2,n)
661    call shift_field(u10,nxfield,ny,1,1,2,n)
662    call shift_field(v10,nxfield,ny,1,1,2,n)
663    call shift_field(tt2,nxfield,ny,1,1,2,n)
664    call shift_field(td2,nxfield,ny,1,1,2,n)
665    call shift_field(lsprec,nxfield,ny,1,1,2,n)
666    call shift_field(convprec,nxfield,ny,1,1,2,n)
667    call shift_field(sshf,nxfield,ny,1,1,2,n)
668    call shift_field(ssr,nxfield,ny,1,1,2,n)
669    call shift_field(hmix,nxfield,ny,1,1,2,n)
670    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
671    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
672    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
673    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
674    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
675! IP & SEC adding GFS Clouds 20181205
676    call shift_field(clwch,nxfield,ny,nuvzmax,nuvz,2,n)
677  endif
678
679  do i=0,nxmin1
680    do j=0,nymin1
681! Convert precip. from mm/s -> mm/hour
682      convprec(i,j,1,n)=convprec(i,j,1,n)*3600.
683      lsprec(i,j,1,n)=lsprec(i,j,1,n)*3600.
684      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
685    end do
686  end do
687
688  if ((.not.hflswitch).or.(.not.strswitch)) then
689!  write(*,*) 'WARNING: No flux data contained in GRIB file ',
690!    +  wfname(indj)
691
692! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
693!***************************************************************************
694
695    do i=0,nxmin1
696      do j=0,nymin1
697        hlev1=30.0                     ! HEIGHT OF FIRST MODEL SIGMA LAYER
698        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
699        fflev1=sqrt(ulev1(i,j)**2+vlev1(i,j)**2)
700        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
701             tt2(i,j,1,n),tlev1(i,j),ff10m,fflev1, &
702             surfstr(i,j,1,n),sshf(i,j,1,n))
703        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
704        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
705      end do
706    end do
707  endif
708
709  if(iumax.ne.nuvz) stop 'READWIND: NUVZ NOT CONSISTENT'
710  if(iumax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
711
712  return
713888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
714  write(*,*) ' #### ',wfname(indj),'                    #### '
715  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
716  stop 'Execution terminated'
717999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
718  write(*,*) ' #### ',wfname(indj),'                    #### '
719  write(*,*) ' #### CANNOT BE OPENED !!!                    #### '
720  stop 'Execution terminated'
721
722end subroutine readwind_gfs
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG