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

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

move license from headers to a different file

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