source: flexpart.git/src/readwind_ecmwf.f90

bugfixes+enhancements
Last change on this file was dba4221, checked in by Pirmin Kaufmann <pirmin.kaufmann@…>, 18 months ago

Bugfixes:

  • options/SPECIES/SPECIES_009: corrected wrong number format, replaced comma by decimal point
  • options/SPECIES/SPECIES_028: corrected wrong number format, moved sign of exponent to after the E
  • options/SPECIES/specoverview.f90: added namelist parameters that appear in SPECIES files but were missing here
  • src/FLEXPART.f90: replaced compiler-specific command line argument routines by standard Fortran intrinsic routines
  • src/FLEXPART_MPI.f90: ditto
  • src/gridcheck_ecmwf.f90: corrected handling of vertical levels when input files do not contain uppermost layers
  • src/gridcheck_nests.f90: ditto
  • src/readwind_ecmwf.f90: corrected handling of vertical levels when input files do not contain uppermost layers
  • readwind_ecmwf_mpi.f90: ditto

Code enhancements:

  • options/OUTGRID: added comments describing contents
  • options/SPECIES/SPECIES_*: aligned comments
  • options/SPECIES/specoverview.f90: removed commented lines, rectified lines indenting
  • src/FLEXPART.f90: rectified lines indenting, updated date in version string
  • src/FLEXPART_MPI.f90: ditto, and realigned code with src/FLEXPART.f90
  • src/gridcheck_*.f90: added code to write out name of file before it is opened (helps a lot when an input file causes troubles)
  • src/par_mod.f90: added comment explaining relevance of nuvzmax for GRIB input
  • src/readreleases.f90: write out warning if too few particles are used to randomize release
  • src/readspecies.f90: write out name of SPECIES file before it is read
  • src/readwind_*.f90: write out name of input file before opening it
  • src/writeheader_txt.f90: removed wrong comment
  • 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!
