source: flexpart.git/src/readwind_gfs.f90 @ 6ecb30a

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 6ecb30a was 6ecb30a, checked in by Espen Sollum ATMOS <eso@…>, 7 years ago

Merged changes from CTBTO project

  • Property mode set to 100644
File size: 25.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
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
111
112  hflswitch=.false.
113  strswitch=.false.
114  levdiff2=nlev_ec-nwz+1
115  iumax=0
116  iwmax=0
117
118
119  ! OPENING OF DATA FILE (GRIB CODE)
120
121  !HSO
1225   call grib_open_file(ifile,path(3)(1:length(3)) &
123         //trim(wfname(indj)),'r',iret)
124  if (iret.ne.GRIB_SUCCESS) then
125    goto 888   ! ERROR DETECTED
126  endif
127  !turn on support for multi fields messages
128  call grib_multi_support_on
129
130  numpt=0
131  numpu=0
132  numpv=0
133  numpw=0
134  numprh=0
135  ifield=0
13610   ifield=ifield+1
137  !
138  ! GET NEXT FIELDS
139  !
140  call grib_new_from_file(ifile,igrib,iret)
141  if (iret.eq.GRIB_END_OF_FILE)  then
142    goto 50    ! EOF DETECTED
143  elseif (iret.ne.GRIB_SUCCESS) then
144    goto 888   ! ERROR DETECTED
145  endif
146
147  !first see if we read GRIB1 or GRIB2
148  call grib_get_int(igrib,'editionNumber',gribVer,iret)
149!  call grib_check(iret,gribFunction,gribErrorMsg)
150
151  if (gribVer.eq.1) then ! GRIB Edition 1
152
153  !read the grib1 identifiers
154  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
155!  call grib_check(iret,gribFunction,gribErrorMsg)
156  call grib_get_int(igrib,'indicatorOfTypeOfLevel',isec1(7),iret)
157!  call grib_check(iret,gribFunction,gribErrorMsg)
158  call grib_get_int(igrib,'level',isec1(8),iret)
159!  call grib_check(iret,gribFunction,gribErrorMsg)
160
161  else ! GRIB Edition 2
162
163  !read the grib2 identifiers
164  call grib_get_int(igrib,'discipline',discipl,iret)
165!  call grib_check(iret,gribFunction,gribErrorMsg)
166  call grib_get_int(igrib,'parameterCategory',parCat,iret)
167!  call grib_check(iret,gribFunction,gribErrorMsg)
168  call grib_get_int(igrib,'parameterNumber',parNum,iret)
169!  call grib_check(iret,gribFunction,gribErrorMsg)
170  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
171!  call grib_check(iret,gribFunction,gribErrorMsg)
172  call grib_get_int(igrib,'scaledValueOfFirstFixedSurface', &
173       valSurf,iret)
174!  call grib_check(iret,gribFunction,gribErrorMsg)
175
176  !convert to grib1 identifiers
177  isec1(6)=-1
178  isec1(7)=-1
179  isec1(8)=-1
180  if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.100)) then ! T
181    isec1(6)=11          ! indicatorOfParameter
182    isec1(7)=100         ! indicatorOfTypeOfLevel
183    isec1(8)=valSurf/100 ! level, convert to hPa
184  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.100)) then ! U
185    isec1(6)=33          ! indicatorOfParameter
186    isec1(7)=100         ! indicatorOfTypeOfLevel
187    isec1(8)=valSurf/100 ! level, convert to hPa
188  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.100)) then ! V
189    isec1(6)=34          ! indicatorOfParameter
190    isec1(7)=100         ! indicatorOfTypeOfLevel
191    isec1(8)=valSurf/100 ! level, convert to hPa
192  elseif ((parCat.eq.2).and.(parNum.eq.8).and.(typSurf.eq.100)) then ! W
193    isec1(6)=39          ! indicatorOfParameter
194    isec1(7)=100         ! indicatorOfTypeOfLevel
195    isec1(8)=valSurf/100 ! level, convert to hPa
196  elseif ((parCat.eq.1).and.(parNum.eq.1).and.(typSurf.eq.100)) then ! RH
197    isec1(6)=52          ! 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.103)) then ! RH2
201    isec1(6)=52          ! indicatorOfParameter
202    isec1(7)=105         ! indicatorOfTypeOfLevel
203    isec1(8)=2
204  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! T2
205    isec1(6)=11          ! indicatorOfParameter
206    isec1(7)=105         ! indicatorOfTypeOfLevel
207    isec1(8)=2
208  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! U10
209    isec1(6)=33          ! indicatorOfParameter
210    isec1(7)=105         ! indicatorOfTypeOfLevel
211    isec1(8)=10
212  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! V10
213    isec1(6)=34          ! indicatorOfParameter
214    isec1(7)=105         ! indicatorOfTypeOfLevel
215    isec1(8)=10
216  elseif ((parCat.eq.3).and.(parNum.eq.1).and.(typSurf.eq.101)) then ! SLP
217    isec1(6)=2           ! indicatorOfParameter
218    isec1(7)=102         ! indicatorOfTypeOfLevel
219    isec1(8)=0
220  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then ! SP
221    isec1(6)=1           ! indicatorOfParameter
222    isec1(7)=1           ! indicatorOfTypeOfLevel
223    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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    isec1(8)=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=-179.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) &
302         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
303    if(abs(yaux-yaux0).gt.eps) &
304         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
305  endif
306  !HSO end of edits
307
308  i179=nint(179./dx)
309  if (dx.lt.0.7) then
310    i180=nint(180./dx)+1    ! 0.5 deg data
311  else
312    i180=nint(179./dx)+1    ! 1 deg data
313  endif
314  i181=i180+1
315
316  if (isec1(6).ne.-1) then
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            do ii=1,nuvz
324              if ((isec1(8)*100.0).eq.akz(ii)) numpt=ii
325            end do
326        endif
327        help=zsec4(nxfield*(ny-j-1)+i+1)
328        if(i.le.i180) then
329          tth(i179+i,j,numpt,n)=help
330        else
331          tth(i-i181,j,numpt,n)=help
332        endif
333      endif
334      if((isec1(6).eq.033).and.(isec1(7).eq.100)) then
335  ! U VELOCITY
336         if((i.eq.0).and.(j.eq.0)) then
337            do ii=1,nuvz
338              if ((isec1(8)*100.0).eq.akz(ii)) numpu=ii
339            end do
340        endif
341        help=zsec4(nxfield*(ny-j-1)+i+1)
342        if(i.le.i180) then
343          uuh(i179+i,j,numpu)=help
344        else
345          uuh(i-i181,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            do ii=1,nuvz
352              if ((isec1(8)*100.0).eq.akz(ii)) numpv=ii
353            end do
354        endif
355        help=zsec4(nxfield*(ny-j-1)+i+1)
356        if(i.le.i180) then
357          vvh(i179+i,j,numpv)=help
358        else
359          vvh(i-i181,j,numpv)=help
360        endif
361      endif
362      if((isec1(6).eq.052).and.(isec1(7).eq.100)) then
363  ! RELATIVE HUMIDITY -> CONVERT TO SPECIFIC HUMIDITY LATER
364         if((i.eq.0).and.(j.eq.0)) then
365            do ii=1,nuvz
366              if ((isec1(8)*100.0).eq.akz(ii)) numprh=ii
367            end do
368        endif
369        help=zsec4(nxfield*(ny-j-1)+i+1)
370        if(i.le.i180) then
371          qvh(i179+i,j,numprh,n)=help
372        else
373          qvh(i-i181,j,numprh,n)=help
374        endif
375      endif
376      if((isec1(6).eq.001).and.(isec1(7).eq.001)) then
377  ! SURFACE PRESSURE
378        help=zsec4(nxfield*(ny-j-1)+i+1)
379        if(i.le.i180) then
380          ps(i179+i,j,1,n)=help
381        else
382          ps(i-i181,j,1,n)=help
383        endif
384      endif
385      if((isec1(6).eq.039).and.(isec1(7).eq.100)) then
386  ! W VELOCITY
387         if((i.eq.0).and.(j.eq.0)) then
388            do ii=1,nuvz
389              if ((isec1(8)*100.0).eq.akz(ii)) numpw=ii
390            end do
391        endif
392        help=zsec4(nxfield*(ny-j-1)+i+1)
393        if(i.le.i180) then
394          wwh(i179+i,j,numpw)=help
395        else
396          wwh(i-i181,j,numpw)=help
397        endif
398      endif
399      if((isec1(6).eq.066).and.(isec1(7).eq.001)) then
400  ! SNOW DEPTH
401        help=zsec4(nxfield*(ny-j-1)+i+1)
402        if(i.le.i180) then
403          sd(i179+i,j,1,n)=help
404        else
405          sd(i-i181,j,1,n)=help
406        endif
407      endif
408      if((isec1(6).eq.002).and.(isec1(7).eq.102)) then
409  ! MEAN SEA LEVEL PRESSURE
410        help=zsec4(nxfield*(ny-j-1)+i+1)
411        if(i.le.i180) then
412          msl(i179+i,j,1,n)=help
413        else
414          msl(i-i181,j,1,n)=help
415        endif
416      endif
417      if((isec1(6).eq.071).and.(isec1(7).eq.244)) then
418  ! TOTAL CLOUD COVER
419        help=zsec4(nxfield*(ny-j-1)+i+1)
420        if(i.le.i180) then
421          tcc(i179+i,j,1,n)=help
422        else
423          tcc(i-i181,j,1,n)=help
424        endif
425      endif
426      if((isec1(6).eq.033).and.(isec1(7).eq.105).and. &
427           (isec1(8).eq.10)) then
428  ! 10 M U VELOCITY
429        help=zsec4(nxfield*(ny-j-1)+i+1)
430        if(i.le.i180) then
431          u10(i179+i,j,1,n)=help
432        else
433          u10(i-i181,j,1,n)=help
434        endif
435      endif
436      if((isec1(6).eq.034).and.(isec1(7).eq.105).and. &
437           (isec1(8).eq.10)) then
438  ! 10 M V VELOCITY
439        help=zsec4(nxfield*(ny-j-1)+i+1)
440        if(i.le.i180) then
441          v10(i179+i,j,1,n)=help
442        else
443          v10(i-i181,j,1,n)=help
444        endif
445      endif
446      if((isec1(6).eq.011).and.(isec1(7).eq.105).and. &
447           (isec1(8).eq.02)) then
448  ! 2 M TEMPERATURE
449        help=zsec4(nxfield*(ny-j-1)+i+1)
450        if(i.le.i180) then
451          tt2(i179+i,j,1,n)=help
452        else
453          tt2(i-i181,j,1,n)=help
454        endif
455      endif
456      if((isec1(6).eq.017).and.(isec1(7).eq.105).and. &
457           (isec1(8).eq.02)) then
458  ! 2 M DEW POINT TEMPERATURE
459        help=zsec4(nxfield*(ny-j-1)+i+1)
460        if(i.le.i180) then
461          td2(i179+i,j,1,n)=help
462        else
463          td2(i-i181,j,1,n)=help
464        endif
465      endif
466      if((isec1(6).eq.062).and.(isec1(7).eq.001)) then
467  ! LARGE SCALE PREC.
468        help=zsec4(nxfield*(ny-j-1)+i+1)
469        if(i.le.i180) then
470          lsprec(i179+i,j,1,n)=help
471        else
472          lsprec(i-i181,j,1,n)=help
473        endif
474      endif
475      if((isec1(6).eq.063).and.(isec1(7).eq.001)) then
476  ! CONVECTIVE PREC.
477        help=zsec4(nxfield*(ny-j-1)+i+1)
478        if(i.le.i180) then
479          convprec(i179+i,j,1,n)=help
480        else
481          convprec(i-i181,j,1,n)=help
482        endif
483      endif
484      if((isec1(6).eq.007).and.(isec1(7).eq.001)) then
485  ! TOPOGRAPHY
486        help=zsec4(nxfield*(ny-j-1)+i+1)
487        if(i.le.i180) then
488          oro(i179+i,j)=help
489          excessoro(i179+i,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
490        else
491          oro(i-i181,j)=help
492          excessoro(i-i181,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
493        endif
494      endif
495      if((isec1(6).eq.081).and.(isec1(7).eq.001)) then
496  ! LAND SEA MASK
497        help=zsec4(nxfield*(ny-j-1)+i+1)
498        if(i.le.i180) then
499          lsm(i179+i,j)=help
500        else
501          lsm(i-i181,j)=help
502        endif
503      endif
504      if((isec1(6).eq.221).and.(isec1(7).eq.001)) then
505  ! MIXING HEIGHT
506        help=zsec4(nxfield*(ny-j-1)+i+1)
507        if(i.le.i180) then
508          hmix(i179+i,j,1,n)=help
509        else
510          hmix(i-i181,j,1,n)=help
511        endif
512      endif
513      if((isec1(6).eq.052).and.(isec1(7).eq.105).and. &
514           (isec1(8).eq.02)) then
515  ! 2 M RELATIVE HUMIDITY
516        help=zsec4(nxfield*(ny-j-1)+i+1)
517        if(i.le.i180) then
518          qvh2(i179+i,j)=help
519        else
520          qvh2(i-i181,j)=help
521        endif
522      endif
523      if((isec1(6).eq.011).and.(isec1(7).eq.107)) then
524  ! TEMPERATURE LOWEST SIGMA LEVEL
525        help=zsec4(nxfield*(ny-j-1)+i+1)
526        if(i.le.i180) then
527          tlev1(i179+i,j)=help
528        else
529          tlev1(i-i181,j)=help
530        endif
531      endif
532      if((isec1(6).eq.033).and.(isec1(7).eq.107)) then
533  ! U VELOCITY LOWEST SIGMA LEVEL
534        help=zsec4(nxfield*(ny-j-1)+i+1)
535        if(i.le.i180) then
536          ulev1(i179+i,j)=help
537        else
538          ulev1(i-i181,j)=help
539        endif
540      endif
541      if((isec1(6).eq.034).and.(isec1(7).eq.107)) then
542  ! V VELOCITY LOWEST SIGMA LEVEL
543        help=zsec4(nxfield*(ny-j-1)+i+1)
544        if(i.le.i180) then
545          vlev1(i179+i,j)=help
546        else
547          vlev1(i-i181,j)=help
548        endif
549      endif
550
551    end do
552  end do
553
554  endif
555
556  if((isec1(6).eq.33).and.(isec1(7).eq.100)) then
557  ! NCEP ISOBARIC LEVELS
558    iumax=iumax+1
559  endif
560
561  call grib_release(igrib)
562  goto 10                      !! READ NEXT LEVEL OR PARAMETER
563  !
564  ! CLOSING OF INPUT DATA FILE
565  !
566
567  !HSO close grib file
56850   continue
569  call grib_close_file(ifile)
570
571  ! SENS. HEAT FLUX
572  sshf(:,:,1,n)=0.0     ! not available from gfs.tccz.pgrbfxx files
573  hflswitch=.false.    ! Heat flux not available
574  ! SOLAR RADIATIVE FLUXES
575  ssr(:,:,1,n)=0.0      ! not available from gfs.tccz.pgrbfxx files
576  ! EW SURFACE STRESS
577  ewss=0.0         ! not available from gfs.tccz.pgrbfxx files
578  ! NS SURFACE STRESS
579  nsss=0.0         ! not available from gfs.tccz.pgrbfxx files
580  strswitch=.false.    ! stress not available
581
582  ! CONVERT TP TO LSP (GRIB2 only)
583  if (gribVer.eq.2) then
584    do j=0,nymin1
585    do i=0,nxfield-1
586     if(i.le.i180) then
587     if (convprec(i179+i,j,1,n).lt.lsprec(i179+i,j,1,n)) then ! neg precip would occur
588         lsprec(i179+i,j,1,n)= &
589              lsprec(i179+i,j,1,n)-convprec(i179+i,j,1,n)
590     else
591         lsprec(i179+i,j,1,n)=0
592     endif
593     else
594     if (convprec(i-i181,j,1,n).lt.lsprec(i-i181,j,1,n)) then
595          lsprec(i-i181,j,1,n)= &
596               lsprec(i-i181,j,1,n)-convprec(i-i181,j,1,n)
597     else
598          lsprec(i-i181,j,1,n)=0
599     endif
600     endif
601    enddo
602    enddo
603  endif
604  !HSO end edits
605
606
607  ! TRANSFORM RH TO SPECIFIC HUMIDITY
608
609  do j=0,ny-1
610    do i=0,nxfield-1
611      do k=1,nuvz
612        help=qvh(i,j,k,n)
613        temp=tth(i,j,k,n)
614        plev1=akm(k)+bkm(k)*ps(i,j,1,n)
615        elev=ew(temp)*help/100.0
616        qvh(i,j,k,n)=xmwml*(elev/(plev1-((1.0-xmwml)*elev)))
617      end do
618    end do
619  end do
620
621  ! CALCULATE 2 M DEW POINT FROM 2 M RELATIVE HUMIDITY
622  ! USING BOLTON'S (1980) FORMULA
623  ! BECAUSE td2 IS NOT AVAILABLE FROM NCEP GFS DATA
624
625  do j=0,ny-1
626    do i=0,nxfield-1
627        help=qvh2(i,j)
628        temp=tt2(i,j,1,n)
629        elev=ew(temp)/100.*help/100.   !vapour pressure in hPa
630        td2(i,j,1,n)=243.5/(17.67/log(elev/6.112)-1)+273.
631        if (help.le.0.) td2(i,j,1,n)=tt2(i,j,1,n)
632    end do
633  end do
634
635  if(levdiff2.eq.0) then
636    iwmax=nlev_ec+1
637    do i=0,nxmin1
638      do j=0,nymin1
639        wwh(i,j,nlev_ec+1)=0.
640      end do
641    end do
642  endif
643
644
645  ! For global fields, assign the leftmost data column also to the rightmost
646  ! data column; if required, shift whole grid by nxshift grid points
647  !*************************************************************************
648
649  if (xglobal) then
650    call shift_field_0(ewss,nxfield,ny)
651    call shift_field_0(nsss,nxfield,ny)
652    call shift_field_0(oro,nxfield,ny)
653    call shift_field_0(excessoro,nxfield,ny)
654    call shift_field_0(lsm,nxfield,ny)
655    call shift_field_0(ulev1,nxfield,ny)
656    call shift_field_0(vlev1,nxfield,ny)
657    call shift_field_0(tlev1,nxfield,ny)
658    call shift_field_0(qvh2,nxfield,ny)
659    call shift_field(ps,nxfield,ny,1,1,2,n)
660    call shift_field(sd,nxfield,ny,1,1,2,n)
661    call shift_field(msl,nxfield,ny,1,1,2,n)
662    call shift_field(tcc,nxfield,ny,1,1,2,n)
663    call shift_field(u10,nxfield,ny,1,1,2,n)
664    call shift_field(v10,nxfield,ny,1,1,2,n)
665    call shift_field(tt2,nxfield,ny,1,1,2,n)
666    call shift_field(td2,nxfield,ny,1,1,2,n)
667    call shift_field(lsprec,nxfield,ny,1,1,2,n)
668    call shift_field(convprec,nxfield,ny,1,1,2,n)
669    call shift_field(sshf,nxfield,ny,1,1,2,n)
670    call shift_field(ssr,nxfield,ny,1,1,2,n)
671    call shift_field(hmix,nxfield,ny,1,1,2,n)
672    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
673    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
674    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
675    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
676    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
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