source: trunk/src/readwind_nests.f90 @ 20

Last change on this file since 20 was 4, checked in by mlanger, 11 years ago
File size: 16.9 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_nests(indj,n,uuhn,vvhn,wwhn)
23  !                           i   i  o    o    o
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads the wind fields for the nested model domains.       *
27  !     It is similar to subroutine readwind, which reads the mother domain.   *
28  !                                                                            *
29  !     Authors: A. Stohl, G. Wotawa                                           *
30  !                                                                            *
31  !     8 February 1999                                                        *
32  !                                                                            *
33  !     Last update: 17 October 2000, A. Stohl                                 *
34  !                                                                            *
35  !*****************************************************************************
36  !  Changes, Bernd C. Krueger, Feb. 2001:                                     *
37  !        Variables tthn and qvhn (on eta coordinates) in common block        *
38  !  CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with ECMWF grib_api    *
39  !  CHANGE: 03/12/2008, Harald Sodemann, update to f90 with ECMWF grib_api    *
40  !*****************************************************************************
41
42  use grib_api
43  use par_mod
44  use com_mod
45
46  implicit none
47
48  !HSO  parameters for grib_api
49  integer :: ifile
50  integer :: iret
51  integer :: igrib
52  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
53  integer :: gotGrid
54  !HSO  end
55
56  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
57  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
58  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
59  integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,l
60
61  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
62
63  ! dimension of isec2 at least (22+n), where n is the number of parallels or
64  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
65
66  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
67  ! coordinate parameters
68
69  integer :: isec1(56),isec2(22+nxmaxn+nymaxn)
70  real(kind=4) :: zsec4(jpunp)
71  real(kind=8) :: xauxin,yauxin
72  real(kind=4) :: xaux,yaux,xaux0,yaux0
73  real :: ewss(0:nxmaxn-1,0:nymaxn-1),nsss(0:nxmaxn-1,0:nymaxn-1)
74  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
75
76  logical :: hflswitch,strswitch
77
78  !HSO  grib api error messages
79  character(len=24) :: gribErrorMsg = 'Error reading grib file'
80  character(len=20) :: gribFunction = 'readwind_nests'
81
82  do l=1,numbnests
83    hflswitch=.false.
84    strswitch=.false.
85    levdiff2=nlev_ec-nwz+1
86    iumax=0
87    iwmax=0
88
89    ifile=0
90    igrib=0
91    iret=0
92
93  !
94  ! OPENING OF DATA FILE (GRIB CODE)
95  !
96
975   call grib_open_file(ifile,path(numpath+2*(l-1)+1) &
98         (1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,indj)),'r')
99  if (iret.ne.GRIB_SUCCESS) then
100    goto 888   ! ERROR DETECTED
101  endif
102  !turn on support for multi fields messages */
103  !call grib_multi_support_on
104
105    gotGrid=0
106    ifield=0
10710   ifield=ifield+1
108  !
109  ! GET NEXT FIELDS
110  !
111  call grib_new_from_file(ifile,igrib,iret)
112  if (iret.eq.GRIB_END_OF_FILE)  then
113    goto 50    ! EOF DETECTED
114  elseif (iret.ne.GRIB_SUCCESS) then
115    goto 888   ! ERROR DETECTED
116  endif
117
118  !first see if we read GRIB1 or GRIB2
119  call grib_get_int(igrib,'editionNumber',gribVer,iret)
120  call grib_check(iret,gribFunction,gribErrorMsg)
121
122  if (gribVer.eq.1) then ! GRIB Edition 1
123
124  !print*,'GRiB Edition 1'
125  !read the grib2 identifiers
126  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
127  call grib_check(iret,gribFunction,gribErrorMsg)
128  call grib_get_int(igrib,'level',isec1(8),iret)
129  call grib_check(iret,gribFunction,gribErrorMsg)
130
131  !change code for etadot to code for omega
132  if (isec1(6).eq.77) then
133    isec1(6)=135
134  endif
135
136  else
137
138  !print*,'GRiB Edition 2'
139  !read the grib2 identifiers
140  call grib_get_int(igrib,'discipline',discipl,iret)
141  call grib_check(iret,gribFunction,gribErrorMsg)
142  call grib_get_int(igrib,'parameterCategory',parCat,iret)
143  call grib_check(iret,gribFunction,gribErrorMsg)
144  call grib_get_int(igrib,'parameterNumber',parNum,iret)
145  call grib_check(iret,gribFunction,gribErrorMsg)
146  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
147  call grib_check(iret,gribFunction,gribErrorMsg)
148  call grib_get_int(igrib,'level',valSurf,iret)
149  call grib_check(iret,gribFunction,gribErrorMsg)
150
151  !print*,discipl,parCat,parNum,typSurf,valSurf
152
153  !convert to grib1 identifiers
154  isec1(6)=-1
155  isec1(7)=-1
156  isec1(8)=-1
157  isec1(8)=valSurf     ! level
158  if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
159    isec1(6)=130         ! indicatorOfParameter
160  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
161    isec1(6)=131         ! indicatorOfParameter
162  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
163    isec1(6)=132         ! indicatorOfParameter
164  elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
165    isec1(6)=133         ! indicatorOfParameter
166  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
167    isec1(6)=134         ! indicatorOfParameter
168  elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
169    isec1(6)=135         ! indicatorOfParameter
170  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
171    isec1(6)=151         ! indicatorOfParameter
172  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
173    isec1(6)=165         ! indicatorOfParameter
174  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
175    isec1(6)=166         ! indicatorOfParameter
176  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
177    isec1(6)=167         ! indicatorOfParameter
178  elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
179    isec1(6)=168         ! indicatorOfParameter
180  elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
181    isec1(6)=141         ! indicatorOfParameter
182  elseif ((parCat.eq.6).and.(parNum.eq.1)) then ! CC
183    isec1(6)=164         ! indicatorOfParameter
184  elseif ((parCat.eq.1).and.(parNum.eq.9)) then ! LSP
185    isec1(6)=142         ! indicatorOfParameter
186  elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
187    isec1(6)=143         ! indicatorOfParameter
188  elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
189    isec1(6)=146         ! indicatorOfParameter
190  elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
191    isec1(6)=176         ! indicatorOfParameter
192  elseif ((parCat.eq.2).and.(parNum.eq.17)) then ! EWSS
193    isec1(6)=180         ! indicatorOfParameter
194  elseif ((parCat.eq.2).and.(parNum.eq.18)) then ! NSSS
195    isec1(6)=181         ! indicatorOfParameter
196  elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
197    isec1(6)=129         ! indicatorOfParameter
198  elseif ((parCat.eq.3).and.(parNum.eq.7)) then ! SDO
199    isec1(6)=160         ! indicatorOfParameter
200  elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
201       (typSurf.eq.1)) then ! LSM
202    isec1(6)=172         ! indicatorOfParameter
203  else
204    print*,'***ERROR: undefined GRiB2 message found!',discipl, &
205         parCat,parNum,typSurf
206  endif
207
208  endif
209
210  !HSO  get the size and data of the values array
211  if (isec1(6).ne.-1) then
212    call grib_get_real4_array(igrib,'values',zsec4,iret)
213    call grib_check(iret,gribFunction,gribErrorMsg)
214  endif
215
216  !HSO  get the required fields from section 2 in a gribex compatible manner
217  if(ifield.eq.1) then
218  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
219       isec2(2),iret)
220  call grib_check(iret,gribFunction,gribErrorMsg)
221  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
222       isec2(3),iret)
223  call grib_check(iret,gribFunction,gribErrorMsg)
224  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
225       isec2(12))
226  call grib_check(iret,gribFunction,gribErrorMsg)
227  ! CHECK GRID SPECIFICATIONS
228  if(isec2(2).ne.nxn(l)) stop &
229       'READWIND: NX NOT CONSISTENT FOR A NESTING LEVEL'
230  if(isec2(3).ne.nyn(l)) stop &
231       'READWIND: NY NOT CONSISTENT FOR A NESTING LEVEL'
232  if(isec2(12)/2-1.ne.nlev_ec) stop 'READWIND: VERTICAL DISCRET&
233       &IZATION NOT CONSISTENT FOR A NESTING LEVEL'
234  endif ! ifield
235
236  !HSO  get the second part of the grid dimensions only from GRiB1 messages
237  if ((gribVer.eq.1).and.(gotGrid.eq.0)) then
238    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
239         xauxin,iret)
240    call grib_check(iret,gribFunction,gribErrorMsg)
241    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
242         yauxin,iret)
243    call grib_check(iret,gribFunction,gribErrorMsg)
244    xaux=xauxin
245    yaux=yauxin
246    xaux0=xlon0n(l)
247    yaux0=ylat0n(l)
248    if(xaux.lt.0.) xaux=xaux+360.
249    if(yaux.lt.0.) yaux=yaux+360.
250    if(xaux0.lt.0.) xaux0=xaux0+360.
251    if(yaux0.lt.0.) yaux0=yaux0+360.
252    if(xaux.ne.xaux0) &
253         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT FOR A NES&
254         &TING LEVEL'
255    if(yaux.ne.yaux0) &
256         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT FOR A NEST&
257         &ING LEVEL'
258    gotGrid=1
259  endif
260
261    do j=0,nyn(l)-1
262      do i=0,nxn(l)-1
263        k=isec1(8)
264        if(isec1(6).eq.130) tthn(i,j,nlev_ec-k+2,n,l)= &!! TEMPERATURE
265             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
266        if(isec1(6).eq.131) uuhn(i,j,nlev_ec-k+2,l)= &!! U VELOCITY
267             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
268        if(isec1(6).eq.132) vvhn(i,j,nlev_ec-k+2,l)= &!! V VELOCITY
269             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
270        if(isec1(6).eq.133) then                         !! SPEC. HUMIDITY
271          qvhn(i,j,nlev_ec-k+2,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
272          if (qvhn(i,j,nlev_ec-k+2,n,l) .lt. 0.) &
273               qvhn(i,j,nlev_ec-k+2,n,l) = 0.
274  !          this is necessary because the gridded data may contain
275  !          spurious negative values
276        endif
277        if(isec1(6).eq.134) psn(i,j,1,n,l)= &!! SURF. PRESS.
278             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
279        if(isec1(6).eq.135) wwhn(i,j,nlev_ec-k+1,l)= &!! W VELOCITY
280             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
281        if(isec1(6).eq.141) sdn(i,j,1,n,l)= &!! SNOW DEPTH
282             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
283        if(isec1(6).eq.151) msln(i,j,1,n,l)= &!! SEA LEVEL PRESS.
284             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
285        if(isec1(6).eq.164) tccn(i,j,1,n,l)= &!! CLOUD COVER
286             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
287        if(isec1(6).eq.165) u10n(i,j,1,n,l)= &!! 10 M U VELOCITY
288             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
289        if(isec1(6).eq.166) v10n(i,j,1,n,l)= &!! 10 M V VELOCITY
290             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
291        if(isec1(6).eq.167) tt2n(i,j,1,n,l)= &!! 2 M TEMPERATURE
292             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
293        if(isec1(6).eq.168) td2n(i,j,1,n,l)= &!! 2 M DEW POINT
294             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
295        if(isec1(6).eq.142) then                         !! LARGE SCALE PREC.
296          lsprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
297          if (lsprecn(i,j,1,n,l).lt.0.) lsprecn(i,j,1,n,l)=0.
298        endif
299        if(isec1(6).eq.143) then                         !! CONVECTIVE PREC.
300          convprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
301          if (convprecn(i,j,1,n,l).lt.0.) convprecn(i,j,1,n,l)=0.
302        endif
303        if(isec1(6).eq.146) sshfn(i,j,1,n,l)= &!! SENS. HEAT FLUX
304             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
305        if((isec1(6).eq.146).and. &
306             (zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) hflswitch=.true.    ! Heat flux available
307        if(isec1(6).eq.176) then                         !! SOLAR RADIATION
308          ssrn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
309          if (ssrn(i,j,1,n,l).lt.0.) ssrn(i,j,1,n,l)=0.
310        endif
311        if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
312             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
313        if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
314             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
315        if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
316             (zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
317        if(isec1(6).eq.129) oron(i,j,l)= &!! ECMWF OROGRAPHY
318             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
319        if(isec1(6).eq.160) excessoron(i,j,l)= &!! STANDARD DEVIATION OF OROGRAPHY
320             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
321        if(isec1(6).eq.172) lsmn(i,j,l)= &!! ECMWF LAND SEA MASK
322             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
323        if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
324        if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
325
326      end do
327    end do
328
329  call grib_release(igrib)
330  goto 10                      !! READ NEXT LEVEL OR PARAMETER
331  !
332  ! CLOSING OF INPUT DATA FILE
333  !
33450   call grib_close_file(ifile)
335
336  !error message if no fields found with correct first longitude in it
337  if (gotGrid.eq.0) then
338    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
339         'messages'
340    stop
341  endif
342
343  if(levdiff2.eq.0) then
344    iwmax=nlev_ec+1
345    do i=0,nxn(l)-1
346      do j=0,nyn(l)-1
347        wwhn(i,j,nlev_ec+1,l)=0.
348      end do
349    end do
350  endif
351
352  do i=0,nxn(l)-1
353    do j=0,nyn(l)-1
354      surfstrn(i,j,1,n,l)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
355    end do
356  end do
357
358  if ((.not.hflswitch).or.(.not.strswitch)) then
359    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
360         wfnamen(l,indj)
361
362  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
363  ! As ECMWF has increased the model resolution, such that now the first model
364  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
365  ! (3rd model level in FLEXPART) for the profile method
366  !***************************************************************************
367
368    do i=0,nxn(l)-1
369      do j=0,nyn(l)-1
370        plev1=akz(3)+bkz(3)*psn(i,j,1,n,l)
371        pmean=0.5*(psn(i,j,1,n,l)+plev1)
372        tv=tthn(i,j,3,n,l)*(1.+0.61*qvhn(i,j,3,n,l))
373        fu=-r_air*tv/ga/pmean
374        hlev1=fu*(plev1-psn(i,j,1,n,l))   ! HEIGTH OF FIRST MODEL LAYER
375        ff10m= sqrt(u10n(i,j,1,n,l)**2+v10n(i,j,1,n,l)**2)
376        fflev1=sqrt(uuhn(i,j,3,l)**2+vvhn(i,j,3,l)**2)
377        call pbl_profile(psn(i,j,1,n,l),td2n(i,j,1,n,l),hlev1, &
378             tt2n(i,j,1,n,l),tthn(i,j,3,n,l),ff10m,fflev1, &
379             surfstrn(i,j,1,n,l),sshfn(i,j,1,n,l))
380        if(sshfn(i,j,1,n,l).gt.200.) sshfn(i,j,1,n,l)=200.
381        if(sshfn(i,j,1,n,l).lt.-400.) sshfn(i,j,1,n,l)=-400.
382      end do
383    end do
384  endif
385
386
387  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
388  ! level at the ground
389  ! Specific humidity is taken the same as at one level above
390  ! Temperature is taken as 2 m temperature
391  !**************************************************************************
392
393    do i=0,nxn(l)-1
394      do j=0,nyn(l)-1
395        uuhn(i,j,1,l)=u10n(i,j,1,n,l)
396        vvhn(i,j,1,l)=v10n(i,j,1,n,l)
397        qvhn(i,j,1,n,l)=qvhn(i,j,2,n,l)
398        tthn(i,j,1,n,l)=tt2n(i,j,1,n,l)
399      end do
400    end do
401
402    if(iumax.ne.nuvz-1) stop &
403         'READWIND: NUVZ NOT CONSISTENT FOR A NESTING LEVEL'
404    if(iwmax.ne.nwz) stop &
405         'READWIND: NWZ NOT CONSISTENT FOR A NESTING LEVEL'
406
407  end do
408
409  return
410888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
411  write(*,*) ' #### ',wfnamen(l,indj),' FOR NESTING LEVEL  #### '
412  write(*,*) ' #### ',l,' IS NOT GRIB FORMAT !!!           #### '
413  stop 'Execution terminated'
414
415
416999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
417  write(*,*) ' #### ',wfnamen(l,indj),'                    #### '
418  write(*,*) ' #### CANNOT BE OPENED FOR NESTING LEVEL ',l,'####'
419
420end subroutine readwind_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG