source: flexpart.git/src/readwind.f90 @ 38b7917

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

Initial code to handle cloud water in nested wind fields (it is not completely implemented in this commit)

  • Property mode set to 100644
File size: 20.7 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!ZHG test the grib fields that have lcwc without using them
109!  readcloud=.false.
110
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! ESO Cloud water is in a) fields CLWC and CIWC, *or* b) field QC
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! ESO qc(=clwc+ciwc):
198    elseif ((parCat.eq.201).and.(parNum.eq.31).and.(typSurf.eq.105)) then ! qc
199      isec1(6)=201031         ! indicatorOfParameter
200    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
201      isec1(6)=134         ! indicatorOfParameter
202    elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
203      isec1(6)=135         ! indicatorOfParameter
204    elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually eta dot
205      isec1(6)=135         ! indicatorOfParameter
206    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
207      isec1(6)=151         ! indicatorOfParameter
208    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
209      isec1(6)=165         ! indicatorOfParameter
210    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
211      isec1(6)=166         ! indicatorOfParameter
212    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
213      isec1(6)=167         ! indicatorOfParameter
214    elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
215      isec1(6)=168         ! indicatorOfParameter
216    elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
217      isec1(6)=141         ! indicatorOfParameter
218      conversion_factor=1000.
219    elseif ((parCat.eq.6).and.(parNum.eq.1) .or. parId .eq. 164) then ! CC
220      isec1(6)=164         ! indicatorOfParameter
221    elseif ((parCat.eq.1).and.(parNum.eq.9) .or. parId .eq. 142) then ! LSP
222      isec1(6)=142         ! indicatorOfParameter
223    elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
224      isec1(6)=143         ! indicatorOfParameter
225      conversion_factor=1000.
226    elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
227      isec1(6)=146         ! indicatorOfParameter
228    elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
229      isec1(6)=176         ! indicatorOfParameter
230!    elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS --wrong
231    elseif ((parCat.eq.2).and.(parNum.eq.38) .or. parId .eq. 180) then ! EWSS --correct
232      isec1(6)=180         ! indicatorOfParameter
233!    elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS --wrong
234    elseif ((parCat.eq.2).and.(parNum.eq.37) .or. parId .eq. 181) then ! NSSS --correct
235      isec1(6)=181         ! indicatorOfParameter
236    elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
237      isec1(6)=129         ! indicatorOfParameter
238    elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO
239      isec1(6)=160         ! indicatorOfParameter
240    elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
241         (typSurf.eq.1)) then ! LSM
242      isec1(6)=172         ! indicatorOfParameter
243    else
244      print*,'***WARNING: undefined GRiB2 message found!',discipl, &
245           parCat,parNum,typSurf
246    endif
247    if(parId .ne. isec1(6) .and. parId .ne. 77) then
248      write(*,*) 'parId',parId, 'isec1(6)',isec1(6)
249!    stop
250    endif
251
252  endif
253
254!HSO  get the size and data of the values array
255  if (isec1(6).ne.-1) then
256    call grib_get_real4_array(igrib,'values',zsec4,iret)
257    call grib_check(iret,gribFunction,gribErrorMsg)
258  endif
259
260!HSO  get the required fields from section 2 in a gribex compatible manner
261  if (ifield.eq.1) then
262    call grib_get_int(igrib,'numberOfPointsAlongAParallel',isec2(2),iret)
263    call grib_check(iret,gribFunction,gribErrorMsg)
264    call grib_get_int(igrib,'numberOfPointsAlongAMeridian',isec2(3),iret)
265    call grib_check(iret,gribFunction,gribErrorMsg)
266    call grib_get_int(igrib,'numberOfVerticalCoordinateValues',isec2(12))
267    call grib_check(iret,gribFunction,gribErrorMsg)
268! CHECK GRID SPECIFICATIONS
269    if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
270    if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
271    if(isec2(12)/2-1.ne.nlev_ec) &
272         stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
273  endif ! ifield
274
275!HSO  get the second part of the grid dimensions only from GRiB1 messages
276  if (isec1(6) .eq. 167 .and. (gotGrid.eq.0)) then
277    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
278         xauxin,iret)
279    call grib_check(iret,gribFunction,gribErrorMsg)
280    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
281         yauxin,iret)
282    call grib_check(iret,gribFunction,gribErrorMsg)
283    if (xauxin.gt.180.) xauxin=xauxin-360.0
284    if (xauxin.lt.-180.) xauxin=xauxin+360.0
285
286    xaux=xauxin+real(nxshift)*dx
287    yaux=yauxin
288    if (xaux.gt.180.) xaux=xaux-360.0
289    if(abs(xaux-xlon0).gt.eps) &
290         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
291    if(abs(yaux-ylat0).gt.eps) &
292         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
293    gotGrid=1
294  endif ! gotGrid
295
296  do j=0,nymin1
297    do i=0,nxfield-1
298      k=isec1(8)
299      if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
300           zsec4(nxfield*(ny-j-1)+i+1)
301      if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
302           zsec4(nxfield*(ny-j-1)+i+1)
303      if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
304           zsec4(nxfield*(ny-j-1)+i+1)
305      if(isec1(6).eq.133) then                      !! SPEC. HUMIDITY
306        qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
307        if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
308             qvh(i,j,nlev_ec-k+2,n) = 0.
309!        this is necessary because the gridded data may contain
310!        spurious negative values
311      endif
312      if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
313           zsec4(nxfield*(ny-j-1)+i+1)
314
315      if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
316           zsec4(nxfield*(ny-j-1)+i+1)
317      if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
318           zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
319      if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
320           zsec4(nxfield*(ny-j-1)+i+1)
321      if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
322           zsec4(nxfield*(ny-j-1)+i+1)
323      if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
324           zsec4(nxfield*(ny-j-1)+i+1)
325      if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
326           zsec4(nxfield*(ny-j-1)+i+1)
327      if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
328           zsec4(nxfield*(ny-j-1)+i+1)
329      if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
330           zsec4(nxfield*(ny-j-1)+i+1)
331      if(isec1(6).eq.142) then                      !! LARGE SCALE PREC.
332        lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
333        if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
334      endif
335      if(isec1(6).eq.143) then                      !! CONVECTIVE PREC.
336        convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
337        if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
338      endif
339      if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
340           zsec4(nxfield*(ny-j-1)+i+1)
341      if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
342           hflswitch=.true.    ! Heat flux available
343      if(isec1(6).eq.176) then                      !! SOLAR RADIATION
344        ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
345        if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
346      endif
347      if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
348           zsec4(nxfield*(ny-j-1)+i+1)
349      if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
350           zsec4(nxfield*(ny-j-1)+i+1)
351      if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
352           (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
353!sec        strswitch=.true.
354      if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
355           zsec4(nxfield*(ny-j-1)+i+1)/ga
356      if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
357           zsec4(nxfield*(ny-j-1)+i+1)
358      if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
359           zsec4(nxfield*(ny-j-1)+i+1)
360      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
361      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
362!ZHG READING CLOUD FIELDS ASWELL
363! ESO TODO: add check for if one of clwc/ciwc missing (error),
364! also if all 3 cw fields present, use qc and disregard the others
365      if(isec1(6).eq.246) then  !! CLWC  Cloud liquid water content [kg/kg]
366        clwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
367        readclouds=.true.
368        sumclouds=.false.
369      endif
370      if(isec1(6).eq.247) then  !! CIWC  Cloud ice water content
371        ciwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
372      endif
373!ZHG end
374!ESO read qc (=clwc+ciwc)
375      if(isec1(6).eq.201031) then  !! QC  Cloud liquid water content [kg/kg]
376        clwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
377        readclouds=.true.
378        sumclouds=.true.
379      endif
380
381    end do
382  end do
383
384  call grib_release(igrib)
385  goto 10                      !! READ NEXT LEVEL OR PARAMETER
386!
387! CLOSING OF INPUT DATA FILE
388!
389
39050 call grib_close_file(ifile)
391
392!error message if no fields found with correct first longitude in it
393  if (gotGrid.eq.0) then
394    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
395         'messages'
396    stop
397  endif
398
399  if(levdiff2.eq.0) then
400    iwmax=nlev_ec+1
401    do i=0,nxmin1
402      do j=0,nymin1
403        wwh(i,j,nlev_ec+1)=0.
404      end do
405    end do
406  endif
407
408! For global fields, assign the leftmost data column also to the rightmost
409! data column; if required, shift whole grid by nxshift grid points
410!*************************************************************************
411
412  if (xglobal) then
413    call shift_field_0(ewss,nxfield,ny)
414    call shift_field_0(nsss,nxfield,ny)
415    call shift_field_0(oro,nxfield,ny)
416    call shift_field_0(excessoro,nxfield,ny)
417    call shift_field_0(lsm,nxfield,ny)
418    call shift_field(ps,nxfield,ny,1,1,2,n)
419    call shift_field(sd,nxfield,ny,1,1,2,n)
420    call shift_field(msl,nxfield,ny,1,1,2,n)
421    call shift_field(tcc,nxfield,ny,1,1,2,n)
422    call shift_field(u10,nxfield,ny,1,1,2,n)
423    call shift_field(v10,nxfield,ny,1,1,2,n)
424    call shift_field(tt2,nxfield,ny,1,1,2,n)
425    call shift_field(td2,nxfield,ny,1,1,2,n)
426    call shift_field(lsprec,nxfield,ny,1,1,2,n)
427    call shift_field(convprec,nxfield,ny,1,1,2,n)
428    call shift_field(sshf,nxfield,ny,1,1,2,n)
429    call shift_field(ssr,nxfield,ny,1,1,2,n)
430    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
431    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
432    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
433    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
434    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
435!ZHG
436    call shift_field(clwch,nxfield,ny,nuvzmax,nuvz,2,n)
437    if (.not.sumclouds) call shift_field(ciwch,nxfield,ny,nuvzmax,nuvz,2,n)
438!ZHG end
439
440  endif
441
442  do i=0,nxmin1
443    do j=0,nymin1
444      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
445    end do
446  end do
447
448  if ((.not.hflswitch).or.(.not.strswitch)) then
449    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
450         wfname(indj)
451
452! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
453! As ECMWF has increased the model resolution, such that now the first model
454! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
455! (3rd model level in FLEXPART) for the profile method
456!***************************************************************************
457
458    do i=0,nxmin1
459      do j=0,nymin1
460        plev1=akz(3)+bkz(3)*ps(i,j,1,n)
461        pmean=0.5*(ps(i,j,1,n)+plev1)
462        tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
463        fu=-r_air*tv/ga/pmean
464        hlev1=fu*(plev1-ps(i,j,1,n))   ! HEIGTH OF FIRST MODEL LAYER
465        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
466        fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
467        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
468             tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
469             surfstr(i,j,1,n),sshf(i,j,1,n))
470        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
471        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
472      end do
473    end do
474  endif
475
476
477! Assign 10 m wind to model level at eta=1.0 to have one additional model
478! level at the ground
479! Specific humidity is taken the same as at one level above
480! Temperature is taken as 2 m temperature
481!**************************************************************************
482
483  do i=0,nxmin1
484    do j=0,nymin1
485      uuh(i,j,1)=u10(i,j,1,n)
486      vvh(i,j,1)=v10(i,j,1,n)
487      qvh(i,j,1,n)=qvh(i,j,2,n)
488      tth(i,j,1,n)=tt2(i,j,1,n)
489    end do
490  end do
491
492  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
493  if(iwmax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
494
495  return
496
497888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
498  write(*,*) ' #### ',wfname(indj),'                    #### '
499  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
500  stop 'Execution terminated'
501
502end subroutine readwind
503
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG