source: flexpart.git/src/readwind_gfs.f90 @ 5d74ed9

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

BUG: var numpclwch in readwind_gfs was not defined

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