source: flexpart.git/src/readwind_ecmwf_mpi.f90 @ 92fab65

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

add SPDX-License-Identifier to all .f90 files

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