source: branches/petra/src/readwind.f90 @ 36

Last change on this file since 36 was 36, checked in by pesei, 9 years ago

Implement switch for incremental deposition, see ticket:113 and many small changes, see changelog.txt

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