source: flexpart.git/src/readwind_mpi.f90 @ 8a65cb0

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

Added code, makefile for dev branch

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