CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
readwind.F90
Go to the documentation of this file.
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 
22 subroutine readwind_ecmwf(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  !**********************************************************************
49  ! Changes Arnold, D. and Morton, D. (2015): *
50  ! -- description of local and common variables *
51  !**********************************************************************
52 
53  use grib_api
54  use par_mod
55  use com_mod
56 
57  implicit none
58 
59  !***********************************************************************
60  ! Subroutine Parameters: *
61  ! input *
62  ! indj indicates number of the wind field to be read in *
63  ! n temporal index for meteorological fields (1 to 3)*
64  ! uuh,vvh, wwh wind components in ECMWF model levels *
65  integer :: indj,n
66  real(kind=4) :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
67  real(kind=4) :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
68  real(kind=4) :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
69  !***********************************************************************
70 
71  !***********************************************************************
72  ! Local variables *
73  ! *
74  ! HSO variables for grib_api: *
75  ! ifile grib file to be opened and read in *
76  ! iret is a return code for successful or not open *
77  ! igrib grib edition number (whether GRIB1 or GRIB2) *
78  ! gribVer where info on igrib is kept, 1 for GRIB1 and *
79  ! 2 for GRIB2 *
80  ! parCat parameter category ( = number) , how FLEXPART *
81  ! identifies fields *
82  ! parNum parameter number by product discipline and *
83  ! parameter category *
84  ! typSurf type of first fixed surface *
85  ! valSurf level *
86  ! discipl discipline of processed data contained in a *
87  ! GRIB message *
88  integer :: ifile
89  integer :: iret
90  integer :: igrib
91  integer :: gribver,parcat,parnum,typsurf,valsurf,discipl
92 #if defined WITH_CTBTO_PATCHES
93  integer :: parid
94 #endif
95  integer :: gotgrid
96  !
97  ! Variables and arrays for grib decoding: *
98  ! dimension of isec2 at least (22+n), where n is the number of *
99  ! parallels or meridians in a quasi-regular (reduced) Gaussian or *
100  ! lat/long grid *
101  ! dimension of zsec2 at least (10+nn), where nn is the number of *
102  ! vertical coordinate parameters *
103  ! *
104  ! isec1 grib definition section (version, center, ... *
105  ! isec2 grid description section *
106  ! zsec4 the binary data section *
107  ! xaux
108  ! yaux
109  ! xauxin
110  ! yauxin
111  ! xaux0 auxiliary variable for xlon0 *
112  ! yaux0 auxiliary variable for xlat0 *
113  ! nsss NS surface stress *
114  ! ewss EW surface stress *
115  ! plev1 pressure of the first model layer *
116  ! pmean mean pressure between ?? *
117  ! tv virtual temperature *
118  ! fu for conversion from pressure to meters *
119  ! hlev1 height of the first model layer *
120  ! ff10m wind speed at 10 m *
121  ! fflev1 wind speed at the first model layere *
122  ! hflswitch logical variable to check existence of flux data *
123  ! strswitch logical variable to check existence of stress *
124  ! data *
125 
126  integer :: isec1(56),isec2(22+nxmax+nymax)
127  real(kind=4) :: zsec4(jpunp)
128  real(kind=4) :: xaux,yaux,xaux0,yaux0
129 #if defined WITH_CTBTO_PATCHES
130  real(kind=4) :: conversion_factor
131 #endif
132  real(kind=8) :: xauxin,yauxin
133  real(kind=4) :: nsss(0:nxmax-1,0:nymax-1),ewss(0:nxmax-1,0:nymax-1)
134  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
135  logical :: hflswitch,strswitch
136 
137  ! Other variables:
138  ! i,j,k loop control indices in each direction *
139  ! levdiff2 number of model levels - number of model levels *
140  ! for the staggered vertical wind +1 *
141  ! ifield index to control field read in *
142  ! iumax
143  ! iwmax
144  integer :: i,j,k,levdiff2,ifield,iumax,iwmax
145  !***********************************************************************
146 
147  !***********************************************************************
148  ! Local constants *
149  real,parameter :: eps=1.e-4
150  !HSO grib api error messages:
151  character(len=24) :: griberrormsg = 'Error reading grib file'
152  character(len=20) :: gribfunction = 'readwind'
153  !***********************************************************************
154 
155  !***********************************************************************
156  ! Global variables *
157  ! from par_mod and com_mod *
158  ! wfname File name of data to be read in *
159  ! nx,ny,nuvz,nwz expected field dimensions *
160  ! nxfield same as nx but for limited area fields *
161  ! nxmin1,nymin1 nx-1, ny-1, respectively *
162  ! nxmax,nymax maximum dimension of wind fields in x and y *
163  ! direction, respectively *
164  ! nuvzmax,nwzmax maximum dimension of (u,v) and (w) wind fields in z
165  ! direction (for fields on eta levels) *
166  ! nlev_ec number of vertical levels ecmwf model *
167  ! u10,v10 wind fields at 10 m *
168  ! tt2,td2 temperature and dew point temperature at 2 m *
169  ! tth,qvh temperature and specific humidity on original *
170  ! eta levels *
171  ! ps surface pressure *
172  ! sd snow depth (but not used afterwards!) *
173  ! msl mean sea level pressure *
174  ! ttc total cloud cover *
175  ! lsprec large scale precipitation *
176  ! convprec convective precipitation *
177  ! sshf sensible heat flux *
178  ! ssr solar radiation *
179  ! surfstr surface stress *
180  ! oro orography *
181  ! excessoro standard deviation of orography *
182  ! lsm land sea mask *
183  ! r_air individual gas constant for dry air [J/kg/K] *
184  ! ga gravity acceleration of earth [m/s**2] *
185  !***********************************************************************
186 
187 !-----------------------------------------------------------------------------
188 
189 
190  hflswitch=.false.
191  strswitch=.false.
192  levdiff2=nlev_ec-nwz+1
193  iumax=0
194  iwmax=0
195 
196  !
197  ! OPENING OF DATA FILE (GRIB CODE)
198  !
199 5 call grib_open_file(ifile,path(3)(1:length(3)) &
200  //trim(wfname(indj)),'r',iret)
201  if (iret.ne.grib_success) then
202  goto 888 ! ERROR DETECTED
203  endif
204  !turn on support for multi fields messages */
205  !call grib_multi_support_on
206 
207  gotgrid=0
208  ifield=0
209 10 ifield=ifield+1
210  !
211  ! GET NEXT FIELDS
212  !
213  call grib_new_from_file(ifile,igrib,iret)
214  if (iret.eq.grib_end_of_file) then
215  goto 50 ! EOF DETECTED
216  elseif (iret.ne.grib_success) then
217  goto 888 ! ERROR DETECTED
218  endif
219 
220 #if defined WITH_CTBTO_PATCHES
221  conversion_factor=1.0
222 #endif
223 
224  !first see if we read GRIB1 or GRIB2
225  call grib_get_int(igrib,'editionNumber',gribver,iret)
226  call grib_check(iret,gribfunction,griberrormsg)
227 
228  if (gribver.eq.1) then ! GRIB Edition 1
229 
230  !print*,'GRiB Edition 1'
231  !read the grib2 identifiers
232  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
233  call grib_check(iret,gribfunction,griberrormsg)
234  call grib_get_int(igrib,'level',isec1(8),iret)
235  call grib_check(iret,gribfunction,griberrormsg)
236 
237  !change code for etadot to code for omega
238 #if defined WITH_CTBTO_PATCHES
239  if (isec1(6).eq.77.or.isec1(6).eq.87) then
240 #else
241  if (isec1(6).eq.77) then
242 #endif
243  isec1(6)=135
244  endif
245 
246  else
247 
248 #if defined WITH_CTBTO_PATCHES
249  call grib2check(igrib,isec1,conversion_factor)
250 #else
251  !print*,'GRiB Edition 2'
252  !read the grib2 identifiers
253  call grib_get_int(igrib,'discipline',discipl,iret)
254  call grib_check(iret,gribfunction,griberrormsg)
255  call grib_get_int(igrib,'parameterCategory',parcat,iret)
256  call grib_check(iret,gribfunction,griberrormsg)
257  call grib_get_int(igrib,'parameterNumber',parnum,iret)
258  call grib_check(iret,gribfunction,griberrormsg)
259  call grib_get_int(igrib,'typeOfFirstFixedSurface',typsurf,iret)
260  call grib_check(iret,gribfunction,griberrormsg)
261  call grib_get_int(igrib,'level',valsurf,iret)
262  call grib_check(iret,gribfunction,griberrormsg)
263 
264  !print*,discipl,parCat,parNum,typSurf,valSurf
265 
266  !convert to grib1 identifiers
267  isec1(6)=-1
268  isec1(7)=-1
269  isec1(8)=-1
270  isec1(8)=valsurf ! level
271  if ((parcat.eq.0).and.(parnum.eq.0).and.(typsurf.eq.105)) then ! T
272  isec1(6)=130 ! indicatorOfParameter
273  elseif ((parcat.eq.2).and.(parnum.eq.2).and.(typsurf.eq.105)) then ! U
274  isec1(6)=131 ! indicatorOfParameter
275  elseif ((parcat.eq.2).and.(parnum.eq.3).and.(typsurf.eq.105)) then ! V
276  isec1(6)=132 ! indicatorOfParameter
277  elseif ((parcat.eq.1).and.(parnum.eq.0).and.(typsurf.eq.105)) then ! Q
278  isec1(6)=133 ! indicatorOfParameter
279  elseif ((parcat.eq.3).and.(parnum.eq.0).and.(typsurf.eq.1)) then !SP
280  isec1(6)=134 ! indicatorOfParameter
281  elseif ((parcat.eq.2).and.(parnum.eq.32)) then ! W, actually eta dot
282  isec1(6)=135 ! indicatorOfParameter
283  elseif ((parcat.eq.3).and.(parnum.eq.0).and.(typsurf.eq.101)) then !SLP
284  isec1(6)=151 ! indicatorOfParameter
285  elseif ((parcat.eq.2).and.(parnum.eq.2).and.(typsurf.eq.103)) then ! 10U
286  isec1(6)=165 ! indicatorOfParameter
287  elseif ((parcat.eq.2).and.(parnum.eq.3).and.(typsurf.eq.103)) then ! 10V
288  isec1(6)=166 ! indicatorOfParameter
289  elseif ((parcat.eq.0).and.(parnum.eq.0).and.(typsurf.eq.103)) then ! 2T
290  isec1(6)=167 ! indicatorOfParameter
291  elseif ((parcat.eq.0).and.(parnum.eq.6).and.(typsurf.eq.103)) then ! 2D
292  isec1(6)=168 ! indicatorOfParameter
293  elseif ((parcat.eq.1).and.(parnum.eq.11).and.(typsurf.eq.1)) then ! SD
294  isec1(6)=141 ! indicatorOfParameter
295  elseif ((parcat.eq.6).and.(parnum.eq.1)) then ! CC
296  isec1(6)=164 ! indicatorOfParameter
297  elseif ((parcat.eq.1).and.(parnum.eq.9)) then ! LSP
298  isec1(6)=142 ! indicatorOfParameter
299  elseif ((parcat.eq.1).and.(parnum.eq.10)) then ! CP
300  isec1(6)=143 ! indicatorOfParameter
301  elseif ((parcat.eq.0).and.(parnum.eq.11).and.(typsurf.eq.1)) then ! SHF
302  isec1(6)=146 ! indicatorOfParameter
303  elseif ((parcat.eq.4).and.(parnum.eq.9).and.(typsurf.eq.1)) then ! SR
304  isec1(6)=176 ! indicatorOfParameter
305  elseif ((parcat.eq.2).and.(parnum.eq.17)) then ! EWSS
306  isec1(6)=180 ! indicatorOfParameter
307  elseif ((parcat.eq.2).and.(parnum.eq.18)) then ! NSSS
308  isec1(6)=181 ! indicatorOfParameter
309  elseif ((parcat.eq.3).and.(parnum.eq.4)) then ! ORO
310  isec1(6)=129 ! indicatorOfParameter
311  elseif ((parcat.eq.3).and.(parnum.eq.7)) then ! SDO
312  isec1(6)=160 ! indicatorOfParameter
313  elseif ((discipl.eq.2).and.(parcat.eq.0).and.(parnum.eq.0).and. &
314  (typsurf.eq.1)) then ! LSM
315  isec1(6)=172 ! indicatorOfParameter
316  else
317  print*,'***ERROR: undefined GRiB2 message found!',discipl, &
318  parcat,parnum,typsurf
319  endif
320 
321 #endif
322 
323  endif
324 
325  !HSO get the size and data of the values array
326  if (isec1(6).ne.-1) then
327  call grib_get_real4_array(igrib,'values',zsec4,iret)
328 #if defined WITH_CTBTO_PATCHES
329  zsec4=zsec4/conversion_factor
330 #endif
331  call grib_check(iret,gribfunction,griberrormsg)
332  endif
333 
334  !HSO get the required fields from section 2 in a gribex compatible manner
335  if (ifield.eq.1) then
336  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
337  isec2(2),iret)
338  call grib_check(iret,gribfunction,griberrormsg)
339  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
340  isec2(3),iret)
341  call grib_check(iret,gribfunction,griberrormsg)
342  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
343  isec2(12))
344  call grib_check(iret,gribfunction,griberrormsg)
345  ! CHECK GRID SPECIFICATIONS
346  if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
347  if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
348  if(isec2(12)/2-1.ne.nlev_ec) &
349  stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
350  endif ! ifield
351 
352  !HSO get the second part of the grid dimensions only from GRiB1 messages
353 #if defined WITH_CTBTO_PATCHES
354  if (gotgrid.eq.0) then
355 #else
356  if ((gribver.eq.1).and.(gotgrid.eq.0)) then
357 #endif
358  call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
359  xauxin,iret)
360  call grib_check(iret,gribfunction,griberrormsg)
361  call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
362  yauxin,iret)
363  call grib_check(iret,gribfunction,griberrormsg)
364  xaux=xauxin+real(nxshift)*dx
365  yaux=yauxin
366  xaux0=xlon0
367  yaux0=ylat0
368  if(xaux.lt.0.) xaux=xaux+360.
369  if(yaux.lt.0.) yaux=yaux+360.
370  if(xaux0.lt.0.) xaux0=xaux0+360.
371  if(yaux0.lt.0.) yaux0=yaux0+360.
372  if(abs(xaux-xaux0).gt.eps) &
373  stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
374  if(abs(yaux-yaux0).gt.eps) &
375  stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
376  gotgrid=1
377  endif ! gotGrid
378 
379  do j=0,nymin1
380  do i=0,nxfield-1
381  k=isec1(8)
382  if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
383  zsec4(nxfield*(ny-j-1)+i+1)
384  if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
385  zsec4(nxfield*(ny-j-1)+i+1)
386  if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
387  zsec4(nxfield*(ny-j-1)+i+1)
388  if(isec1(6).eq.133) then !! SPEC. HUMIDITY
389  qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
390  if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
391  qvh(i,j,nlev_ec-k+2,n) = 0.
392  ! this is necessary because the gridded data may contain
393  ! spurious negative values
394  endif
395  if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
396  zsec4(nxfield*(ny-j-1)+i+1)
397 
398  if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
399  zsec4(nxfield*(ny-j-1)+i+1)
400  if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
401  zsec4(nxfield*(ny-j-1)+i+1)
402  if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
403  zsec4(nxfield*(ny-j-1)+i+1)
404  if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
405  zsec4(nxfield*(ny-j-1)+i+1)
406  if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
407  zsec4(nxfield*(ny-j-1)+i+1)
408  if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
409  zsec4(nxfield*(ny-j-1)+i+1)
410  if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
411  zsec4(nxfield*(ny-j-1)+i+1)
412  if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
413  zsec4(nxfield*(ny-j-1)+i+1)
414  if(isec1(6).eq.142) then !! LARGE SCALE PREC.
415  lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
416  if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
417  endif
418  if(isec1(6).eq.143) then !! CONVECTIVE PREC.
419  convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
420  if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
421  endif
422  if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
423  zsec4(nxfield*(ny-j-1)+i+1)
424  if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
425  hflswitch=.true. ! Heat flux available
426  if(isec1(6).eq.176) then !! SOLAR RADIATION
427  ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
428  if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
429  endif
430  if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
431  zsec4(nxfield*(ny-j-1)+i+1)
432  if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
433  zsec4(nxfield*(ny-j-1)+i+1)
434  if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
435  (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true. ! stress available
436  !sec strswitch=.true.
437  if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
438  zsec4(nxfield*(ny-j-1)+i+1)/ga
439  if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
440  zsec4(nxfield*(ny-j-1)+i+1)
441  if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
442  zsec4(nxfield*(ny-j-1)+i+1)
443  if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
444  if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
445 
446  end do
447  end do
448 
449  call grib_release(igrib)
450  goto 10 !! READ NEXT LEVEL OR PARAMETER
451  !
452  ! CLOSING OF INPUT DATA FILE
453  !
454 
455 50 call grib_close_file(ifile)
456 
457  !error message if no fields found with correct first longitude in it
458  if (gotgrid.eq.0) then
459  print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
460  'messages'
461  stop
462  endif
463 
464  if(levdiff2.eq.0) then
465  iwmax=nlev_ec+1
466  do i=0,nxmin1
467  do j=0,nymin1
468  wwh(i,j,nlev_ec+1)=0.
469  end do
470  end do
471  endif
472 
473  ! For global fields, assign the leftmost data column also to the rightmost
474  ! data column; if required, shift whole grid by nxshift grid points
475  !*************************************************************************
476 
477  if (xglobal) then
478  call shift_field_0(ewss,nxfield,ny)
479  call shift_field_0(nsss,nxfield,ny)
480  call shift_field_0(oro,nxfield,ny)
481  call shift_field_0(excessoro,nxfield,ny)
482  call shift_field_0(lsm,nxfield,ny)
483  call shift_field(ps,nxfield,ny,1,1,2,n)
484  call shift_field(sd,nxfield,ny,1,1,2,n)
485  call shift_field(msl,nxfield,ny,1,1,2,n)
486  call shift_field(tcc,nxfield,ny,1,1,2,n)
487  call shift_field(u10,nxfield,ny,1,1,2,n)
488  call shift_field(v10,nxfield,ny,1,1,2,n)
489  call shift_field(tt2,nxfield,ny,1,1,2,n)
490  call shift_field(td2,nxfield,ny,1,1,2,n)
491  call shift_field(lsprec,nxfield,ny,1,1,2,n)
492  call shift_field(convprec,nxfield,ny,1,1,2,n)
493  call shift_field(sshf,nxfield,ny,1,1,2,n)
494  call shift_field(ssr,nxfield,ny,1,1,2,n)
495  call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
496  call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
497  call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
498  call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
499  call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
500  endif
501 
502  do i=0,nxmin1
503  do j=0,nymin1
504  surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
505  end do
506  end do
507 
508  if ((.not.hflswitch).or.(.not.strswitch)) then
509  write(*,*) 'WARNING: No flux data contained in GRIB file ', &
510  wfname(indj)
511 
512  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
513  ! As ECMWF has increased the model resolution, such that now the first model
514  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
515  ! (3rd model level in FLEXPART) for the profile method
516  !***************************************************************************
517 
518  do i=0,nxmin1
519  do j=0,nymin1
520  plev1=akz(3)+bkz(3)*ps(i,j,1,n)
521  pmean=0.5*(ps(i,j,1,n)+plev1)
522  tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
523  fu=-r_air*tv/ga/pmean
524  hlev1=fu*(plev1-ps(i,j,1,n)) ! HEIGTH OF FIRST MODEL LAYER
525  ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
526  fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
527  call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
528  tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
529  surfstr(i,j,1,n),sshf(i,j,1,n))
530  if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
531  if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
532  end do
533  end do
534  endif
535 
536 
537  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
538  ! level at the ground
539  ! Specific humidity is taken the same as at one level above
540  ! Temperature is taken as 2 m temperature
541  !**************************************************************************
542 
543  do i=0,nxmin1
544  do j=0,nymin1
545  uuh(i,j,1)=u10(i,j,1,n)
546  vvh(i,j,1)=v10(i,j,1,n)
547  qvh(i,j,1,n)=qvh(i,j,2,n)
548  tth(i,j,1,n)=tt2(i,j,1,n)
549  end do
550  end do
551 
552  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
553  if(iwmax.ne.nwz) stop 'READWIND: NWZ NOT CONSISTENT'
554 
555  return
556 888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
557  write(*,*) ' #### ',wfname(indj),' #### '
558  write(*,*) ' #### IS NOT GRIB FORMAT !!! #### '
559  stop 'Execution terminated'
560 999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
561  write(*,*) ' #### ',wfname(indj),' #### '
562  write(*,*) ' #### CANNOT BE OPENED !!! #### '
563  stop 'Execution terminated'
564 
565 end subroutine readwind_ecmwf
subroutine readwind_ecmwf(indj, n, uuh, vvh, wwh)
Definition: readwind.F90:22
subroutine shift_field(field, nxf, nyf, nzfmax, nzf, nmax, n)
Definition: shift_field.f90:22
subroutine pbl_profile(ps, td2m, zml1, t2m, tml1, u10m, uml1, stress, hf)
Definition: pbl_profile.f90:22
subroutine shift_field_0(field, nxf, nyf)
subroutine grib2check(igrib, isec1, conversion_factor)
Definition: grib2check.F90:1