source: flexpart.git/src/readwind_ecmwf.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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