105  write(*,*)   'Reading: '//path(3)(1:length(3)) &
106       //trim(wfname(indj))
1075 call grib_open_file(ifile,path(3)(1:length(3)) &
108       //trim(wfname(indj)),'r',iret)
109  if (iret.ne.GRIB_SUCCESS) then
110    goto 888   ! ERROR DETECTED
111  endif
112!turn on support for multi fields messages */
113!call grib_multi_support_on
114
115  gotGrid=0
116  ifield=0
11710 ifield=ifield+1
118!
119! GET NEXT FIELDS
120!
121  call grib_new_from_file(ifile,igrib,iret)
122  if (iret.eq.GRIB_END_OF_FILE)  then
123    goto 50    ! EOF DETECTED
124  elseif (iret.ne.GRIB_SUCCESS) then
125    goto 888   ! ERROR DETECTED
126  endif
127
128!first see if we read GRIB1 or GRIB2
129  call grib_get_int(igrib,'editionNumber',gribVer,iret)
130  call grib_check(iret,gribFunction,gribErrorMsg)
131
132  if (gribVer.eq.1) then ! GRIB Edition 1
133
134!print*,'GRiB Edition 1'
135!read the grib2 identifiers
136    call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
137    call grib_check(iret,gribFunction,gribErrorMsg)
138    call grib_get_int(igrib,'level',isec1(8),iret)
139    call grib_check(iret,gribFunction,gribErrorMsg)
140
141!change code for etadot to code for omega
142    if (isec1(6).eq.77) then
143      isec1(6)=135
144    endif
145
146    conversion_factor=1.
147
148  else
149
150!print*,'GRiB Edition 2'
151!read the grib2 identifiers
152    call grib_get_int(igrib,'discipline',discipl,iret)
153    call grib_check(iret,gribFunction,gribErrorMsg)
154    call grib_get_int(igrib,'parameterCategory',parCat,iret)
155    call grib_check(iret,gribFunction,gribErrorMsg)
156    call grib_get_int(igrib,'parameterNumber',parNum,iret)
157    call grib_check(iret,gribFunction,gribErrorMsg)
158    call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
159    call grib_check(iret,gribFunction,gribErrorMsg)
160    call grib_get_int(igrib,'level',valSurf,iret)
161    call grib_check(iret,gribFunction,gribErrorMsg)
162    call grib_get_int(igrib,'paramId',parId,iret)
163    call grib_check(iret,gribFunction,gribErrorMsg)
164
165!print*,discipl,parCat,parNum,typSurf,valSurf
166
167!convert to grib1 identifiers
168    isec1(6)=-1
169    isec1(7)=-1
170    isec1(8)=-1
171    isec1(8)=valSurf     ! level
172    conversion_factor=1.
173    if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
174      isec1(6)=130         ! indicatorOfParameter
175    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
176      isec1(6)=131         ! indicatorOfParameter
177    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
178      isec1(6)=132         ! indicatorOfParameter
179    elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
180      isec1(6)=133         ! indicatorOfParameter
181! ESO Cloud water is in a) fields CLWC and CIWC, *or* b) field QC
182    elseif ((parCat.eq.1).and.(parNum.eq.83).and.(typSurf.eq.105)) then ! clwc
183      isec1(6)=246         ! indicatorOfParameter
184    elseif ((parCat.eq.1).and.(parNum.eq.84).and.(typSurf.eq.105)) then ! ciwc
185      isec1(6)=247         ! indicatorOfParameter
186! ESO qc(=clwc+ciwc):
187    elseif ((parCat.eq.201).and.(parNum.eq.31).and.(typSurf.eq.105)) then ! qc
188      isec1(6)=201031         ! indicatorOfParameter
189    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
190      isec1(6)=134         ! indicatorOfParameter
191    elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
192      isec1(6)=135         ! indicatorOfParameter
193    elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually eta dot
194      isec1(6)=135         ! indicatorOfParameter
195    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
196      isec1(6)=151         ! indicatorOfParameter
197    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
198      isec1(6)=165         ! indicatorOfParameter
199    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
200      isec1(6)=166         ! indicatorOfParameter
201    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
202      isec1(6)=167         ! indicatorOfParameter
203    elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
204      isec1(6)=168         ! indicatorOfParameter
205    elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
206      isec1(6)=141         ! indicatorOfParameter
207      conversion_factor=1000.
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 ! CP
213      isec1(6)=143         ! indicatorOfParameter
214      conversion_factor=1000.
215    elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
216      isec1(6)=146         ! indicatorOfParameter
217    elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
218      isec1(6)=176         ! indicatorOfParameter
219!    elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS --wrong
220    elseif ((parCat.eq.2).and.(parNum.eq.38) .or. parId .eq. 180) then ! EWSS --correct
221      isec1(6)=180         ! indicatorOfParameter
222!    elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS --wrong
223    elseif ((parCat.eq.2).and.(parNum.eq.37) .or. parId .eq. 181) then ! NSSS --correct
224      isec1(6)=181         ! indicatorOfParameter
225    elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
226      isec1(6)=129         ! indicatorOfParameter
227    elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO
228      isec1(6)=160         ! indicatorOfParameter
229    elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
230         (typSurf.eq.1)) then ! LSM
231      isec1(6)=172         ! indicatorOfParameter
232    elseif (parNum.eq.152) then
233      isec1(6)=152         ! avoid warning for lnsp   
234    else
235      print*,'***WARNING: undefined GRiB2 message found!',discipl, &
236           parCat,parNum,typSurf
237    endif
238    if(parId .ne. isec1(6) .and. parId .ne. 77) then
239      write(*,*) 'parId',parId, 'isec1(6)',isec1(6)
240!    stop
241    endif
242
243  endif
244
245!HSO  get the size and data of the values array
246  if (isec1(6).ne.-1) then
247    call grib_get_real4_array(igrib,'values',zsec4,iret)
248    call grib_check(iret,gribFunction,gribErrorMsg)
249  endif
250
251!HSO  get the required fields from section 2 in a gribex compatible manner
252  if (ifield.eq.1) then
253    call grib_get_int(igrib,'numberOfPointsAlongAParallel',isec2(2),iret)
254    call grib_check(iret,gribFunction,gribErrorMsg)
255    call grib_get_int(igrib,'numberOfPointsAlongAMeridian',isec2(3),iret)
256    call grib_check(iret,gribFunction,gribErrorMsg)
257    call grib_get_int(igrib,'numberOfVerticalCoordinateValues',isec2(12))
258    call grib_check(iret,gribFunction,gribErrorMsg)
259! CHECK GRID SPECIFICATIONS
260    if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
261    if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
262    if(isec2(12)/2-1.ne.nlev_ec) &
263         stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
264  endif ! ifield
265
266!HSO  get the second part of the grid dimensions only from GRiB1 messages
267  if (isec1(6) .eq. 167 .and. (gotGrid.eq.0)) then
268    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
269         xauxin,iret)
270    call grib_check(iret,gribFunction,gribErrorMsg)
271    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
272         yauxin,iret)
273    call grib_check(iret,gribFunction,gribErrorMsg)
274    if (xauxin.gt.180.) xauxin=xauxin-360.0
275    if (xauxin.lt.-180.) xauxin=xauxin+360.0
276
277    xaux=xauxin+real(nxshift)*dx
278    yaux=yauxin
279    if (xaux.gt.180.) xaux=xaux-360.0
280    if(abs(xaux-xlon0).gt.eps) &
281         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
282    if(abs(yaux-ylat0).gt.eps) &
283         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
284    gotGrid=1
285  endif ! gotGrid
286
287  do j=0,nymin1
288    do i=0,nxfield-1
289      k=isec1(8)
290      if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
291           zsec4(nxfield*(ny-j-1)+i+1)
292      if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
293           zsec4(nxfield*(ny-j-1)+i+1)
294      if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
295           zsec4(nxfield*(ny-j-1)+i+1)
296      if(isec1(6).eq.133) then                      !! SPEC. HUMIDITY
297        qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
298        if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
299             qvh(i,j,nlev_ec-k+2,n) = 0.
300!        this is necessary because the gridded data may contain
301!        spurious negative values
302      endif
303      if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
304           zsec4(nxfield*(ny-j-1)+i+1)
305
306      if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
307           zsec4(nxfield*(ny-j-1)+i+1)
308      if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
309           zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
310      if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
311           zsec4(nxfield*(ny-j-1)+i+1)
312      if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
313           zsec4(nxfield*(ny-j-1)+i+1)
314      if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
315           zsec4(nxfield*(ny-j-1)+i+1)
316      if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
317           zsec4(nxfield*(ny-j-1)+i+1)
318      if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
319           zsec4(nxfield*(ny-j-1)+i+1)
320      if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
321           zsec4(nxfield*(ny-j-1)+i+1)
322      if(isec1(6).eq.142) then                      !! LARGE SCALE PREC.
323        lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
324        if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
325      endif
326      if(isec1(6).eq.143) then                      !! CONVECTIVE PREC.
327        convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)/conversion_factor
328        if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
329      endif
330      if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
331           zsec4(nxfield*(ny-j-1)+i+1)
332      if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
333           hflswitch=.true.    ! Heat flux available
334      if(isec1(6).eq.176) then                      !! SOLAR RADIATION
335        ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
336        if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
337      endif
338      if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
339           zsec4(nxfield*(ny-j-1)+i+1)
340      if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
341           zsec4(nxfield*(ny-j-1)+i+1)
342      if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
343           (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
344!sec        strswitch=.true.
345      if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
346           zsec4(nxfield*(ny-j-1)+i+1)/ga
347      if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
348           zsec4(nxfield*(ny-j-1)+i+1)
349      if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
350           zsec4(nxfield*(ny-j-1)+i+1)
351      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
352      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
353!ZHG READING CLOUD FIELDS ASWELL
354! ESO TODO: add check for if one of clwc/ciwc missing (error),
355! also if all 3 cw fields present, use qc and disregard the others
356      if(isec1(6).eq.246) then  !! CLWC  Cloud liquid water content [kg/kg]
357        clwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
358        readclouds=.true.
359        sumclouds=.false.
360      endif
361      if(isec1(6).eq.247) then  !! CIWC  Cloud ice water content
362        ciwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
363      endif
364!ZHG end
365!ESO read qc (=clwc+ciwc)
366      if(isec1(6).eq.201031) then  !! QC  Cloud liquid water content [kg/kg]
367        clwch(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
368        readclouds=.true.
369        sumclouds=.true.
370      endif
371
372    end do
373  end do
374
375  call grib_release(igrib)
376  goto 10                      !! READ NEXT LEVEL OR PARAMETER
377!
378! CLOSING OF INPUT DATA FILE
379!
380
38150 call grib_close_file(ifile)
382
383!error message if no fields found with correct first longitude in it
384  if (gotGrid.eq.0) then
385    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
386         'messages'
387    stop
388  endif
389
390!  if(levdiff2.eq.0) then
391    iwmax=iwmax+1
392    do i=0,nxmin1
393      do j=0,nymin1
394        wwh(i,j,iwmax)=0.
395      end do
396    end do
397!  endif
398
399! For global fields, assign the leftmost data column also to the rightmost
400! data column; if required, shift whole grid by nxshift grid points
401!*************************************************************************
402
403  if (xglobal) then
404    call shift_field_0(ewss,nxfield,ny)
405    call shift_field_0(nsss,nxfield,ny)
406    call shift_field_0(oro,nxfield,ny)
407    call shift_field_0(excessoro,nxfield,ny)
408    call shift_field_0(lsm,nxfield,ny)
409    call shift_field(ps,nxfield,ny,1,1,2,n)
410    call shift_field(sd,nxfield,ny,1,1,2,n)
411    call shift_field(msl,nxfield,ny,1,1,2,n)
412    call shift_field(tcc,nxfield,ny,1,1,2,n)
413    call shift_field(u10,nxfield,ny,1,1,2,n)
414    call shift_field(v10,nxfield,ny,1,1,2,n)
415    call shift_field(tt2,nxfield,ny,1,1,2,n)
416    call shift_field(td2,nxfield,ny,1,1,2,n)
417    call shift_field(lsprec,nxfield,ny,1,1,2,n)
418    call shift_field(convprec,nxfield,ny,1,1,2,n)
419    call shift_field(sshf,nxfield,ny,1,1,2,n)
420    call shift_field(ssr,nxfield,ny,1,1,2,n)
421    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
422    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
423    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
424    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
425    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
426!ZHG
427    call shift_field(clwch,nxfield,ny,nuvzmax,nuvz,2,n)
428    if (.not.sumclouds) call shift_field(ciwch,nxfield,ny,nuvzmax,nuvz,2,n)
429!ZHG end
430
431  endif
432
433  do i=0,nxmin1
434    do j=0,nymin1
435      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
436    end do
437  end do
438
439  if ((.not.hflswitch).or.(.not.strswitch)) then
440    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
441         wfname(indj)
442
443! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
444! As ECMWF has increased the model resolution, such that now the first model
445! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
446! (3rd model level in FLEXPART) for the profile method
447!***************************************************************************
448
449    do i=0,nxmin1
450      do j=0,nymin1
451        plev1=akz(3)+bkz(3)*ps(i,j,1,n)
452        pmean=0.5*(ps(i,j,1,n)+plev1)
453        tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
454        fu=-r_air*tv/ga/pmean
455        hlev1=fu*(plev1-ps(i,j,1,n))   ! HEIGTH OF FIRST MODEL LAYER
456        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
457        fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
458        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
459             tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
460             surfstr(i,j,1,n),sshf(i,j,1,n))
461        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
462        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
463      end do
464    end do
465  endif
466
467
468! Assign 10 m wind to model level at eta=1.0 to have one additional model
469! level at the ground
470! Specific humidity is taken the same as at one level above
471! Temperature is taken as 2 m temperature
472!**************************************************************************
473
474  do i=0,nxmin1
475    do j=0,nymin1
476      uuh(i,j,1)=u10(i,j,1,n)
477      vvh(i,j,1)=v10(i,j,1,n)
478      qvh(i,j,1,n)=qvh(i,j,2,n)
479      tth(i,j,1,n)=tt2(i,j,1,n)
480    end do
481  end do
482
483  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
484  if(iwmax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
485
486  return
487
488888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
489  write(*,*) ' #### ',wfname(indj),'                    #### '
490  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
491  stop 'Execution terminated'
492
493end subroutine readwind_ecmwf
494
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG