source: flexpart.git/src/readwind.f90 @ 5f9d14a

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

Updated wet depo scheme

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