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

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since db91eb7 was db91eb7, checked in by Sabine <sabine.eckhardt@…>, 5 years ago

reading clouds in NCEP

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