source: flexpart.git/src/readwind_ecmwf.f90 @ 553d0a7

scaling-bug
Last change on this file since 553d0a7 was 553d0a7, checked in by Pirmin Kaufmann <pirmin.kaufmann@…>, 18 months ago

Removed conversion factor for CP in GRIB2. Changed conversion factor for SDE (snow depth) and added recognition of SD (water equivalent) parameterNumber in GRIB2.

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