source: trunk/src/gridcheck_nests_emos.f90 @ 28

Last change on this file since 28 was 4, checked in by mlanger, 11 years ago
  • Property svn:executable set to *
File size: 10.1 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 gridcheck_nests
23
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine checks the grid specification for the nested model        *
27  !     domains. It is similar to subroutine gridcheck, which checks the       *
28  !     mother domain.                                                         *
29  !                                                                            *
30  !     Authors: A. Stohl, G. Wotawa                                           *
31  !                                                                            *
32  !     8 February 1999                                                        *
33  !                                                                            *
34  !*****************************************************************************
35
36  use par_mod
37  use com_mod
38
39  implicit none
40
41  integer :: i,j,k,l,ifn,ifield,iumax,iwmax,numskip,nlev_ecn
42  integer :: nuvzn,nwzn
43  real :: akmn(nwzmax),bkmn(nwzmax),akzn(nuvzmax),bkzn(nuvzmax)
44  real :: xaux1,xaux2,yaux1,yaux2
45
46  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
47
48  ! dimension of isec2 at least (22+n), where n is the number of parallels or
49  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
50
51  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
52  ! coordinate parameters
53
54  integer :: isec0(2),isec1(56),isec2(22+nxmaxn+nymaxn),isec3(2)
55  integer :: isec4(64),inbuff(jpack),ilen,ierr,iword,lunit
56  !integer iswap
57  real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp)
58  character(len=1) :: yoper = 'D'
59
60  xresoln(0)=1.       ! resolution enhancement for mother grid
61  yresoln(0)=1.       ! resolution enhancement for mother grid
62
63  ! Loop about all nesting levels
64  !******************************
65
66  do l=1,numbnests
67
68    iumax=0
69    iwmax=0
70
71    if(ideltas.gt.0) then
72      ifn=1
73    else
74      ifn=numbwf
75    endif
76  !
77  ! OPENING OF DATA FILE (GRIB CODE)
78  !
795   call pbopen(lunit,path(numpath+2*(l-1)+1) &
80         (1:length(numpath+2*(l-1)+1))//wfnamen(l,ifn),'r',ierr)
81    if(ierr.lt.0) goto 999
82
83    ifield=0
8410    ifield=ifield+1
85  !
86  ! GET NEXT FIELDS
87  !
88      call pbgrib(lunit,inbuff,jpack,ilen,ierr)
89      if(ierr.eq.-1) goto 30    ! EOF DETECTED
90      if(ierr.lt.-1) goto 999   ! ERROR DETECTED
91
92      ierr=1
93
94  ! Check whether we are on a little endian or on a big endian computer
95  !********************************************************************
96
97  !    if (inbuff(1).eq.1112101447) then         ! little endian, swap bytes
98  !      iswap=1+ilen/4
99  !      call swap32(inbuff,iswap)
100  !    else if (inbuff(1).ne.1196575042) then    ! big endian
101  !      stop 'subroutine gridcheck: corrupt GRIB data'
102  !    endif
103
104      call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, &
105           zsec4,jpunp,inbuff,jpack,iword,yoper,ierr)
106      if (ierr.ne.0) goto 999   ! ERROR DETECTED
107
108      if(ifield.eq.1) then
109        nxn(l)=isec2(2)
110        nyn(l)=isec2(3)
111        xaux1=real(isec2(5))/1000.
112        xaux2=real(isec2(8))/1000.
113        if(xaux1.gt.180) xaux1=xaux1-360.0
114        if(xaux2.gt.180) xaux2=xaux2-360.0
115        if(xaux1.lt.-180) xaux1=xaux1+360.0
116        if(xaux2.lt.-180) xaux2=xaux2+360.0
117        if (xaux2.lt.xaux1) xaux2=xaux2+360.
118        yaux1=real(isec2(7))/1000.
119        yaux2=real(isec2(4))/1000.
120        xlon0n(l)=xaux1
121        ylat0n(l)=yaux1
122        dxn(l)=(xaux2-xaux1)/real(nxn(l)-1)
123        dyn(l)=(yaux2-yaux1)/real(nyn(l)-1)
124        nlev_ecn=isec2(12)/2-1
125      endif
126
127
128      k=isec1(8)
129      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
130      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
131
132      if(isec1(6).eq.129) then
133      do j=0,nyn(l)-1
134        do i=0,nxn(l)-1
135          oron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
136        end do
137      end do
138      endif
139      if(isec1(6).eq.172) then
140      do j=0,nyn(l)-1
141        do i=0,nxn(l)-1
142          lsmn(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
143        end do
144      end do
145      endif
146      if(isec1(6).eq.160) then
147      do j=0,nyn(l)-1
148        do i=0,nxn(l)-1
149          excessoron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
150        end do
151      end do
152      endif
153
154      goto 10                 !! READ NEXT LEVEL OR PARAMETER
155  !
156  ! CLOSING OF INPUT DATA FILE
157  !
15830   call pbclose(lunit,ierr)     !! FINISHED READING / CLOSING GRIB FILE
159
160    nuvzn=iumax
161    nwzn=iwmax
162    if(nuvzn.eq.nlev_ec) nwzn=nlev_ecn+1
163
164    if (nxn(l).gt.nxmaxn) then
165  write(*,*) 'FLEXPART error: Too many grid points in x direction.'
166  write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)'
167  write(*,*) 'for nesting level ',l
168  write(*,*) 'Or change parameter settings in file par_mod.'
169  write(*,*) nxn(l),nxmaxn
170      stop
171    endif
172
173    if (nyn(l).gt.nymaxn) then
174  write(*,*) 'FLEXPART error: Too many grid points in y direction.'
175  write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)'
176  write(*,*) 'for nesting level ',l
177  write(*,*) 'Or change parameter settings in file par_mod.'
178  write(*,*) nyn(l),nymaxn
179      stop
180    endif
181
182    if ((nuvzn.gt.nuvzmax).or.(nwzn.gt.nwzmax)) then
183  write(*,*) 'FLEXPART error: Nested wind fields have too many'// &
184       'vertical levels.'
185  write(*,*) 'Problem was encountered for nesting level ',l
186      stop
187    endif
188
189
190  ! Output of grid info
191  !********************
192
193  write(*,'(a,i2)') 'Nested domain #: ',l
194  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Longitude range: ', &
195       xlon0n(l),' to ',xlon0n(l)+(nxn(l)-1)*dxn(l), &
196       '   Grid distance: ',dxn(l)
197  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Latitude range: ', &
198       ylat0n(l),' to ',ylat0n(l)+(nyn(l)-1)*dyn(l), &
199       '   Grid distance: ',dyn(l)
200  write(*,*)
201
202  ! Determine, how much the resolutions in the nests are enhanced as
203  ! compared to the mother grid
204  !*****************************************************************
205
206    xresoln(l)=dx/dxn(l)
207    yresoln(l)=dy/dyn(l)
208
209  ! Determine the mother grid coordinates of the corner points of the
210  ! nested grids
211  ! Convert first to geographical coordinates, then to grid coordinates
212  !********************************************************************
213
214    xaux1=xlon0n(l)
215    xaux2=xlon0n(l)+real(nxn(l)-1)*dxn(l)
216    yaux1=ylat0n(l)
217    yaux2=ylat0n(l)+real(nyn(l)-1)*dyn(l)
218
219    xln(l)=(xaux1-xlon0)/dx
220    xrn(l)=(xaux2-xlon0)/dx
221    yln(l)=(yaux1-ylat0)/dy
222    yrn(l)=(yaux2-ylat0)/dy
223
224
225    if ((xln(l).lt.0.).or.(yln(l).lt.0.).or. &
226         (xrn(l).gt.real(nxmin1)).or.(yrn(l).gt.real(nymin1))) then
227      write(*,*) 'Nested domain does not fit into mother domain'
228      write(*,*) 'For global mother domain fields, you can shift'
229      write(*,*) 'shift the mother domain into x-direction'
230      write(*,*) 'by setting nxshift (file par_mod) to a'
231      write(*,*) 'positive value. Execution is terminated.'
232      stop
233    endif
234
235
236  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
237  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
238
239    numskip=nlev_ecn-nuvzn ! number of ecmwf model layers not used by FLEXPART
240    do i=1,nwzn
241      j=numskip+i
242      k=nlev_ecn+1+numskip+i
243      akmn(nwzn-i+1)=zsec2(10+j)
244      bkmn(nwzn-i+1)=zsec2(10+k)
245    end do
246
247  !
248  ! CALCULATION OF AKZ, BKZ
249  ! AKZ,BKZ: model discretization parameters at the center of each model
250  !     layer
251  !
252  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
253  ! i.e. ground level
254  !*****************************************************************************
255
256    akzn(1)=0.
257    bkzn(1)=1.0
258    do i=1,nuvzn
259      akzn(i+1)=0.5*(akmn(i+1)+akmn(i))
260      bkzn(i+1)=0.5*(bkmn(i+1)+bkmn(i))
261    end do
262    nuvzn=nuvzn+1
263
264  ! Check, whether the heights of the model levels of the nested
265  ! wind fields are consistent with those of the mother domain.
266  ! If not, terminate model run.
267  !*************************************************************
268
269    do i=1,nuvz
270      if ((akzn(i).ne.akz(i)).or.(bkzn(i).ne.bkz(i))) then
271  write(*,*) 'FLEXPART error: The wind fields of nesting level',l
272  write(*,*) 'are not consistent with the mother domain:'
273  write(*,*) 'Differences in vertical levels detected.'
274        stop
275      endif
276    end do
277
278    do i=1,nwz
279      if ((akmn(i).ne.akm(i)).or.(bkmn(i).ne.bkm(i))) then
280  write(*,*) 'FLEXPART error: The wind fields of nesting level',l
281  write(*,*) 'are not consistent with the mother domain:'
282  write(*,*) 'Differences in vertical levels detected.'
283        stop
284      endif
285    end do
286
287  end do
288
289  return
290
291999   write(*,*)
292  write(*,*) ' ###########################################'// &
293       '###### '
294  write(*,*) '                FLEXPART SUBROUTINE GRIDCHECK:'
295  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfnamen(l,ifn)
296  write(*,*) ' FOR NESTING LEVEL ',k
297  write(*,*) ' ###########################################'// &
298       '###### '
299  stop
300
301end subroutine gridcheck_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG