source: trunk/src/gridcheck_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: 14.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 gridcheck
23
24  !**********************************************************************
25  !                                                                     *
26  !             FLEXPART MODEL SUBROUTINE GRIDCHECK                     *
27  !                                                                     *
28  !**********************************************************************
29  !                                                                     *
30  !             AUTHOR:      G. WOTAWA                                  *
31  !             DATE:        1997-08-06                                 *
32  !             LAST UPDATE: 1997-10-10                                 *
33  !                                                                     *
34  !             Update:      1999-02-08, global fields allowed, A. Stohl*
35  !                                                                     *
36  !**********************************************************************
37  !                                                                     *
38  ! DESCRIPTION:                                                        *
39  !                                                                     *
40  ! THIS SUBROUTINE DETERMINES THE GRID SPECIFICATIONS (LOWER LEFT      *
41  ! LONGITUDE, LOWER LEFT LATITUDE, NUMBER OF GRID POINTS, GRID DIST-   *
42  ! ANCE AND VERTICAL DISCRETIZATION OF THE ECMWF MODEL) FROM THE       *
43  ! GRIB HEADER OF THE FIRST INPUT FILE. THE CONSISTANCY (NO CHANGES    *
44  ! WITHIN ONE FLEXPART RUN) IS CHECKED IN THE ROUTINE "READWIND" AT    *
45  ! ANY CALL.                                                           *
46  !                                                                     *
47  ! XLON0                geographical longitude of lower left gridpoint *
48  ! YLAT0                geographical latitude of lower left gridpoint  *
49  ! NX                   number of grid points x-direction              *
50  ! NY                   number of grid points y-direction              *
51  ! DX                   grid distance x-direction                      *
52  ! DY                   grid distance y-direction                      *
53  ! NUVZ                 number of grid points for horizontal wind      *
54  !                      components in z direction                      *
55  ! NWZ                  number of grid points for vertical wind        *
56  !                      component in z direction                       *
57  ! sizesouth, sizenorth give the map scale (i.e. number of virtual grid*
58  !                      points of the polar stereographic grid):       *
59  !                      used to check the CFL criterion                *
60  ! UVHEIGHT(1)-         heights of gridpoints where u and v are        *
61  ! UVHEIGHT(NUVZ)       given                                          *
62  ! WHEIGHT(1)-          heights of gridpoints where w is given         *
63  ! WHEIGHT(NWZ)                                                        *
64  !                                                                     *
65  !**********************************************************************
66 
67  use par_mod
68  use com_mod
69  use conv_mod
70
71  implicit none
72
73  integer :: ix,jy,i,ifn,ifield,j,k,iumax,iwmax,numskip
74  real :: xaux1,xaux2,yaux1,yaux2,sizesouth,sizenorth,xauxa,pint
75
76  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
77
78  ! dimension of isec2 at least (22+n), where n is the number of parallels or
79  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
80
81  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
82  ! coordinate parameters
83
84  integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2)
85  integer :: isec4(64),inbuff(jpack),ilen,ierr,iword,lunit
86  !integer iswap
87  real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp)
88  character(len=1) :: opt, yoper = 'D'
89
90
91  iumax=0
92  iwmax=0
93
94  if(ideltas.gt.0) then
95    ifn=1
96  else
97    ifn=numbwf
98  endif
99  !
100  ! OPENING OF DATA FILE (GRIB CODE)
101  !
1025   call pbopen(lunit,path(3)(1:length(3))//wfname(ifn),'r',ierr)
103  if(ierr.lt.0) goto 999
104
105  ifield=0
10610   ifield=ifield+1
107  !
108  ! GET NEXT FIELDS
109  !
110  call pbgrib(lunit,inbuff,jpack,ilen,ierr)
111  if(ierr.eq.-1) goto 30    ! EOF DETECTED
112  if(ierr.lt.-1) goto 999   ! ERROR DETECTED
113
114  ierr=1
115
116  ! Check whether we are on a little endian or on a big endian computer
117  !********************************************************************
118
119  !if (inbuff(1).eq.1112101447) then         ! little endian, swap bytes
120  !  iswap=1+ilen/4
121  !  call swap32(inbuff,iswap)
122  !else if (inbuff(1).ne.1196575042) then    ! big endian
123  !  stop 'subroutine gridcheck: corrupt GRIB data'
124  !endif
125
126  call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, &
127       zsec4,jpunp,inbuff,jpack,iword,yoper,ierr)
128  if (ierr.ne.0) goto 999   ! ERROR DETECTED
129
130  if(ifield.eq.1) then
131    nxfield=isec2(2)
132    ny=isec2(3)
133    xaux1=real(isec2(5))/1000.
134    xaux2=real(isec2(8))/1000.
135    if (xaux1.gt.180) xaux1=xaux1-360.0
136    if (xaux2.gt.180) xaux2=xaux2-360.0
137    if (xaux1.lt.-180) xaux1=xaux1+360.0
138    if (xaux2.lt.-180) xaux2=xaux2+360.0
139    if (xaux2.lt.xaux1) xaux2=xaux2+360.
140    yaux1=real(isec2(7))/1000.
141    yaux2=real(isec2(4))/1000.
142    xlon0=xaux1
143    ylat0=yaux1
144    dx=(xaux2-xaux1)/real(nxfield-1)
145    dy=(yaux2-yaux1)/real(ny-1)
146    dxconst=180./(dx*r_earth*pi)
147    dyconst=180./(dy*r_earth*pi)
148    nlev_ec=isec2(12)/2-1
149
150
151  ! Check whether fields are global
152  ! If they contain the poles, specify polar stereographic map
153  ! projections using the stlmbr- and stcm2p-calls
154  !***********************************************************
155
156    xauxa=abs(xaux2+dx-360.-xaux1)
157    if (xauxa.lt.0.001) then
158      nx=nxfield+1                 ! field is cyclic
159      xglobal=.true.
160      if (abs(nxshift).ge.nx) &
161           stop 'nxshift in file par_mod is too large'
162      xlon0=xlon0+real(nxshift)*dx
163    else
164      nx=nxfield
165      xglobal=.false.
166      if (nxshift.ne.0) &
167           stop 'nxshift (par_mod) must be zero for non-global domain'
168    endif
169    nxmin1=nx-1
170    nymin1=ny-1
171    if (xlon0.gt.180.) xlon0=xlon0-360.
172    xauxa=abs(yaux1+90.)
173    if (xglobal.and.xauxa.lt.0.001) then
174      sglobal=.true.               ! field contains south pole
175  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
176  ! map scale
177      sizesouth=6.*(switchsouth+90.)/dy
178      call stlmbr(southpolemap,-90.,0.)
179      call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, &
180           sizesouth,switchsouth,180.)
181      switchsouthg=(switchsouth-ylat0)/dy
182    else
183      sglobal=.false.
184      switchsouthg=999999.
185    endif
186    xauxa=abs(yaux2-90.)
187    if (xglobal.and.xauxa.lt.0.001) then
188      nglobal=.true.               ! field contains north pole
189  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
190  ! map scale
191      sizenorth=6.*(90.-switchnorth)/dy
192      call stlmbr(northpolemap,90.,0.)
193      call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, &
194           sizenorth,switchnorth,180.)
195      switchnorthg=(switchnorth-ylat0)/dy
196    else
197      nglobal=.false.
198      switchnorthg=999999.
199    endif
200  endif
201
202  if (nxshift.lt.0) stop 'nxshift (par_mod) must not be negative'
203  if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large'
204
205  k=isec1(8)
206  if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
207  if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
208
209  if(isec1(6).eq.129) then
210    do jy=0,ny-1
211      do ix=0,nxfield-1
212  !       write (*,*) 'ich stop!',nxfield,ny,jy,ix,ga
213  !       stop
214        oro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)/ga
215      end do
216    end do
217  endif
218  if(isec1(6).eq.172) then
219    do jy=0,ny-1
220      do ix=0,nxfield-1
221        lsm(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
222      end do
223    end do
224  endif
225  if(isec1(6).eq.160) then
226    do jy=0,ny-1
227      do ix=0,nxfield-1
228        excessoro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
229      end do
230    end do
231  endif
232
233  goto 10                      !! READ NEXT LEVEL OR PARAMETER
234  !
235  ! CLOSING OF INPUT DATA FILE
236  !
23730   call pbclose(lunit,ierr)     !! FINISHED READING / CLOSING GRIB FILE
238
239  nuvz=iumax
240  nwz =iwmax
241  if(nuvz.eq.nlev_ec) nwz=nlev_ec+1
242
243  if (nx.gt.nxmax) then
244   write(*,*) 'FLEXPART error: Too many grid points in x direction.'
245    write(*,*) 'Reduce resolution of wind fields.'
246    write(*,*) 'Or change parameter settings in file par_mod.'
247    write(*,*) nx,nxmax
248    stop
249  endif
250
251  if (ny.gt.nymax) then
252   write(*,*) 'FLEXPART error: Too many grid points in y direction.'
253    write(*,*) 'Reduce resolution of wind fields.'
254    write(*,*) 'Or change parameter settings in file par_mod.'
255    write(*,*) ny,nymax
256    stop
257  endif
258
259  if (nuvz+1.gt.nuvzmax) then
260    write(*,*) 'FLEXPART error: Too many u,v grid points in z '// &
261         'direction.'
262    write(*,*) 'Reduce resolution of wind fields.'
263    write(*,*) 'Or change parameter settings in file par_mod.'
264    write(*,*) nuvz+1,nuvzmax
265    stop
266  endif
267
268  if (nwz.gt.nwzmax) then
269    write(*,*) 'FLEXPART error: Too many w grid points in z '// &
270         'direction.'
271    write(*,*) 'Reduce resolution of wind fields.'
272    write(*,*) 'Or change parameter settings in file par_mod.'
273    write(*,*) nwz,nwzmax
274    stop
275  endif
276
277  ! If desired, shift all grids by nxshift grid cells
278  !**************************************************
279
280  if (xglobal) then
281    call shift_field_0(oro,nxfield,ny)
282    call shift_field_0(lsm,nxfield,ny)
283    call shift_field_0(excessoro,nxfield,ny)
284  endif
285
286  ! Output of grid info
287  !********************
288
289  write(*,*)
290  write(*,*)
291  write(*,'(a,2i7)') '# of vertical levels in ECMWF data: ', &
292       nuvz+1,nwz
293  write(*,*)
294  write(*,'(a)') 'Mother domain:'
295  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Longitude range: ', &
296       xlon0,' to ',xlon0+(nx-1)*dx,'   Grid distance: ',dx
297  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Latitude range: ', &
298       ylat0,' to ',ylat0+(ny-1)*dy,'   Grid distance: ',dy
299  write(*,*)
300
301
302  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
303  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
304
305  numskip=nlev_ec-nuvz  ! number of ecmwf model layers not used
306                        ! by trajectory model
307  do i=1,nwz
308    j=numskip+i
309    k=nlev_ec+1+numskip+i
310    akm(nwz-i+1)=zsec2(10+j)
311    bkm(nwz-i+1)=zsec2(10+k)
312  end do
313
314  !
315  ! CALCULATION OF AKZ, BKZ
316  ! AKZ,BKZ: model discretization parameters at the center of each model
317  !     layer
318  !
319  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
320  ! i.e. ground level
321  !*****************************************************************************
322
323  akz(1)=0.
324  bkz(1)=1.0
325  do i=1,nuvz
326     akz(i+1)=0.5*(akm(i+1)+akm(i))
327     bkz(i+1)=0.5*(bkm(i+1)+bkm(i))
328  end do
329  nuvz=nuvz+1
330
331  ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled
332  ! upon the transformation to z levels. In order to save computer memory, this is
333  ! not done anymore in the standard version. However, this option can still be
334  ! switched on by replacing the following lines with those below, that are
335  ! currently commented out. For this, similar changes are necessary in
336  ! verttransform.f and verttranform_nests.f
337  !*****************************************************************************
338
339  nz=nuvz
340  if (nz.gt.nzmax) stop 'nzmax too small'
341  do i=1,nuvz
342    aknew(i)=akz(i)
343    bknew(i)=bkz(i)
344  end do
345
346  ! Switch on following lines to use doubled vertical resolution
347  !*************************************************************
348  !nz=nuvz+nwz-1
349  !if (nz.gt.nzmax) stop 'nzmax too small'
350  !do 100 i=1,nwz
351  !  aknew(2*(i-1)+1)=akm(i)
352  !00     bknew(2*(i-1)+1)=bkm(i)
353  !do 110 i=2,nuvz
354  !  aknew(2*(i-1))=akz(i)
355  !10     bknew(2*(i-1))=bkz(i)
356  ! End doubled vertical resolution
357
358
359  ! Determine the uppermost level for which the convection scheme shall be applied
360  ! by assuming that there is no convection above 50 hPa (for standard SLP)
361  !*****************************************************************************
362
363  do i=1,nuvz-2
364    pint=akz(i)+bkz(i)*101325.
365    if (pint.lt.5000.) goto 96
366  end do
36796   nconvlev=i
368  if (nconvlev.gt.nconvlevmax-1) then
369    nconvlev=nconvlevmax-1
370    write(*,*) 'Attention, convection only calculated up to ', &
371         akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa'
372  endif
373
374  return
375
376999   write(*,*)
377  write(*,*) ' ###########################################'// &
378       '###### '
379  write(*,*) '       TRAJECTORY MODEL SUBROUTINE GRIDCHECK:'
380  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn)
381  write(*,*) ' ###########################################'// &
382       '###### '
383  write(*,*)
384  write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND   !!!'
385  write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE...     !!!'
386  write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!'
387  write(*,'(a)') '!!! THE "X" KEY...                   !!!'
388  write(*,*)
389  read(*,'(a)') opt
390  if(opt.eq.'X') then
391    stop
392  else
393    goto 5
394  endif
395
396end subroutine gridcheck
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG