source: branches/petra/src/readwind_nests.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: 18.4 KB
Line 
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_nests(indj,n,uuhn,vvhn,wwhn)
23  !                           i   i  o    o    o
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads the wind fields for the nested model domains.       *
27  !     It is similar to subroutine readwind, which reads the mother domain.   *
28  !                                                                            *
29  !     Authors: A. Stohl, G. Wotawa                                           *
30  !                                                                            *
31  !     8 February 1999                                                        *
32  !                                                                            *
33  !     Last update: 17 October 2000, A. Stohl                                 *
34  !                                                                            *
35  !*****************************************************************************
36  !  Changes, Bernd C. Krueger, Feb. 2001:                                     *
37  !        Variables tthn and qvhn (on eta coordinates) in common block        *
38  !  CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with ECMWF grib_api    *
39  !  CHANGE: 03/12/2008, Harald Sodemann, update to f90 with ECMWF grib_api    *
40  ! 2/2015, PS: add missing paramter iret in call to grib subr
41  !
42  !*****************************************************************************
43
44  use grib_api
45  use par_mod
46  use com_mod
47
48  implicit none
49
50  !HSO  parameters for grib_api
51  integer :: ifile
52  integer :: iret
53  integer :: igrib
54  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
55  integer :: parId !!added by mc for making it consistent with new readwind.f90
56  integer :: gotGrid
57  !HSO  end
58
59  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
60  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
61  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
62  integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,l
63
64  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
65
66  ! dimension of isec2 at least (22+n), where n is the number of parallels or
67  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
68
69  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
70  ! coordinate parameters
71
72  integer :: isec1(56),isec2(22+nxmaxn+nymaxn)
73  real(kind=4) :: zsec4(jpunp)
74  real(kind=4) :: xaux,yaux
75  real(kind=8) :: xauxin,yauxin
76  real,parameter :: eps=1.e-4
77  real :: ewss(0:nxmaxn-1,0:nymaxn-1),nsss(0:nxmaxn-1,0:nymaxn-1)
78  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
79  real :: conversion_factor !added by mc to make it consistent with new gridchek.f90
80
81  logical :: hflswitch,strswitch
82
83  !HSO  grib api error messages
84  character(len=24) :: gribErrorMsg = 'Error reading grib file'
85  character(len=20) :: gribFunction = 'readwind_nests'
86
87  do l=1,numbnests
88    hflswitch=.false.
89    strswitch=.false.
90    levdiff2=nlev_ec-nwz+1
91    iumax=0
92    iwmax=0
93
94    ifile=0
95    igrib=0
96    iret=0
97
98  !
99  ! OPENING OF DATA FILE (GRIB CODE)
100  !
101
1025   call grib_open_file(ifile,path(numpath+2*(l-1)+1) &
103         (1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,indj)),'r')
104  if (iret.ne.GRIB_SUCCESS) then
105    goto 888   ! ERROR DETECTED
106  endif
107  !turn on support for multi fields messages */
108  !call grib_multi_support_on
109
110    gotGrid=0
111    ifield=0
11210   ifield=ifield+1
113  !
114  ! GET NEXT FIELDS
115  !
116  call grib_new_from_file(ifile,igrib,iret)
117  if (iret.eq.GRIB_END_OF_FILE)  then
118    goto 50    ! EOF DETECTED
119  elseif (iret.ne.GRIB_SUCCESS) then
120    goto 888   ! ERROR DETECTED
121  endif
122
123  !first see if we read GRIB1 or GRIB2
124  call grib_get_int(igrib,'editionNumber',gribVer,iret)
125  call grib_check(iret,gribFunction,gribErrorMsg)
126
127  if (gribVer.eq.1) then ! GRIB Edition 1
128
129  !print*,'GRiB Edition 1'
130  !read the grib2 identifiers
131  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
132  call grib_check(iret,gribFunction,gribErrorMsg)
133  call grib_get_int(igrib,'level',isec1(8),iret)
134  call grib_check(iret,gribFunction,gribErrorMsg)
135
136  !change code for etadot to code for omega
137  if (isec1(6).eq.77) then
138    isec1(6)=135
139  endif
140
141  conversion_factor=1.
142
143
144  else
145
146  !print*,'GRiB Edition 2'
147  !read the grib2 identifiers
148  call grib_get_int(igrib,'discipline',discipl,iret)
149  call grib_check(iret,gribFunction,gribErrorMsg)
150  call grib_get_int(igrib,'parameterCategory',parCat,iret)
151  call grib_check(iret,gribFunction,gribErrorMsg)
152  call grib_get_int(igrib,'parameterNumber',parNum,iret)
153  call grib_check(iret,gribFunction,gribErrorMsg)
154  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
155  call grib_check(iret,gribFunction,gribErrorMsg)
156  call grib_get_int(igrib,'level',valSurf,iret)
157  call grib_check(iret,gribFunction,gribErrorMsg)
158  call grib_get_int(igrib,'paramId',parId,iret) !added by mc to make it consisitent with new readwind.f90
159  call grib_check(iret,gribFunction,gribErrorMsg) !added by mc to make it consisitent with new readwind.f90
160
161  !print*,discipl,parCat,parNum,typSurf,valSurf
162
163  !convert to grib1 identifiers
164  isec1(6)=-1
165  isec1(7)=-1
166  isec1(8)=-1
167  isec1(8)=valSurf     ! level
168   conversion_factor=1.
169  if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
170    isec1(6)=130         ! indicatorOfParameter
171  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
172    isec1(6)=131         ! indicatorOfParameter
173  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
174    isec1(6)=132         ! indicatorOfParameter
175  elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
176    isec1(6)=133         ! indicatorOfParameter
177  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
178    isec1(6)=134         ! indicatorOfParameter
179  elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot !
180    isec1(6)=135         ! indicatorOfParameter
181  elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually eta dot !added by mc to make it consisitent with new readwind.f90
182    isec1(6)=135         ! indicatorOfParameter    !added by mc to make it consisitent with new readwind.f90
183  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
184    isec1(6)=151         ! indicatorOfParameter
185  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
186    isec1(6)=165         ! indicatorOfParameter
187  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
188    isec1(6)=166         ! indicatorOfParameter
189  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
190    isec1(6)=167         ! indicatorOfParameter
191  elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
192    isec1(6)=168         ! indicatorOfParameter
193  elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
194    isec1(6)=141         ! indicatorOfParameter
195    conversion_factor=1000. !added by mc to make it consisitent with new readwind.f90
196  elseif ((parCat.eq.6).and.(parNum.eq.1) .or. parId .eq. 164) then ! CC !added by mc to make it consisitent with new readwind.f90
197    isec1(6)=164         ! indicatorOfParameter
198  elseif ((parCat.eq.1).and.(parNum.eq.9) .or. parId .eq. 142) then ! LSP !added by mc to make it consisitent with new readwind.f90
199    isec1(6)=142         ! indicatorOfParameter
200  elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
201    isec1(6)=143         ! indicatorOfParameter
202    conversion_factor=1000. !added by mc to make it consisitent with new readwind.f90
203  elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
204    isec1(6)=146         ! indicatorOfParameter
205  elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
206    isec1(6)=176         ! indicatorOfParameter
207  elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS !added by mc to make it consisitent with new readwind.f90
208    isec1(6)=180         ! indicatorOfParameter
209  elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS !added by mc to make it consisitent with new readwind.f90
210    isec1(6)=181         ! indicatorOfParameter
211  elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
212    isec1(6)=129         ! indicatorOfParameter
213   elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO !added by mc to make it consisitent with new readwind.f90
214    isec1(6)=160         ! indicatorOfParameter
215  elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
216       (typSurf.eq.1)) then ! LSM
217    isec1(6)=172         ! indicatorOfParameter
218  else
219    print*,'***WARNING: undefined GRiB2 message found!',discipl, &
220         parCat,parNum,typSurf
221  endif
222  if(parId .ne. isec1(6) .and. parId .ne. 77) then !added by mc to make it consisitent with new readwind.f90
223    write(*,*) 'parId',parId, 'isec1(6)',isec1(6)  !
224!    stop
225  endif
226
227  endif
228
229  !HSO  get the size and data of the values array
230  if (isec1(6).ne.-1) then
231    call grib_get_real4_array(igrib,'values',zsec4,iret)
232    call grib_check(iret,gribFunction,gribErrorMsg)
233  endif
234
235  !HSO  get the required fields from section 2 in a gribex compatible manner
236  if(ifield.eq.1) then
237  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
238       isec2(2),iret)
239  call grib_check(iret,gribFunction,gribErrorMsg)
240  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
241       isec2(3),iret)
242  call grib_check(iret,gribFunction,gribErrorMsg)
243  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
244       isec2(12),iret)
245  call grib_check(iret,gribFunction,gribErrorMsg)
246  ! CHECK GRID SPECIFICATIONS
247  if(isec2(2).ne.nxn(l)) stop &
248  'READWIND: NX NOT CONSISTENT FOR A NESTING LEVEL'
249  if(isec2(3).ne.nyn(l)) stop &
250  'READWIND: NY NOT CONSISTENT FOR A NESTING LEVEL'
251  if(isec2(12)/2-1.ne.nlev_ec) stop 'READWIND: VERTICAL DISCRET&
252  IZATION NOT CONSISTENT FOR A NESTING LEVEL'
253  endif ! ifield
254
255  !HSO  get the second part of the grid dimensions only from GRiB1 messages
256 if (isec1(6) .eq. 167 .and. (gotGrid.eq.0)) then ! !added by mc to make it consisitent with new readwind.f90
257    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
258         xauxin,iret)
259    call grib_check(iret,gribFunction,gribErrorMsg)
260    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
261         yauxin,iret)
262    call grib_check(iret,gribFunction,gribErrorMsg)
263    if (xauxin.gt.180.) xauxin=xauxin-360.0
264    if (xauxin.lt.-180.) xauxin=xauxin+360.0
265
266    xaux=xauxin
267    yaux=yauxin
268    if (abs(xaux-xlon0n(l)).gt.eps) &
269    stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT FOR A NESTING LEVEL'
270    if (abs(yaux-ylat0n(l)).gt.eps) &
271    stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT FOR A NESTING LEVEL'
272    gotGrid=1
273  endif
274
275    do j=0,nyn(l)-1
276      do i=0,nxn(l)-1
277        k=isec1(8)
278        if(isec1(6).eq.130) tthn(i,j,nlev_ec-k+2,n,l)= &!! TEMPERATURE
279             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
280        if(isec1(6).eq.131) uuhn(i,j,nlev_ec-k+2,l)= &!! U VELOCITY
281             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
282        if(isec1(6).eq.132) vvhn(i,j,nlev_ec-k+2,l)= &!! V VELOCITY
283             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
284        if(isec1(6).eq.133) then                         !! SPEC. HUMIDITY
285          qvhn(i,j,nlev_ec-k+2,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
286          if (qvhn(i,j,nlev_ec-k+2,n,l) .lt. 0.) &
287               qvhn(i,j,nlev_ec-k+2,n,l) = 0.
288  !          this is necessary because the gridded data may contain
289  !          spurious negative values
290        endif
291        if(isec1(6).eq.134) psn(i,j,1,n,l)= &!! SURF. PRESS.
292             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
293        if(isec1(6).eq.135) wwhn(i,j,nlev_ec-k+1,l)= &!! W VELOCITY
294             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
295        if(isec1(6).eq.141) sdn(i,j,1,n,l)= &!! SNOW DEPTH
296             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/conversion_factor !added by mc to make it consisitent with new readwind.f90!
297        if(isec1(6).eq.151) msln(i,j,1,n,l)= &!! SEA LEVEL PRESS.
298             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
299        if(isec1(6).eq.164) tccn(i,j,1,n,l)= &!! CLOUD COVER
300             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
301        if(isec1(6).eq.165) u10n(i,j,1,n,l)= &!! 10 M U VELOCITY
302             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
303        if(isec1(6).eq.166) v10n(i,j,1,n,l)= &!! 10 M V VELOCITY
304             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
305        if(isec1(6).eq.167) tt2n(i,j,1,n,l)= &!! 2 M TEMPERATURE
306             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
307        if(isec1(6).eq.168) td2n(i,j,1,n,l)= &!! 2 M DEW POINT
308             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
309        if(isec1(6).eq.142) then                         !! LARGE SCALE PREC.
310          lsprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
311          if (lsprecn(i,j,1,n,l).lt.0.) lsprecn(i,j,1,n,l)=0.
312        endif
313        if(isec1(6).eq.143) then                         !! CONVECTIVE PREC.
314          convprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/conversion_factor !added by mc to make it consisitent with new readwind.f90
315          if (convprecn(i,j,1,n,l).lt.0.) convprecn(i,j,1,n,l)=0.
316        endif
317        if(isec1(6).eq.146) sshfn(i,j,1,n,l)= &!! SENS. HEAT FLUX
318             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
319        if((isec1(6).eq.146).and. &
320             (zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) hflswitch=.true.    ! Heat flux available
321        if(isec1(6).eq.176) then                         !! SOLAR RADIATION
322          ssrn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
323          if (ssrn(i,j,1,n,l).lt.0.) ssrn(i,j,1,n,l)=0.
324        endif
325        if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
326             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
327        if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
328             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
329        if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
330             (zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
331        if(isec1(6).eq.129) oron(i,j,l)= &!! ECMWF OROGRAPHY
332             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
333        if(isec1(6).eq.160) excessoron(i,j,l)= &!! STANDARD DEVIATION OF OROGRAPHY
334             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
335        if(isec1(6).eq.172) lsmn(i,j,l)= &!! ECMWF LAND SEA MASK
336             zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
337        if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
338        if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
339
340      end do
341    end do
342
343  call grib_release(igrib)
344  goto 10                      !! READ NEXT LEVEL OR PARAMETER
345  !
346  ! CLOSING OF INPUT DATA FILE
347  !
34850   call grib_close_file(ifile)
349
350  !error message if no fields found with correct first longitude in it
351  if (gotGrid.eq.0) then
352    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
353         'messages'
354    stop
355  endif
356
357  if(levdiff2.eq.0) then
358    iwmax=nlev_ec+1
359    do i=0,nxn(l)-1
360      do j=0,nyn(l)-1
361        wwhn(i,j,nlev_ec+1,l)=0.
362      end do
363    end do
364  endif
365
366  do i=0,nxn(l)-1
367    do j=0,nyn(l)-1
368      surfstrn(i,j,1,n,l)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
369    end do
370  end do
371
372  if ((.not.hflswitch).or.(.not.strswitch)) then
373    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
374         wfnamen(l,indj)
375
376  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
377  ! As ECMWF has increased the model resolution, such that now the first model
378  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
379  ! (3rd model level in FLEXPART) for the profile method
380  !***************************************************************************
381
382    do i=0,nxn(l)-1
383      do j=0,nyn(l)-1
384        plev1=akz(3)+bkz(3)*psn(i,j,1,n,l)
385        pmean=0.5*(psn(i,j,1,n,l)+plev1)
386        tv=tthn(i,j,3,n,l)*(1.+0.61*qvhn(i,j,3,n,l))
387        fu=-r_air*tv/ga/pmean
388        hlev1=fu*(plev1-psn(i,j,1,n,l))   ! HEIGTH OF FIRST MODEL LAYER
389        ff10m= sqrt(u10n(i,j,1,n,l)**2+v10n(i,j,1,n,l)**2)
390        fflev1=sqrt(uuhn(i,j,3,l)**2+vvhn(i,j,3,l)**2)
391        call pbl_profile(psn(i,j,1,n,l),td2n(i,j,1,n,l),hlev1, &
392             tt2n(i,j,1,n,l),tthn(i,j,3,n,l),ff10m,fflev1, &
393             surfstrn(i,j,1,n,l),sshfn(i,j,1,n,l))
394        if(sshfn(i,j,1,n,l).gt.200.) sshfn(i,j,1,n,l)=200.
395        if(sshfn(i,j,1,n,l).lt.-400.) sshfn(i,j,1,n,l)=-400.
396      end do
397    end do
398  endif
399
400
401  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
402  ! level at the ground
403  ! Specific humidity is taken the same as at one level above
404  ! Temperature is taken as 2 m temperature
405  !**************************************************************************
406
407    do i=0,nxn(l)-1
408      do j=0,nyn(l)-1
409        uuhn(i,j,1,l)=u10n(i,j,1,n,l)
410        vvhn(i,j,1,l)=v10n(i,j,1,n,l)
411        qvhn(i,j,1,n,l)=qvhn(i,j,2,n,l)
412        tthn(i,j,1,n,l)=tt2n(i,j,1,n,l)
413      end do
414    end do
415
416    if(iumax.ne.nuvz-1) stop &
417         'READWIND: NUVZ NOT CONSISTENT FOR A NESTING LEVEL'
418    if(iwmax.ne.nwz) stop &
419         'READWIND: NWZ NOT CONSISTENT FOR A NESTING LEVEL'
420
421  end do
422
423  return
424888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
425  write(*,*) ' #### ',wfnamen(l,indj),' FOR NESTING LEVEL  #### '
426  write(*,*) ' #### ',l,' IS NOT GRIB FORMAT !!!           #### '
427  stop 'Execution terminated'
428
429
430999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
431  write(*,*) ' #### ',wfnamen(l,indj),'                    #### '
432  write(*,*) ' #### CANNOT BE OPENED FOR NESTING LEVEL ',l,'####'
433
434end subroutine readwind_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG