source: flexpart.git/src/readwind_gfs.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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