source: trunk/src/readwind_emos.f90 @ 28

Last change on this file since 28 was 4, checked in by mlanger, 11 years ago
  • Property svn:executable set to *
File size: 13.4 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  !                                                                     *
34  !**********************************************************************
35  !  Changes, Bernd C. Krueger, Feb. 2001:
36  !   Variables tth and qvh (on eta coordinates) in common block
37  !**********************************************************************
38  !                                                                     *
39  ! DESCRIPTION:                                                        *
40  !                                                                     *
41  ! READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE   *
42  ! INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE          *
43  !                                                                     *
44  ! INPUT:                                                              *
45  ! indj               indicates number of the wind field to be read in *
46  ! n                  temporal index for meteorological fields (1 to 3)*
47  !                                                                     *
48  ! IMPORTANT VARIABLES FROM COMMON BLOCK:                              *
49  !                                                                     *
50  ! wfname             File name of data to be read in                  *
51  ! nx,ny,nuvz,nwz     expected field dimensions                        *
52  ! nlev_ec            number of vertical levels ecmwf model            *
53  ! uu,vv,ww           wind fields                                      *
54  ! tt,qv              temperature and specific humidity                *
55  ! ps                 surface pressure                                 *
56  !                                                                     *
57  !**********************************************************************
58
59  use par_mod
60  use com_mod
61
62  implicit none
63
64  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
65  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
66  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
67  integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,lunit
68
69  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
70
71  ! dimension of isec2 at least (22+n), where n is the number of parallels or
72  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
73
74  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
75  ! coordinate parameters
76
77  integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2)
78  integer :: isec4(64),inbuff(jpack),ilen,ierr,iword
79  !integer iswap
80  real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp)
81  real :: xaux,yaux,xaux0,yaux0
82  real,parameter :: eps=1.e-4
83  real :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1)
84  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
85
86  character(len=1) :: yoper = 'D'
87  logical :: hflswitch,strswitch
88
89  hflswitch=.false.
90  strswitch=.false.
91  levdiff2=nlev_ec-nwz+1
92  iumax=0
93  iwmax=0
94
95  !
96  ! OPENING OF DATA FILE (GRIB CODE)
97  !
985   call pbopen(lunit,path(3)(1:length(3))//wfname(indj),'r',ierr)
99  if(ierr.lt.0) goto 999
100
101  ifield=0
10210   ifield=ifield+1
103  !
104  ! GET NEXT FIELDS
105  !
106  call pbgrib(lunit,inbuff,jpack,ilen,ierr)
107  if(ierr.eq.-1) goto 50    ! EOF DETECTED
108  if(ierr.lt.-1) goto 888   ! ERROR DETECTED
109
110  ierr=1
111
112  ! Check whether we are on a little endian or on a big endian computer
113  !********************************************************************
114
115  !if (inbuff(1).eq.1112101447) then         ! little endian, swap bytes
116  !  iswap=1+ilen/4
117  !  call swap32(inbuff,iswap)
118  !else if (inbuff(1).ne.1196575042) then    ! big endian
119  !  stop 'subroutine gridcheck: corrupt GRIB data'
120  !endif
121
122
123  call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, &
124       zsec4,jpunp,inbuff,jpack,iword,yoper,ierr)
125  if (ierr.ne.0) goto 888   ! ERROR DETECTED
126
127  if(ifield.eq.1) then
128
129  ! CHECK GRID SPECIFICATIONS
130
131    if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
132    if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
133    if(isec2(12)/2-1.ne.nlev_ec) &
134         stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
135    xaux=real(isec2(5))/1000.+real(nxshift)*dx
136    yaux=real(isec2(7))/1000.
137    xaux0=xlon0
138    yaux0=ylat0
139    if(xaux.lt.0.) xaux=xaux+360.
140    if(yaux.lt.0.) yaux=yaux+360.
141    if(xaux0.lt.0.) xaux0=xaux0+360.
142    if(yaux0.lt.0.) yaux0=yaux0+360.
143    if(abs(xaux-xaux0).gt.eps) &
144         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
145    if(abs(yaux-yaux0).gt.eps) &
146         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
147  endif
148
149  do j=0,nymin1
150    do i=0,nxfield-1
151      k=isec1(8)
152      if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
153           zsec4(nxfield*(ny-j-1)+i+1)
154      if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
155           zsec4(nxfield*(ny-j-1)+i+1)
156      if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
157           zsec4(nxfield*(ny-j-1)+i+1)
158      if(isec1(6).eq.133) then                      !! SPEC. HUMIDITY
159        qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
160        if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
161             qvh(i,j,nlev_ec-k+2,n) = 0.
162  !        this is necessary because the gridded data may contain
163  !        spurious negative values
164      endif
165      if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
166           zsec4(nxfield*(ny-j-1)+i+1)
167
168      if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
169           zsec4(nxfield*(ny-j-1)+i+1)
170      if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
171           zsec4(nxfield*(ny-j-1)+i+1)
172      if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
173           zsec4(nxfield*(ny-j-1)+i+1)
174      if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
175           zsec4(nxfield*(ny-j-1)+i+1)
176      if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
177           zsec4(nxfield*(ny-j-1)+i+1)
178      if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
179           zsec4(nxfield*(ny-j-1)+i+1)
180      if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
181           zsec4(nxfield*(ny-j-1)+i+1)
182      if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
183           zsec4(nxfield*(ny-j-1)+i+1)
184      if(isec1(6).eq.142) then                      !! LARGE SCALE PREC.
185        lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
186        if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
187      endif
188      if(isec1(6).eq.143) then                      !! CONVECTIVE PREC.
189        convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
190        if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
191      endif
192      if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
193           zsec4(nxfield*(ny-j-1)+i+1)
194      if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
195           hflswitch=.true.    ! Heat flux available
196      if(isec1(6).eq.176) then                      !! SOLAR RADIATION
197        ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
198        if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
199      endif
200      if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
201           zsec4(nxfield*(ny-j-1)+i+1)
202      if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
203           zsec4(nxfield*(ny-j-1)+i+1)
204      if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
205           (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
206      if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
207           zsec4(nxfield*(ny-j-1)+i+1)/ga
208      if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
209           zsec4(nxfield*(ny-j-1)+i+1)
210      if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
211           zsec4(nxfield*(ny-j-1)+i+1)
212      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
213      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
214    end do
215  end do
216
217  goto 10                      !! READ NEXT LEVEL OR PARAMETER
218  !
219  ! CLOSING OF INPUT DATA FILE
220  !
22150   call pbclose(lunit,ierr)     !! FINNISHED READING / CLOSING GRIB FILE
222
223  if(levdiff2.eq.0) then
224    iwmax=nlev_ec+1
225    do i=0,nxmin1
226      do j=0,nymin1
227        wwh(i,j,nlev_ec+1)=0.
228      end do
229    end do
230  endif
231
232
233  ! For global fields, assign the leftmost data column also to the rightmost
234  ! data column; if required, shift whole grid by nxshift grid points
235  !*************************************************************************
236
237  if (xglobal) then
238    call shift_field_0(ewss,nxfield,ny)
239    call shift_field_0(nsss,nxfield,ny)
240    call shift_field_0(oro,nxfield,ny)
241    call shift_field_0(excessoro,nxfield,ny)
242    call shift_field_0(lsm,nxfield,ny)
243    call shift_field(ps,nxfield,ny,1,1,2,n)
244    call shift_field(sd,nxfield,ny,1,1,2,n)
245    call shift_field(msl,nxfield,ny,1,1,2,n)
246    call shift_field(tcc,nxfield,ny,1,1,2,n)
247    call shift_field(u10,nxfield,ny,1,1,2,n)
248    call shift_field(v10,nxfield,ny,1,1,2,n)
249    call shift_field(tt2,nxfield,ny,1,1,2,n)
250    call shift_field(td2,nxfield,ny,1,1,2,n)
251    call shift_field(lsprec,nxfield,ny,1,1,2,n)
252    call shift_field(convprec,nxfield,ny,1,1,2,n)
253    call shift_field(sshf,nxfield,ny,1,1,2,n)
254    call shift_field(ssr,nxfield,ny,1,1,2,n)
255    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
256    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
257    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
258    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
259    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
260  endif
261
262  do i=0,nxmin1
263    do j=0,nymin1
264      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
265    end do
266  end do
267
268  if ((.not.hflswitch).or.(.not.strswitch)) then
269    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
270         wfname(indj)
271
272  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
273  ! As ECMWF has increased the model resolution, such that now the first model
274  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
275  ! (3rd model level in FLEXPART) for the profile method
276  !***************************************************************************
277
278    do i=0,nxmin1
279      do j=0,nymin1
280        plev1=akz(3)+bkz(3)*ps(i,j,1,n)
281        pmean=0.5*(ps(i,j,1,n)+plev1)
282        tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
283        fu=-r_air*tv/ga/pmean
284        hlev1=fu*(plev1-ps(i,j,1,n))   ! HEIGTH OF FIRST MODEL LAYER
285        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
286        fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
287        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
288             tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
289             surfstr(i,j,1,n),sshf(i,j,1,n))
290        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
291        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
292      end do
293    end do
294  endif
295
296
297  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
298  ! level at the ground
299  ! Specific humidity is taken the same as at one level above
300  ! Temperature is taken as 2 m temperature
301  !**************************************************************************
302
303     do i=0,nxmin1
304        do j=0,nymin1
305           uuh(i,j,1)=u10(i,j,1,n)
306           vvh(i,j,1)=v10(i,j,1,n)
307           qvh(i,j,1,n)=qvh(i,j,2,n)
308           tth(i,j,1,n)=tt2(i,j,1,n)
309        end do
310     end do
311
312  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
313  if(iwmax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
314
315  return
316888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
317  write(*,*) ' #### ',wfname(indj),'                    #### '
318  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
319  stop 'Execution terminated'
320999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
321  write(*,*) ' #### ',wfname(indj),'                    #### '
322  write(*,*) ' #### CANNOT BE OPENED !!!                    #### '
323  stop 'Execution terminated'
324
325end subroutine readwind
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG