source: flexpart.git/src/readwind_mpi.f90 @ d6a0977

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since d6a0977 was d6a0977, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Updates to Henrik's wet depo scheme

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