FLEXPART CTBTO WO8
 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  use class_vtable
57 
58  implicit none
59 
60  !***********************************************************************
61  ! Subroutine Parameters: *
62  ! input *
63  ! indj indicates number of the wind field to be read in *
64  ! n temporal index for meteorological fields (1 to 3)*
65  ! uuh,vvh, wwh wind components in ECMWF model levels *
66  integer :: indj,n
67  real(kind=4) :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
68  real(kind=4) :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
69  real(kind=4) :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
70  !***********************************************************************
71 
72  !***********************************************************************
73  ! Local variables *
74  ! *
75  ! HSO variables for grib_api: *
76  ! ifile grib file to be opened and read in *
77  ! iret is a return code for successful or not open *
78  ! igrib grib edition number (whether GRIB1 or GRIB2) *
79  ! gribVer where info on igrib is kept, 1 for GRIB1 and *
80  ! 2 for GRIB2 *
81  ! parCat parameter category ( = number) , how FLEXPART *
82  ! identifies fields *
83  ! parNum parameter number by product discipline and *
84  ! parameter category *
85  ! typSurf type of first fixed surface *
86  ! valSurf level *
87  ! discipl discipline of processed data contained in a *
88  ! GRIB message *
89  integer :: ifile
90  integer :: iret
91  integer :: igrib
92  integer :: gribver,parcat,parnum,typsurf,valsurf,discipl
93 #if defined WITH_CTBTO_PATCHES
94  integer :: parid
95 #endif
96  integer :: gotgrid
97 
98 
99  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
100  !!!! Vtable related variables
101 
102  ! Paths to Vtables (current implementation will assume they are in the cwd
103  ! This is a crappy place for these parameters. Need to move them.
104  character(LEN=255), parameter :: vtable_ecmwf_grib1_path = &
105  "Vtable_ecmwf_grib1", &
106  vtable_ecmwf_grib2_path = &
107  "Vtable_ecmwf_grib2", &
108  vtable_ecmwf_grib1_2_path = &
109  "Vtable_ecmwf_grib1_2"
110 
111  integer :: gribfile_type
112  integer :: current_grib_level ! This "was" isec1(8) in previous version
113  character(len=255) :: gribfile_name
114  character(len=255) :: vtable_path
115  character(len=15) :: fpname
116 
117  type(vtable) :: my_vtable ! unallocated
118  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
119 
120 
121 
122  !
123  ! Variables and arrays for grib decoding: *
124  ! dimension of isec2 at least (22+n), where n is the number of *
125  ! parallels or meridians in a quasi-regular (reduced) Gaussian or *
126  ! lat/long grid *
127  ! dimension of zsec2 at least (10+nn), where nn is the number of *
128  ! vertical coordinate parameters *
129  ! *
130  ! isec1 grib definition section (version, center, ... *
131  ! isec2 grid description section *
132  ! zsec4 the binary data section *
133  ! xaux
134  ! yaux
135  ! xauxin
136  ! yauxin
137  ! xaux0 auxiliary variable for xlon0 *
138  ! yaux0 auxiliary variable for xlat0 *
139  ! nsss NS surface stress *
140  ! ewss EW surface stress *
141  ! plev1 pressure of the first model layer *
142  ! pmean mean pressure between ?? *
143  ! tv virtual temperature *
144  ! fu for conversion from pressure to meters *
145  ! hlev1 height of the first model layer *
146  ! ff10m wind speed at 10 m *
147  ! fflev1 wind speed at the first model layere *
148  ! hflswitch logical variable to check existence of flux data *
149  ! strswitch logical variable to check existence of stress *
150  ! data *
151 
152  integer :: isec1(56),isec2(22+nxmax+nymax)
153  real(kind=4) :: zsec4(jpunp)
154  real(kind=4) :: xaux,yaux,xaux0,yaux0
155 #if defined WITH_CTBTO_PATCHES
156  real(kind=4) :: conversion_factor
157 #endif
158  real(kind=8) :: xauxin,yauxin
159  real(kind=4) :: nsss(0:nxmax-1,0:nymax-1),ewss(0:nxmax-1,0:nymax-1)
160  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
161  logical :: hflswitch,strswitch
162 
163  ! Other variables:
164  ! i,j,k loop control indices in each direction *
165  ! levdiff2 number of model levels - number of model levels *
166  ! for the staggered vertical wind +1 *
167  ! ifield index to control field read in *
168  ! iumax
169  ! iwmax
170  integer :: i,j,k,levdiff2,ifield,iumax,iwmax
171  !***********************************************************************
172 
173  !***********************************************************************
174  ! Local constants *
175  real,parameter :: eps=1.e-4
176  !HSO grib api error messages:
177  character(len=24) :: griberrormsg = 'Error reading grib file'
178  character(len=20) :: gribfunction = 'readwind'
179  !***********************************************************************
180 
181  !***********************************************************************
182  ! Global variables *
183  ! from par_mod and com_mod *
184  ! wfname File name of data to be read in *
185  ! nx,ny,nuvz,nwz expected field dimensions *
186  ! nxfield same as nx but for limited area fields *
187  ! nxmin1,nymin1 nx-1, ny-1, respectively *
188  ! nxmax,nymax maximum dimension of wind fields in x and y *
189  ! direction, respectively *
190  ! nuvzmax,nwzmax maximum dimension of (u,v) and (w) wind fields in z
191  ! direction (for fields on eta levels) *
192  ! nlev_ec number of vertical levels ecmwf model *
193  ! u10,v10 wind fields at 10 m *
194  ! tt2,td2 temperature and dew point temperature at 2 m *
195  ! tth,qvh temperature and specific humidity on original *
196  ! eta levels *
197  ! ps surface pressure *
198  ! sd snow depth (but not used afterwards!) *
199  ! msl mean sea level pressure *
200  ! ttc total cloud cover *
201  ! lsprec large scale precipitation *
202  ! convprec convective precipitation *
203  ! sshf sensible heat flux *
204  ! ssr solar radiation *
205  ! surfstr surface stress *
206  ! oro orography *
207  ! excessoro standard deviation of orography *
208  ! lsm land sea mask *
209  ! r_air individual gas constant for dry air [J/kg/K] *
210  ! ga gravity acceleration of earth [m/s**2] *
211  !***********************************************************************
212 
213 !-----------------------------------------------------------------------------
214 
215 
216  hflswitch=.false.
217  strswitch=.false.
218  levdiff2=nlev_ec-nwz+1
219  iumax=0
220  iwmax=0
221 
222  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
223  !!!!!!!!!!!!!!!!!!! VTABLE code
224  !!!!!!! Vtable choice
225  gribfile_name = path(3)(1:length(3))//trim(wfname(indj))
226  print *, 'gribfile_name: ', gribfile_name
227 
228  gribfile_type = vtable_detect_gribfile_type( gribfile_name )
229 
230  print *, 'gribfile_type: ', gribfile_type
231 
232  if (gribfile_type .eq. vtable_gribfile_type_ecmwf_grib1) then
233  vtable_path = vtable_ecmwf_grib1_path
234  else if (gribfile_type .eq. vtable_gribfile_type_ecmwf_grib2) then
235  vtable_path = vtable_ecmwf_grib2_path
236  else if (gribfile_type .eq. vtable_gribfile_type_ecmwf_grib1_2) then
237  vtable_path = vtable_ecmwf_grib1_2_path
238  else
239  print *, 'Unsupported gribfile_type: ', gribfile_type
240  stop
241  endif
242 
243 
244  ! Load the Vtable into 'my_vtable'
245  print *, 'Loading Vtable: ', vtable_path
246  call vtable_load_by_name(vtable_path, my_vtable)
247  print *, 'Vtable Initialized: ', my_vtable%initialized
248  print *, 'Vtable num_entries: ', my_vtable%num_entries
249  !!!!!!!!!!!!!!!!!!! VTABLE code
250  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
251 
252 
253 
254 
255 
256  !
257  ! OPENING OF DATA FILE (GRIB CODE)
258  !
259 5 call grib_open_file(ifile,path(3)(1:length(3)) &
260  //trim(wfname(indj)),'r',iret)
261  if (iret.ne.grib_success) then
262  goto 888 ! ERROR DETECTED
263  endif
264  !turn on support for multi fields messages */
265  !call grib_multi_support_on
266 
267  gotgrid=0
268  ifield=0
269 10 ifield=ifield+1
270  !
271  ! GET NEXT FIELDS
272  !
273  call grib_new_from_file(ifile,igrib,iret)
274  if (iret.eq.grib_end_of_file) then
275  goto 50 ! EOF DETECTED
276  elseif (iret.ne.grib_success) then
277  goto 888 ! ERROR DETECTED
278  endif
279 
280  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
281  !!!!!!!!!!!!!!!!!!! VTABLE code
282  ! Get the fpname
283  fpname = vtable_get_fpname(igrib, my_vtable)
284  !print *, 'fpname: ', trim(fpname)
285 
286 
287  !!!!!!!!!!!!!!!!!!! VTABLE code
288  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
289 
290 #if defined WITH_CTBTO_PATCHES
291  conversion_factor=1.0
292 #endif
293 
294  !first see if we read GRIB1 or GRIB2
295  call grib_get_int(igrib,'editionNumber',gribver,iret)
296  call grib_check(iret,gribfunction,griberrormsg)
297 
298  !!!!!!! DJM - candidate for consolidation
299  if (gribver.eq.1) then ! GRIB Edition 1
300 
301  call grib_get_int(igrib,'level', current_grib_level,iret)
302  call grib_check(iret,gribfunction,griberrormsg)
303 
304 
305  else ! GRIB Edition 2
306 
307 #if defined WITH_CTBTO_PATCHES
308  !print *, 'fpname: ', fpname
309  call grib2check(igrib, fpname, conversion_factor)
310 #endif
311 
312  call grib_get_int(igrib,'level', current_grib_level,iret)
313  call grib_check(iret,gribfunction,griberrormsg)
314 
315 
316  endif ! END IF Grib Version selection
317 
318 
319 
320  !HSO get the size and data of the values array
321  if (trim(fpname) .ne. 'None') then
322  call grib_get_real4_array(igrib,'values',zsec4,iret)
323 #if defined WITH_CTBTO_PATCHES
324  zsec4=zsec4/conversion_factor
325 #endif
326  call grib_check(iret,gribfunction,griberrormsg)
327  endif
328 
329  !HSO get the required fields from section 2 in a gribex compatible manner
330  if (ifield.eq.1) then
331  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
332  isec2(2),iret)
333  call grib_check(iret,gribfunction,griberrormsg)
334  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
335  isec2(3),iret)
336  call grib_check(iret,gribfunction,griberrormsg)
337  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
338  isec2(12))
339  call grib_check(iret,gribfunction,griberrormsg)
340  ! CHECK GRID SPECIFICATIONS
341  if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
342  if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
343  if(isec2(12)/2-1.ne.nlev_ec) &
344  stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
345  endif ! ifield
346 
347  !HSO get the second part of the grid dimensions only from GRiB1 messages
348 #if defined WITH_CTBTO_PATCHES
349 ! if (gotGrid.eq.0) then
350 #else
351 ! if ((gribVer.eq.1).and.(gotGrid.eq.0)) then
352 #endif
353  call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
354  xauxin,iret)
355  call grib_check(iret,gribfunction,griberrormsg)
356  call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
357  yauxin,iret)
358  call grib_check(iret,gribfunction,griberrormsg)
359  xaux=xauxin+real(nxshift)*dx
360  yaux=yauxin
361  xaux0=xlon0
362  yaux0=ylat0
363  if(xaux.lt.0.) xaux=xaux+360.
364  if(yaux.lt.0.) yaux=yaux+360.
365  if(xaux0.lt.0.) xaux0=xaux0+360.
366  if(yaux0.lt.0.) yaux0=yaux0+360.
367  if(abs(xaux-xaux0).gt.eps) &
368  stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
369  if(abs(yaux-yaux0).gt.eps) &
370  stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
371  gotgrid=1
372 ! endif ! gotGrid
373 
374 
375 k = current_grib_level
376 
377 !! TEMPERATURE
378 if(trim(fpname) .eq. 'TT') then
379  do j=0,nymin1
380  do i=0,nxfield-1
381  tth(i,j,nlev_ec-k+2,n) = zsec4(nxfield*(ny-j-1)+i+1)
382  enddo
383  enddo
384 endif
385 
386 
387 !! U VELOCITY
388 if(trim(fpname) .eq. 'UU') then
389  iumax=max(iumax,nlev_ec-k+1)
390  do j=0,nymin1
391  do i=0,nxfield-1
392  uuh(i,j,nlev_ec-k+2) = zsec4(nxfield*(ny-j-1)+i+1)
393  enddo
394  enddo
395 endif
396 
397 
398 !! V VELOCITY
399 if(trim(fpname) .eq. 'VV') then
400  do j=0,nymin1
401  do i=0,nxfield-1
402  vvh(i,j,nlev_ec-k+2) = zsec4(nxfield*(ny-j-1)+i+1)
403  enddo
404  enddo
405 endif
406 
407 !! W VELOCITY
408 if(trim(fpname) .eq. 'ETADOT') then
409  iwmax=max(iwmax,nlev_ec-k+1)
410  do j=0,nymin1
411  do i=0,nxfield-1
412  wwh(i,j,nlev_ec-k+1) = zsec4(nxfield*(ny-j-1)+i+1)
413  enddo
414  enddo
415 endif
416 
417 !! SPEC HUMIDITY
418 if(trim(fpname) .eq. 'QV') then
419  do j=0,nymin1
420  do i=0,nxfield-1
421  qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
422  if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) qvh(i,j,nlev_ec-k+2,n) = 0.
423  ! this is necessary because the gridded data may contain
424  ! spurious negative values
425  enddo
426  enddo
427 endif
428 
429 !! SURF. PRESS.
430 if(trim(fpname) .eq. 'PS') then
431  do j=0,nymin1
432  do i=0,nxfield-1
433  ps(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
434  enddo
435  enddo
436 endif
437 
438 !! SNOW DEPTH
439 if(trim(fpname) .eq. 'SD') then
440  do j=0,nymin1
441  do i=0,nxfield-1
442  sd(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
443  enddo
444  enddo
445 endif
446 
447 !! SEA LEVEL PRESS.
448 if(trim(fpname) .eq. 'MSL') then
449  do j=0,nymin1
450  do i=0,nxfield-1
451  msl(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
452  enddo
453  enddo
454 endif
455 
456 !! CLOUD COVER
457 if(trim(fpname) .eq. 'TCC') then
458  do j=0,nymin1
459  do i=0,nxfield-1
460  tcc(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
461  enddo
462  enddo
463 endif
464 
465 !! 10 M U VELOCITY
466 if(trim(fpname) .eq. 'U10') then
467  do j=0,nymin1
468  do i=0,nxfield-1
469  u10(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
470  enddo
471  enddo
472 endif
473 
474 !! 10 M V VELOCITY
475 if(trim(fpname) .eq. 'V10') then
476  do j=0,nymin1
477  do i=0,nxfield-1
478  v10(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
479  enddo
480  enddo
481 endif
482 
483 !! 2 M TEMPERATURE
484 if(trim(fpname) .eq. 'T2') then
485  do j=0,nymin1
486  do i=0,nxfield-1
487  tt2(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
488  enddo
489  enddo
490 endif
491 
492 !! 2 M DEW POINT
493 if(trim(fpname) .eq. 'TD2') then
494  do j=0,nymin1
495  do i=0,nxfield-1
496  td2(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
497  enddo
498  enddo
499 endif
500 
501 !! LARGE SCALE PREC.
502 if(trim(fpname) .eq. 'LSPREC') then
503  do j=0,nymin1
504  do i=0,nxfield-1
505  lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
506  if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
507  enddo
508  enddo
509 endif
510 
511 !! CONVECTIVE PREC.
512 if(trim(fpname) .eq. 'CONVPREC') then
513  do j=0,nymin1
514  do i=0,nxfield-1
515  convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
516  if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
517  enddo
518  enddo
519 endif
520 
521 !! SENS. HEAT FLUX
522 if(trim(fpname) .eq. 'SHF') then
523  do j=0,nymin1
524  do i=0,nxfield-1
525  sshf(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
526  if (zsec4(nxfield*(ny-j-1)+i+1).ne.0.) hflswitch=.true. ! Heat flux available
527  enddo
528  enddo
529 endif
530 
531 !! SOLAR RADIATION
532 if(trim(fpname) .eq. 'SR') then
533  do j=0,nymin1
534  do i=0,nxfield-1
535  ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
536  if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
537  enddo
538  enddo
539 endif
540 
541 !! EW SURFACE STRESS
542 if(trim(fpname) .eq. 'EWSS') then
543  do j=0,nymin1
544  do i=0,nxfield-1
545  ewss(i,j)=zsec4(nxfield*(ny-j-1)+i+1)
546  if (zsec4(nxfield*(ny-j-1)+i+1).ne.0.) strswitch=.true. ! stress available
547  enddo
548  enddo
549 endif
550 
551 !! NS SURFACE STRESS
552 if(trim(fpname) .eq. 'NSSS') then
553  do j=0,nymin1
554  do i=0,nxfield-1
555  nsss(i,j)=zsec4(nxfield*(ny-j-1)+i+1)
556  if (zsec4(nxfield*(ny-j-1)+i+1).ne.0.) strswitch=.true. ! stress available
557  enddo
558  enddo
559 endif
560 
561 !! ECMWF OROGRAPHY
562 if(trim(fpname) .eq. 'ORO') then
563  do j=0,nymin1
564  do i=0,nxfield-1
565  oro(i,j)=zsec4(nxfield*(ny-j-1)+i+1)/ga
566  enddo
567  enddo
568 endif
569 
570 !! ECMWF OROGRAPHY
571 if(trim(fpname) .eq. 'EXCESSORO') then
572  do j=0,nymin1
573  do i=0,nxfield-1
574  excessoro(i,j)=zsec4(nxfield*(ny-j-1)+i+1)
575  enddo
576  enddo
577 endif
578 
579 !! ECMWF LAND SEA MASK
580 if(trim(fpname) .eq. 'LSM') then
581  do j=0,nymin1
582  do i=0,nxfield-1
583  lsm(i,j)=zsec4(nxfield*(ny-j-1)+i+1)
584  enddo
585  enddo
586 endif
587 
588 
589  call grib_release(igrib)
590  goto 10 !! READ NEXT LEVEL OR PARAMETER
591  !
592  ! CLOSING OF INPUT DATA FILE
593  !
594 
595 50 call grib_close_file(ifile)
596 
597  !error message if no fields found with correct first longitude in it
598  if (gotgrid.eq.0) then
599  print*,'***ERROR: no fields found with correct first longitude'// &
600  'messages'
601  stop
602  endif
603 
604  if(levdiff2.eq.0) then
605  iwmax=nlev_ec+1
606  do i=0,nxmin1
607  do j=0,nymin1
608  wwh(i,j,nlev_ec+1)=0.
609  end do
610  end do
611  endif
612 
613  ! For global fields, assign the leftmost data column also to the rightmost
614  ! data column; if required, shift whole grid by nxshift grid points
615  !*************************************************************************
616 
617  if (xglobal) then
618  call shift_field_0(ewss,nxfield,ny)
619  call shift_field_0(nsss,nxfield,ny)
620  call shift_field_0(oro,nxfield,ny)
621  call shift_field_0(excessoro,nxfield,ny)
622  call shift_field_0(lsm,nxfield,ny)
623  call shift_field(ps,nxfield,ny,1,1,2,n)
624  call shift_field(sd,nxfield,ny,1,1,2,n)
625  call shift_field(msl,nxfield,ny,1,1,2,n)
626  call shift_field(tcc,nxfield,ny,1,1,2,n)
627  call shift_field(u10,nxfield,ny,1,1,2,n)
628  call shift_field(v10,nxfield,ny,1,1,2,n)
629  call shift_field(tt2,nxfield,ny,1,1,2,n)
630  call shift_field(td2,nxfield,ny,1,1,2,n)
631  call shift_field(lsprec,nxfield,ny,1,1,2,n)
632  call shift_field(convprec,nxfield,ny,1,1,2,n)
633  call shift_field(sshf,nxfield,ny,1,1,2,n)
634  call shift_field(ssr,nxfield,ny,1,1,2,n)
635  call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
636  call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
637  call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
638  call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
639  call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
640  endif
641 
642  do i=0,nxmin1
643  do j=0,nymin1
644  surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
645  end do
646  end do
647 
648  if ((.not.hflswitch).or.(.not.strswitch)) then
649  write(*,*) 'WARNING: No flux data contained in GRIB file ', &
650  wfname(indj)
651 
652  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
653  ! As ECMWF has increased the model resolution, such that now the first model
654  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
655  ! (3rd model level in FLEXPART) for the profile method
656  !***************************************************************************
657 
658  do i=0,nxmin1
659  do j=0,nymin1
660  plev1=akz(3)+bkz(3)*ps(i,j,1,n)
661  pmean=0.5*(ps(i,j,1,n)+plev1)
662  tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
663  fu=-r_air*tv/ga/pmean
664  hlev1=fu*(plev1-ps(i,j,1,n)) ! HEIGTH OF FIRST MODEL LAYER
665  ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
666  fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
667  call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
668  tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
669  surfstr(i,j,1,n),sshf(i,j,1,n))
670  if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
671  if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
672  end do
673  end do
674  endif
675 
676 
677  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
678  ! level at the ground
679  ! Specific humidity is taken the same as at one level above
680  ! Temperature is taken as 2 m temperature
681  !**************************************************************************
682 
683  do i=0,nxmin1
684  do j=0,nymin1
685  uuh(i,j,1)=u10(i,j,1,n)
686  vvh(i,j,1)=v10(i,j,1,n)
687  qvh(i,j,1,n)=qvh(i,j,2,n)
688  tth(i,j,1,n)=tt2(i,j,1,n)
689  end do
690  end do
691 
692  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
693  if(iwmax.ne.nwz) stop 'READWIND: NWZ NOT CONSISTENT'
694 
695  return
696 888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
697  write(*,*) ' #### ',wfname(indj),' #### '
698  write(*,*) ' #### IS NOT GRIB FORMAT !!! #### '
699  stop 'Execution terminated'
700 999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
701  write(*,*) ' #### ',wfname(indj),' #### '
702  write(*,*) ' #### CANNOT BE OPENED !!! #### '
703  stop 'Execution terminated'
704 
705 end subroutine readwind_ecmwf
character(len=15) function, public vtable_get_fpname(igrib, vtable_object)
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 grib2check(igrib, fpname, conversion_factor)
Definition: grib2check.F90:1
subroutine pbl_profile(ps, td2m, zml1, t2m, tml1, u10m, uml1, stress, hf)
Definition: pbl_profile.f90:22
subroutine, public vtable_load_by_name(vtable_name, the_vtable_data)
subroutine shift_field_0(field, nxf, nyf)
integer function, public vtable_detect_gribfile_type(gribfilename)