source: trunk/src/readwind_gfs_emos.f90 @ 38

Last change on this file since 38 was 4, checked in by mlanger, 11 years ago
  • Property svn:executable set to *
File size: 17.8 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(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  !             CAHENGE: 16/11/2005, Caroline Forster, GFS data         *
34  !                                                                     *
35  !**********************************************************************
36  !  Changes, Bernd C. Krueger, Feb. 2001:
37  !   Variables tth and qvh (on eta coordinates) in common block
38  !**********************************************************************
39  !                                                                     *
40  ! DESCRIPTION:                                                        *
41  !                                                                     *
42  ! READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE   *
43  ! INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE          *
44  !                                                                     *
45  ! INPUT:                                                              *
46  ! indj               indicates number of the wind field to be read in *
47  ! n                  temporal index for meteorological fields (1 to 3)*
48  !                                                                     *
49  ! IMPORTANT VARIABLES FROM COMMON BLOCK:                              *
50  !                                                                     *
51  ! wfname             File name of data to be read in                  *
52  ! nx,ny,nuvz,nwz     expected field dimensions                        *
53  ! nlev_ec            number of vertical levels ecmwf model            *
54  ! uu,vv,ww           wind fields                                      *
55  ! tt,qv              temperature and specific humidity                *
56  ! ps                 surface pressure                                 *
57  !                                                                     *
58  !**********************************************************************
59
60  use par_mod
61  use com_mod
62
63  implicit none
64
65  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
66  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
67  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
68  integer :: ii,indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,lunit
69
70  ! NCEP
71
72  integer :: numpt,numpu,numpv,numpw,numprh
73  real :: help, temp, ew
74  real :: elev
75  real :: ulev1(0:nxmax-1,0:nymax-1),vlev1(0:nxmax-1,0:nymax-1)
76  real :: tlev1(0:nxmax-1,0:nymax-1)
77  real :: qvh2(0:nxmax-1,0:nymax-1)
78
79  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
80
81  ! dimension of isec2 at least (22+n), where n is the number of parallels or
82  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
83
84  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
85  ! coordinate parameters
86
87  integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2)
88  integer :: isec4(64),inbuff(jpack),ilen,iswap,ierr,iword
89  real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp)
90  real :: xaux,yaux,xaux0,yaux0
91  real,parameter :: eps=1.e-4
92  real :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1)
93  real :: plev1,hlev1,ff10m,fflev1
94
95  character(len=1) :: yoper = 'D'
96  logical :: hflswitch,strswitch
97
98  hflswitch=.false.
99  strswitch=.false.
100  levdiff2=nlev_ec-nwz+1
101  iumax=0
102  iwmax=0
103
104  !
105  ! OPENING OF DATA FILE (GRIB CODE)
106  !
1075   call pbopen(lunit,path(3)(1:length(3))//wfname(indj),'r',ierr)
108  if(ierr.lt.0) goto 999
109
110  numpt=0
111  numpu=0
112  numpv=0
113  numpw=0
114  numprh=0
115  ifield=0
11610   ifield=ifield+1
117  !
118  ! GET NEXT FIELDS
119  !
120  call pbgrib(lunit,inbuff,jpack,ilen,ierr)
121  if(ierr.eq.-1) goto 50    ! EOF DETECTED
122  if(ierr.lt.-1) goto 888   ! ERROR DETECTED
123
124  ierr=1
125
126  ! Check whether we are on a little endian or on a big endian computer
127  !********************************************************************
128
129  !if (inbuff(1).eq.1112101447) then         ! little endian, swap bytes
130  !  iswap=1+ilen/4
131  !  call swap32(inbuff,iswap)
132  !else if (inbuff(1).ne.1196575042) then    ! big endian
133  !  stop 'subroutine gridcheck: corrupt GRIB data'
134  !endif
135
136
137  call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, &
138       zsec4,jpunp,inbuff,jpack,iword,yoper,ierr)
139  if (ierr.ne.0) goto 10   ! ERROR DETECTED
140
141  if(ifield.eq.1) then
142
143  ! CHECK GRID SPECIFICATIONS
144
145    if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
146    if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
147    xaux=real(isec2(5))/1000.+real(nxshift)*dx
148    if(xaux.eq.0.) xaux=-179.0     ! NCEP DATA
149    yaux=real(isec2(7))/1000.
150    xaux0=xlon0
151    yaux0=ylat0
152    if(xaux.lt.0.) xaux=xaux+360.
153    if(yaux.lt.0.) yaux=yaux+360.
154    if(xaux0.lt.0.) xaux0=xaux0+360.
155    if(yaux0.lt.0.) yaux0=yaux0+360.
156    if(abs(xaux-xaux0).gt.eps) &
157         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
158    if(abs(yaux-yaux0).gt.eps) &
159         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
160  endif
161
162  do j=0,nymin1
163    do i=0,nxfield-1
164      if((isec1(6).eq.011).and.(isec1(7).eq.100)) then
165  ! TEMPERATURE
166         if((i.eq.0).and.(j.eq.0)) then
167            do ii=1,nuvz
168              if ((isec1(8)*100.0).eq.akz(ii)) numpt=ii
169            end do
170        endif
171        help=zsec4(nxfield*(ny-j-1)+i+1)
172        if(i.le.180) then
173          tth(179+i,j,numpt,n)=help
174        else
175          tth(i-181,j,numpt,n)=help
176        endif
177      endif
178      if((isec1(6).eq.033).and.(isec1(7).eq.100)) then
179  ! U VELOCITY
180         if((i.eq.0).and.(j.eq.0)) then
181            do ii=1,nuvz
182              if ((isec1(8)*100.0).eq.akz(ii)) numpu=ii
183            end do
184        endif
185        help=zsec4(nxfield*(ny-j-1)+i+1)
186        if(i.le.180) then
187          uuh(179+i,j,numpu)=help
188        else
189          uuh(i-181,j,numpu)=help
190        endif
191      endif
192      if((isec1(6).eq.034).and.(isec1(7).eq.100)) then
193  ! V VELOCITY
194         if((i.eq.0).and.(j.eq.0)) then
195            do ii=1,nuvz
196              if ((isec1(8)*100.0).eq.akz(ii)) numpv=ii
197            end do
198        endif
199        help=zsec4(nxfield*(ny-j-1)+i+1)
200        if(i.le.180) then
201          vvh(179+i,j,numpv)=help
202        else
203          vvh(i-181,j,numpv)=help
204        endif
205      endif
206      if((isec1(6).eq.052).and.(isec1(7).eq.100)) then
207  ! RELATIVE HUMIDITY -> CONVERT TO SPECIFIC HUMIDITY LATER
208         if((i.eq.0).and.(j.eq.0)) then
209            do ii=1,nuvz
210              if ((isec1(8)*100.0).eq.akz(ii)) numprh=ii
211            end do
212        endif
213        help=zsec4(nxfield*(ny-j-1)+i+1)
214        if(i.le.180) then
215          qvh(179+i,j,numprh,n)=help
216        else
217          qvh(i-181,j,numprh,n)=help
218        endif
219      endif
220      if((isec1(6).eq.001).and.(isec1(7).eq.001)) then
221  ! SURFACE PRESSURE
222        help=zsec4(nxfield*(ny-j-1)+i+1)
223        if(i.le.180) then
224          ps(179+i,j,1,n)=help
225        else
226          ps(i-181,j,1,n)=help
227        endif
228      endif
229      if((isec1(6).eq.039).and.(isec1(7).eq.100)) then
230  ! W VELOCITY
231         if((i.eq.0).and.(j.eq.0)) then
232            do ii=1,nuvz
233              if ((isec1(8)*100.0).eq.akz(ii)) numpw=ii
234            end do
235        endif
236        help=zsec4(nxfield*(ny-j-1)+i+1)
237        if(i.le.180) then
238          wwh(179+i,j,numpw)=help
239        else
240          wwh(i-181,j,numpw)=help
241        endif
242      endif
243      if((isec1(6).eq.066).and.(isec1(7).eq.001)) then
244  ! SNOW DEPTH
245        help=zsec4(nxfield*(ny-j-1)+i+1)
246        if(i.le.180) then
247          sd(179+i,j,1,n)=help
248        else
249          sd(i-181,j,1,n)=help
250        endif
251      endif
252      if((isec1(6).eq.002).and.(isec1(7).eq.102)) then
253  ! MEAN SEA LEVEL PRESSURE
254        help=zsec4(nxfield*(ny-j-1)+i+1)
255        if(i.le.180) then
256          msl(179+i,j,1,n)=help
257        else
258          msl(i-181,j,1,n)=help
259        endif
260      endif
261      if((isec1(6).eq.071).and.(isec1(7).eq.244)) then
262  ! TOTAL CLOUD COVER
263        help=zsec4(nxfield*(ny-j-1)+i+1)
264        if(i.le.180) then
265          tcc(179+i,j,1,n)=help
266        else
267          tcc(i-181,j,1,n)=help
268        endif
269      endif
270      if((isec1(6).eq.033).and.(isec1(7).eq.105).and. &
271           (isec1(8).eq.10)) then
272  ! 10 M U VELOCITY
273        help=zsec4(nxfield*(ny-j-1)+i+1)
274        if(i.le.180) then
275          u10(179+i,j,1,n)=help
276        else
277          u10(i-181,j,1,n)=help
278        endif
279      endif
280      if((isec1(6).eq.034).and.(isec1(7).eq.105).and. &
281           (isec1(8).eq.10)) then
282  ! 10 M V VELOCITY
283        help=zsec4(nxfield*(ny-j-1)+i+1)
284        if(i.le.180) then
285          v10(179+i,j,1,n)=help
286        else
287          v10(i-181,j,1,n)=help
288        endif
289      endif
290      if((isec1(6).eq.011).and.(isec1(7).eq.105).and. &
291           (isec1(8).eq.02)) then
292  ! 2 M TEMPERATURE
293        help=zsec4(nxfield*(ny-j-1)+i+1)
294        if(i.le.180) then
295          tt2(179+i,j,1,n)=help
296        else
297          tt2(i-181,j,1,n)=help
298        endif
299      endif
300      if((isec1(6).eq.017).and.(isec1(7).eq.105).and. &
301           (isec1(8).eq.02)) then
302  ! 2 M DEW POINT TEMPERATURE
303        help=zsec4(nxfield*(ny-j-1)+i+1)
304        if(i.le.180) then
305          td2(179+i,j,1,n)=help
306        else
307          td2(i-181,j,1,n)=help
308        endif
309      endif
310      if((isec1(6).eq.062).and.(isec1(7).eq.001)) then
311  ! LARGE SCALE PREC.
312        help=zsec4(nxfield*(ny-j-1)+i+1)
313        if(i.le.180) then
314          lsprec(179+i,j,1,n)=help
315        else
316          lsprec(i-181,j,1,n)=help
317        endif
318      endif
319      if((isec1(6).eq.063).and.(isec1(7).eq.001)) then
320  ! CONVECTIVE PREC.
321        help=zsec4(nxfield*(ny-j-1)+i+1)
322        if(i.le.180) then
323          convprec(179+i,j,1,n)=help
324        else
325          convprec(i-181,j,1,n)=help
326        endif
327      endif
328  ! SENS. HEAT FLUX
329      sshf(i,j,1,n)=0.0     ! not available from gfs.tccz.pgrbfxx files
330      hflswitch=.false.    ! Heat flux not available
331  ! SOLAR RADIATIVE FLUXES
332      ssr(i,j,1,n)=0.0      ! not available from gfs.tccz.pgrbfxx files
333  ! EW SURFACE STRESS
334      ewss(i,j)=0.0         ! not available from gfs.tccz.pgrbfxx files
335  ! NS SURFACE STRESS
336      nsss(i,j)=0.0         ! not available from gfs.tccz.pgrbfxx files
337      strswitch=.false.    ! stress not available
338      if((isec1(6).eq.007).and.(isec1(7).eq.001)) then
339  ! TOPOGRAPHY
340        help=zsec4(nxfield*(ny-j-1)+i+1)
341        if(i.le.180) then
342          oro(179+i,j)=help
343          excessoro(179+i,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
344        else
345          oro(i-181,j)=help
346          excessoro(i-181,j)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
347        endif
348      endif
349      if((isec1(6).eq.081).and.(isec1(7).eq.001)) then
350  ! LAND SEA MASK
351        help=zsec4(nxfield*(ny-j-1)+i+1)
352        if(i.le.180) then
353          lsm(179+i,j)=help
354        else
355          lsm(i-181,j)=help
356        endif
357      endif
358      if((isec1(6).eq.221).and.(isec1(7).eq.001)) then
359  ! MIXING HEIGHT
360        help=zsec4(nxfield*(ny-j-1)+i+1)
361        if(i.le.180) then
362          hmix(179+i,j,1,n)=help
363        else
364          hmix(i-181,j,1,n)=help
365        endif
366      endif
367      if((isec1(6).eq.052).and.(isec1(7).eq.105).and. &
368           (isec1(8).eq.02)) then
369  ! 2 M RELATIVE HUMIDITY
370        help=zsec4(nxfield*(ny-j-1)+i+1)
371        if(i.le.180) then
372          qvh2(179+i,j)=help
373        else
374          qvh2(i-181,j)=help
375        endif
376      endif
377      if((isec1(6).eq.011).and.(isec1(7).eq.107)) then
378  ! TEMPERATURE LOWEST SIGMA LEVEL
379        help=zsec4(nxfield*(ny-j-1)+i+1)
380        if(i.le.180) then
381          tlev1(179+i,j)=help
382        else
383          tlev1(i-181,j)=help
384        endif
385      endif
386      if((isec1(6).eq.033).and.(isec1(7).eq.107)) then
387  ! U VELOCITY LOWEST SIGMA LEVEL
388        help=zsec4(nxfield*(ny-j-1)+i+1)
389        if(i.le.180) then
390          ulev1(179+i,j)=help
391        else
392          ulev1(i-181,j)=help
393        endif
394      endif
395      if((isec1(6).eq.034).and.(isec1(7).eq.107)) then
396  ! V VELOCITY LOWEST SIGMA LEVEL
397        help=zsec4(nxfield*(ny-j-1)+i+1)
398        if(i.le.180) then
399          vlev1(179+i,j)=help
400        else
401          vlev1(i-181,j)=help
402        endif
403      endif
404
405    end do
406  end do
407
408  if((isec1(6).eq.33).and.(isec1(7).eq.100)) then
409  ! NCEP ISOBARIC LEVELS
410    iumax=iumax+1
411  endif
412
413
414  goto 10                      !! READ NEXT LEVEL OR PARAMETER
415  !
416  ! CLOSING OF INPUT DATA FILE
417  !
41850   call pbclose(lunit,ierr)     !! FINNISHED READING / CLOSING GRIB FILE
419
420  ! TRANSFORM RH TO SPECIFIC HUMIDITY
421
422  do j=0,ny-1
423    do i=0,nxfield-1
424      do k=1,nuvz
425        help=qvh(i,j,k,n)
426        temp=tth(i,j,k,n)
427        plev1=akm(k)+bkm(k)*ps(i,j,1,n)
428        elev=ew(temp)*help/100.0
429        qvh(i,j,k,n)=xmwml*(elev/(plev1-((1.0-xmwml)*elev)))
430      end do
431    end do
432  end do
433
434  ! CALCULATE 2 M DEW POINT FROM 2 M RELATIVE HUMIDITY
435  ! USING BOLTON'S (1980) FORMULA
436  ! BECAUSE td2 IS NOT AVAILABLE FROM NCEP GFS DATA
437
438  do j=0,ny-1
439    do i=0,nxfield-1
440        help=qvh2(i,j)
441        temp=tt2(i,j,1,n)
442        elev=ew(temp)/100.*help/100.   !vapour pressure in hPa
443        td2(i,j,1,n)=243.5/(17.67/log(elev/6.112)-1)+273.
444        if (help.le.0.) td2(i,j,1,n)=tt2(i,j,1,n)
445    end do
446  end do
447
448  if(levdiff2.eq.0) then
449    iwmax=nlev_ec+1
450    do i=0,nxmin1
451      do j=0,nymin1
452        wwh(i,j,nlev_ec+1)=0.
453      end do
454    end do
455  endif
456
457
458  ! For global fields, assign the leftmost data column also to the rightmost
459  ! data column; if required, shift whole grid by nxshift grid points
460  !*************************************************************************
461
462  if (xglobal) then
463    call shift_field_0(ewss,nxfield,ny)
464    call shift_field_0(nsss,nxfield,ny)
465    call shift_field_0(oro,nxfield,ny)
466    call shift_field_0(excessoro,nxfield,ny)
467    call shift_field_0(lsm,nxfield,ny)
468    call shift_field_0(ulev1,nxfield,ny)
469    call shift_field_0(vlev1,nxfield,ny)
470    call shift_field_0(tlev1,nxfield,ny)
471    call shift_field_0(qvh2,nxfield,ny)
472    call shift_field(ps,nxfield,ny,1,1,2,n)
473    call shift_field(sd,nxfield,ny,1,1,2,n)
474    call shift_field(msl,nxfield,ny,1,1,2,n)
475    call shift_field(tcc,nxfield,ny,1,1,2,n)
476    call shift_field(u10,nxfield,ny,1,1,2,n)
477    call shift_field(v10,nxfield,ny,1,1,2,n)
478    call shift_field(tt2,nxfield,ny,1,1,2,n)
479    call shift_field(td2,nxfield,ny,1,1,2,n)
480    call shift_field(lsprec,nxfield,ny,1,1,2,n)
481    call shift_field(convprec,nxfield,ny,1,1,2,n)
482    call shift_field(sshf,nxfield,ny,1,1,2,n)
483    call shift_field(ssr,nxfield,ny,1,1,2,n)
484    call shift_field(hmix,nxfield,ny,1,1,2,n)
485    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
486    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
487    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
488    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
489    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
490  endif
491
492  do i=0,nxmin1
493    do j=0,nymin1
494      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
495    end do
496  end do
497
498  if ((.not.hflswitch).or.(.not.strswitch)) then
499  !  write(*,*) 'WARNING: No flux data contained in GRIB file ',
500  !    +  wfname(indj)
501
502  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
503  !***************************************************************************
504
505    do i=0,nxmin1
506      do j=0,nymin1
507        hlev1=30.0                     ! HEIGHT OF FIRST MODEL SIGMA LAYER
508        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
509        fflev1=sqrt(ulev1(i,j)**2+vlev1(i,j)**2)
510        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
511             tt2(i,j,1,n),tlev1(i,j),ff10m,fflev1, &
512             surfstr(i,j,1,n),sshf(i,j,1,n))
513        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
514        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
515      end do
516    end do
517  endif
518
519
520  if(iumax.ne.nuvz) stop 'READWIND: NUVZ NOT CONSISTENT'
521  if(iumax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
522
523  return
524888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
525  write(*,*) ' #### ',wfname(indj),'                    #### '
526  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
527  stop 'Execution terminated'
528999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
529  write(*,*) ' #### ',wfname(indj),'                    #### '
530  write(*,*) ' #### CANNOT BE OPENED !!!                    #### '
531  stop 'Execution terminated'
532
533end subroutine readwind
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG