source: branches/petra/src/readwind.f90 @ 37

Last change on this file since 37 was 37, checked in by pesei, 9 years ago

Wet dep quick fix and other small changes. Wet depo quick fix not final yet.

File size: 20.4 KB
Line 
1!**********************************************************************
2! Copyright 1998-2015                                                 *
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  !             CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with *
34  !                                 ECMWF grib_api                      *
35  !             CHANGE: 03/12/2008, Harald Sodemann, update to f90 with *
36  !                                 ECMWF grib_api                      *
37  !                                                                     *
38  !**********************************************************************
39  !  Changes, Bernd C. Krueger, Feb. 2001:
40  !   Variables tth and qvh (on eta coordinates) in common block
41  ! 2/2015, PS: add missing paramter iret in call to grib subr
42  ! 3/2015, PS: add some verbosity messages
43  ! 3/2015, PS: add GRIB2 coding for deta/dt dp/deta used by ZAMG
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  parameters for grib_api
74  integer :: ifile
75  integer :: iret
76  integer :: igrib
77  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl,parId
78  integer :: gotGrid
79  !HSO  end
80
81  real(kind=4) :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
82  real(kind=4) :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
83  real(kind=4) :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
84  integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax
85
86  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
87
88  ! dimension of isec2 at least (22+n), where n is the number of parallels or
89  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
90
91  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
92  ! coordinate parameters
93
94  integer :: isec1(56),isec2(22+nxmax+nymax)
95  real(kind=4) :: zsec4(jpunp)
96  real(kind=4) :: xaux,yaux
97  real(kind=8) :: xauxin,yauxin
98  real,parameter :: eps=1.e-4
99  real(kind=4) :: nsss(0:nxmax-1,0:nymax-1),ewss(0:nxmax-1,0:nymax-1)
100  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1,conversion_factor
101
102  logical :: hflswitch,strswitch
103
104  !HSO  grib api error messages
105  character(len=24) :: gribErrorMsg = 'Error reading grib file'
106  character(len=20) :: gribFunction = 'readwind'
107
108  hflswitch=.false.
109  strswitch=.false.
110  levdiff2=nlev_ec-nwz+1
111  iumax=0
112  iwmax=0
113
114  !
115  ! OPENING OF DATA FILE (GRIB CODE)
116  !
117 
1185 continue
119  call verboutput(verbosity,1,'      call grib_open_file '//path(3)(1:length(3))&
120         //trim(wfname(indj)))
121  call grib_open_file(ifile,path(3)(1:length(3)) &
122         //trim(wfname(indj)),'r',iret)
123  if (iret.ne.GRIB_SUCCESS) then
124    goto 888   ! ERROR DETECTED
125  endif
126  !turn on support for multi fields messages */
127  !call grib_multi_support_on
128
129  gotGrid=0
130  ifield=0
13110   ifield=ifield+1
132  !
133  ! GET NEXT FIELDS
134  !
135  call verboutput(verbosity,2,'      call grib_new_from_file')
136  call grib_new_from_file(ifile,igrib,iret)
137  if (iret.eq.GRIB_END_OF_FILE)  then
138    goto 50    ! EOF DETECTED
139  elseif (iret.ne.GRIB_SUCCESS) then
140    goto 888   ! ERROR DETECTED
141  endif
142
143  !first see if we read GRIB1 or GRIB2
144  call verboutput(verbosity,2,'      call grib_get_int')
145  call grib_get_int(igrib,'editionNumber',gribVer,iret)
146  call verboutput(verbosity,2,'      call grib_check')
147  call grib_check(iret,gribFunction,gribErrorMsg)
148
149  if (gribVer.eq.1) then ! GRIB Edition 1
150
151    !print*,'GRiB Edition 1'
152
153    ! read the grib2 identifiers
154
155    call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
156    call grib_check(iret,gribFunction,gribErrorMsg)
157    call grib_get_int(igrib,'level',isec1(8),iret)
158    call grib_check(iret,gribFunction,gribErrorMsg)
159
160    !change code for etadot to code for omega
161    if (isec1(6).eq.77) then
162      isec1(6)=135
163    endif
164
165    conversion_factor=1.
166
167  else ! GRIB 2
168
169    !print*,'GRiB Edition 2'
170
171    ! read the grib2 identifiers
172
173    call grib_get_int(igrib,'discipline',discipl,iret)
174    call grib_check(iret,gribFunction,gribErrorMsg)
175    call grib_get_int(igrib,'parameterCategory',parCat,iret)
176    call grib_check(iret,gribFunction,gribErrorMsg)
177    call grib_get_int(igrib,'parameterNumber',parNum,iret)
178    call grib_check(iret,gribFunction,gribErrorMsg)
179    call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
180    call grib_check(iret,gribFunction,gribErrorMsg)
181    call grib_get_int(igrib,'level',valSurf,iret)
182    call grib_check(iret,gribFunction,gribErrorMsg)
183    call grib_get_int(igrib,'paramId',parId,iret)
184    call grib_check(iret,gribFunction,gribErrorMsg)
185
186    !print*,discipl,parCat,parNum,typSurf,valSurf
187
188    ! convert to grib1 identifiers
189
190    isec1(6)=-1
191    isec1(7)=-1
192    isec1(8)=-1
193    isec1(8)=valSurf     ! level
194    conversion_factor=1.
195    if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
196      isec1(6)=130         ! indicatorOfParameter
197    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
198      isec1(6)=131         ! indicatorOfParameter
199    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
200      isec1(6)=132         ! indicatorOfParameter
201    elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
202      isec1(6)=133         ! indicatorOfParameter
203    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
204      isec1(6)=134         ! indicatorOfParameter
205    elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually detadtdpdeta
206      isec1(6)=135         ! indicatorOfParameter
207    elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually detadtdpdeta
208      isec1(6)=135         ! indicatorOfParameter
209    elseif ((parCat.eq.2).and.(parNum.eq.8)) then ! Omega, actually detadtdpdeta
210      isec1(6)=135         ! indicatorOfParameter
211    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
212      isec1(6)=151         ! indicatorOfParameter
213    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
214      isec1(6)=165         ! indicatorOfParameter
215    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
216      isec1(6)=166         ! indicatorOfParameter
217    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
218      isec1(6)=167         ! indicatorOfParameter
219    elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
220      isec1(6)=168         ! indicatorOfParameter
221    elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
222      isec1(6)=141         ! indicatorOfParameter
223      conversion_factor=1000.
224    elseif ((parCat.eq.6).and.(parNum.eq.1) .or. parId .eq. 164) then ! CC
225      isec1(6)=164         ! indicatorOfParameter
226    elseif ((parCat.eq.1).and.(parNum.eq.9) .or. parId .eq. 142) then ! LSP
227      isec1(6)=142         ! indicatorOfParameter
228    elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
229      isec1(6)=143         ! indicatorOfParameter
230      conversion_factor=1000.
231    elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
232      isec1(6)=146         ! indicatorOfParameter
233    elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
234      isec1(6)=176         ! indicatorOfParameter
235    elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS
236      isec1(6)=180         ! indicatorOfParameter
237    elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS
238      isec1(6)=181         ! indicatorOfParameter
239    elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
240      isec1(6)=129         ! indicatorOfParameter
241    elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO
242      isec1(6)=160         ! indicatorOfParameter
243    elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
244         (typSurf.eq.1)) then ! LSM
245      isec1(6)=172         ! indicatorOfParameter
246    else
247      print*,'***WARNING: undefined GRiB2 message found!',discipl, &
248           parCat,parNum,typSurf
249    endif
250    if(parId .ne. isec1(6) .and. parId .ne. 77) then
251      write(*,*) 'parId',parId, 'isec1(6)',isec1(6)
252  !    stop
253    endif
254
255  endif ! end GRIB 1 / GRIB 2 cases
256
257  !HSO  get the size and data of the values array
258  if (isec1(6).ne.-1) then
259    call grib_get_real4_array(igrib,'values',zsec4,iret)
260    call grib_check(iret,gribFunction,gribErrorMsg)
261  endif
262
263  !HSO  get the required fields from section 2 in a gribex compatible manner
264  if (ifield.eq.1) then
265  call grib_get_int(igrib,'numberOfPointsAlongAParallel',isec2(2),iret)
266  call grib_check(iret,gribFunction,gribErrorMsg)
267  call grib_get_int(igrib,'numberOfPointsAlongAMeridian',isec2(3),iret)
268  call grib_check(iret,gribFunction,gribErrorMsg)
269  call grib_get_int(igrib,'numberOfVerticalCoordinateValues',isec2(12),iret)
270  call grib_check(iret,gribFunction,gribErrorMsg)
271  ! CHECK GRID SPECIFICATIONS
272  if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
273  if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
274  if(isec2(12)/2-1.ne.nlev_ec) &
275  stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
276  endif ! ifield
277
278  !HSO  get the second part of the grid dimensions only from GRiB1 messages
279  if (isec1(6) .eq. 167 .and. (gotGrid.eq.0)) then
280    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
281         xauxin,iret)
282    call grib_check(iret,gribFunction,gribErrorMsg)
283    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
284         yauxin,iret)
285    call grib_check(iret,gribFunction,gribErrorMsg)
286    if (xauxin.gt.180.) xauxin=xauxin-360.0
287    if (xauxin.lt.-180.) xauxin=xauxin+360.0
288
289    xaux=xauxin+real(nxshift)*dx
290    yaux=yauxin
291    if (xaux.gt.180.) xaux=xaux-360.0
292    if(abs(xaux-xlon0).gt.eps) &
293    stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
294    if(abs(yaux-ylat0).gt.eps) &
295    stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
296    gotGrid=1
297  endif ! gotGrid
298
299  do j=0,nymin1
300    do i=0,nxfield-1
301      k=isec1(8)
302      if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
303           zsec4(nxfield*(ny-j-1)+i+1)
304      if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
305           zsec4(nxfield*(ny-j-1)+i+1)
306      if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
307           zsec4(nxfield*(ny-j-1)+i+1)
308      if(isec1(6).eq.133) then                      !! SPEC. HUMIDITY
309        qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
310        if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
311             qvh(i,j,nlev_ec-k+2,n) = 0.
312  !        this is necessary because the gridded data may contain
313  !        spurious negative values
314      endif
315      if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
316           zsec4(nxfield*(ny-j-1)+i+1)
317
318      if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
319           zsec4(nxfield*(ny-j-1)+i+1)
320      if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
321           zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
322      if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
323           zsec4(nxfield*(ny-j-1)+i+1)
324      if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
325           zsec4(nxfield*(ny-j-1)+i+1)
326      if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
327           zsec4(nxfield*(ny-j-1)+i+1)
328      if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
329           zsec4(nxfield*(ny-j-1)+i+1)
330      if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
331           zsec4(nxfield*(ny-j-1)+i+1)
332      if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
333           zsec4(nxfield*(ny-j-1)+i+1)
334      if(isec1(6).eq.142) then                      !! LARGE SCALE PREC.
335        lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
336        if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
337      endif
338      if(isec1(6).eq.143) then                      !! CONVECTIVE PREC.
339        convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
340        if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
341      endif
342      if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
343           zsec4(nxfield*(ny-j-1)+i+1)
344      if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
345           hflswitch=.true.    ! Heat flux available
346      if(isec1(6).eq.176) then                      !! SOLAR RADIATION
347        ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
348        if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
349      endif
350      if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
351           zsec4(nxfield*(ny-j-1)+i+1)
352      if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
353           zsec4(nxfield*(ny-j-1)+i+1)
354      if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
355           (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
356  !sec        strswitch=.true.
357      if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
358           zsec4(nxfield*(ny-j-1)+i+1)/ga
359      if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
360           zsec4(nxfield*(ny-j-1)+i+1)
361      if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
362           zsec4(nxfield*(ny-j-1)+i+1)
363      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
364      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
365
366    end do
367  end do
368
369  call verboutput(verbosity,2,'      call grib_release')
370  call grib_release(igrib)
371  goto 10                      !! READ NEXT LEVEL OR PARAMETER
372  !
373  ! CLOSING OF INPUT DATA FILE
374  !
375
37650   call grib_close_file(ifile)
377
378  !error message if no fields found with correct first longitude in it
379  if (gotGrid.eq.0) then
380    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
381         'messages'
382    stop
383  endif
384
385  if(levdiff2.eq.0) then
386    iwmax=nlev_ec+1
387    do i=0,nxmin1
388      do j=0,nymin1
389        wwh(i,j,nlev_ec+1)=0.
390      end do
391    end do
392  endif
393
394  ! For global fields, assign the leftmost data column also to the rightmost
395  ! data column; if required, shift whole grid by nxshift grid points
396  !*************************************************************************
397
398  if (xglobal) then
399    call verboutput(verbosity,2,'      call shift_field_0')
400    call shift_field_0(ewss,nxfield,ny)
401    call shift_field_0(nsss,nxfield,ny)
402    call shift_field_0(oro,nxfield,ny)
403    call shift_field_0(excessoro,nxfield,ny)
404    call shift_field_0(lsm,nxfield,ny)
405    call shift_field(ps,nxfield,ny,1,1,2,n)
406    call shift_field(sd,nxfield,ny,1,1,2,n)
407    call shift_field(msl,nxfield,ny,1,1,2,n)
408    call shift_field(tcc,nxfield,ny,1,1,2,n)
409    call shift_field(u10,nxfield,ny,1,1,2,n)
410    call shift_field(v10,nxfield,ny,1,1,2,n)
411    call shift_field(tt2,nxfield,ny,1,1,2,n)
412    call shift_field(td2,nxfield,ny,1,1,2,n)
413    call shift_field(lsprec,nxfield,ny,1,1,2,n)
414    call shift_field(convprec,nxfield,ny,1,1,2,n)
415    call shift_field(sshf,nxfield,ny,1,1,2,n)
416    call shift_field(ssr,nxfield,ny,1,1,2,n)
417    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
418    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
419    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
420    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
421    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
422  endif
423
424  do i=0,nxmin1
425    do j=0,nymin1
426      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
427    end do
428  end do
429
430  if ((.not.hflswitch).or.(.not.strswitch)) then
431    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
432         wfname(indj)
433
434  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
435  ! As ECMWF has increased the model resolution, such that now the first model
436  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
437  ! (3rd model level in FLEXPART) for the profile method
438  !***************************************************************************
439
440    do i=0,nxmin1
441      do j=0,nymin1
442        plev1=akz(3)+bkz(3)*ps(i,j,1,n)
443        pmean=0.5*(ps(i,j,1,n)+plev1)
444        tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
445        fu=-r_air*tv/ga/pmean
446        hlev1=fu*(plev1-ps(i,j,1,n))   ! HEIGTH OF FIRST MODEL LAYER
447        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
448        fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
449        call verboutput(verbosity,2,'      call pbl_profile')
450        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
451             tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
452             surfstr(i,j,1,n),sshf(i,j,1,n))
453        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
454        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
455      end do
456    end do
457  endif
458
459
460  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
461  ! level at the ground
462  ! Specific humidity is taken the same as at one level above
463  ! Temperature is taken as 2 m temperature
464  !**************************************************************************
465
466     do i=0,nxmin1
467        do j=0,nymin1
468           uuh(i,j,1)=u10(i,j,1,n)
469           vvh(i,j,1)=v10(i,j,1,n)
470           qvh(i,j,1,n)=qvh(i,j,2,n)
471           tth(i,j,1,n)=tt2(i,j,1,n)
472        end do
473     end do
474
475  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
476  if(iwmax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
477
478  call verboutput(verbosity,1,'      leaving readwind')
479  return
480
481888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
482  write(*,*) ' #### ',wfname(indj),'                    #### '
483  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
484  stop 'Execution terminated'
485999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
486  write(*,*) ' #### ',wfname(indj),'                    #### '
487  write(*,*) ' #### CANNOT BE OPENED !!!                    #### '
488  stop 'Execution terminated'
489
490end subroutine readwind
491
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG