source: branches/jerome/src/gridcheck_gfs.f90 @ 15

Last change on this file since 15 was 15, checked in by jebri, 11 years ago

fixes for flexpart91

File size: 18.9 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  !             CHANGE: 17/11/2005, Caroline Forster, GFS data          *
36  !             CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with *
37  !                                 ECMWF grib_api                      *
38  !             CHANGE: 03/12/2008, Harald Sodemann, update to f90 with *
39  !                                 ECMWF grib_api                      *
40  !             CHANGE:             Jerome Brioude, sort akm            *
41  !**********************************************************************
42  !                                                                     *
43  ! DESCRIPTION:                                                        *
44  !                                                                     *
45  ! THIS SUBROUTINE DETERMINES THE GRID SPECIFICATIONS (LOWER LEFT      *
46  ! LONGITUDE, LOWER LEFT LATITUDE, NUMBER OF GRID POINTS, GRID DIST-   *
47  ! ANCE AND VERTICAL DISCRETIZATION OF THE ECMWF MODEL) FROM THE       *
48  ! GRIB HEADER OF THE FIRST INPUT FILE. THE CONSISTANCY (NO CHANGES    *
49  ! WITHIN ONE FLEXPART RUN) IS CHECKED IN THE ROUTINE "READWIND" AT    *
50  ! ANY CALL.                                                           *
51  !                                                                     *
52  ! XLON0                geographical longitude of lower left gridpoint *
53  ! YLAT0                geographical latitude of lower left gridpoint  *
54  ! NX                   number of grid points x-direction              *
55  ! NY                   number of grid points y-direction              *
56  ! DX                   grid distance x-direction                      *
57  ! DY                   grid distance y-direction                      *
58  ! NUVZ                 number of grid points for horizontal wind      *
59  !                      components in z direction                      *
60  ! NWZ                  number of grid points for vertical wind        *
61  !                      component in z direction                       *
62  ! sizesouth, sizenorth give the map scale (i.e. number of virtual grid*
63  !                      points of the polar stereographic grid):       *
64  !                      used to check the CFL criterion                *
65  ! UVHEIGHT(1)-         heights of gridpoints where u and v are        *
66  ! UVHEIGHT(NUVZ)       given                                          *
67  ! WHEIGHT(1)-          heights of gridpoints where w is given         *
68  ! WHEIGHT(NWZ)                                                        *
69  !                                                                     *
70  !**********************************************************************
71
72  use grib_api
73  use par_mod
74  use com_mod
75  use conv_mod
76  use cmapf_mod, only: stlmbr,stcm2p
77
78  implicit none
79
80  !HSO  parameters for grib_api
81  integer :: ifile
82  integer :: iret
83  integer :: igrib,inc
84  real(kind=4) :: xaux1,xaux2,yaux1,yaux2
85  real(kind=8) :: xaux1in,xaux2in,yaux1in,yaux2in
86  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
87  !HSO  end
88  integer :: ix,jy,i,ifn,ifield,j,k,iumax,iwmax,numskip
89  real :: sizesouth,sizenorth,xauxa,pint
90  real :: akm_usort(nwzmax),akm2(nwzmax)
91  real,parameter :: eps=0.0001
92
93  ! NCEP GFS
94  real :: pres(nwzmax), help
95
96  integer :: i179,i180,i181
97
98  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
99
100  integer :: isec1(8),isec2(3)
101  real(kind=4) :: zsec4(jpunp)
102  character(len=1) :: opt
103
104  !HSO  grib api error messages
105  character(len=24) :: gribErrorMsg = 'Error reading grib file'
106  character(len=20) :: gribFunction = 'gridcheckwind_gfs'
107  !
108  if (numbnests.ge.1) then
109  write(*,*) ' ###########################################'
110  write(*,*) ' FLEXPART ERROR SUBROUTINE GRIDCHECK:'
111  write(*,*) ' NO NESTED WINDFIELDAS ALLOWED FOR GFS!      '
112  write(*,*) ' ###########################################'
113  stop
114  endif
115
116  iumax=0
117  iwmax=0
118
119  if(ideltas.gt.0) then
120    ifn=1
121  else
122    ifn=numbwf
123  endif
124  !
125  ! OPENING OF DATA FILE (GRIB CODE)
126  !
1275   call grib_open_file(ifile,path(3)(1:length(3)) &
128         //trim(wfname(ifn)),'r',iret)
129  if (iret.ne.GRIB_SUCCESS) then
130    goto 999   ! ERROR DETECTED
131  endif
132  !turn on support for multi fields messages
133  call grib_multi_support_on
134
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 30    ! EOF DETECTED
143  elseif (iret.ne.GRIB_SUCCESS) then
144    goto 999   ! 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  !read the grib1 identifiers
154  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
155  call grib_check(iret,gribFunction,gribErrorMsg)
156  call grib_get_int(igrib,'indicatorOfTypeOfLevel',isec1(7),iret)
157  call grib_check(iret,gribFunction,gribErrorMsg)
158  call grib_get_int(igrib,'level',isec1(8),iret)
159  call grib_check(iret,gribFunction,gribErrorMsg)
160
161  !get the size and data of the values array
162  call grib_get_real4_array(igrib,'values',zsec4,iret)
163  call grib_check(iret,gribFunction,gribErrorMsg)
164
165  else ! GRIB Edition 2
166
167  !read the grib2 identifiers
168  call grib_get_int(igrib,'discipline',discipl,iret)
169  call grib_check(iret,gribFunction,gribErrorMsg)
170  call grib_get_int(igrib,'parameterCategory',parCat,iret)
171  call grib_check(iret,gribFunction,gribErrorMsg)
172  call grib_get_int(igrib,'parameterNumber',parNum,iret)
173  call grib_check(iret,gribFunction,gribErrorMsg)
174  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
175  call grib_check(iret,gribFunction,gribErrorMsg)
176  call grib_get_int(igrib,'scaledValueOfFirstFixedSurface', &
177       valSurf,iret)
178  call grib_check(iret,gribFunction,gribErrorMsg)
179
180  !convert to grib1 identifiers
181  isec1(6)=-1
182  isec1(7)=-1
183  isec1(8)=-1
184  if ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.100)) then ! U
185    isec1(6)=33          ! indicatorOfParameter
186    isec1(7)=100         ! indicatorOfTypeOfLevel
187    isec1(8)=valSurf/100 ! level, convert to hPa
188  elseif ((parCat.eq.3).and.(parNum.eq.5).and.(typSurf.eq.1)) then ! TOPO
189    isec1(6)=7           ! indicatorOfParameter
190    isec1(7)=1           ! indicatorOfTypeOfLevel
191    isec1(8)=0
192  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.1) &
193       .and.(discipl.eq.2)) then ! LSM
194    isec1(6)=81          ! indicatorOfParameter
195    isec1(7)=1           ! indicatorOfTypeOfLevel
196    isec1(8)=0
197  endif
198
199  if (isec1(6).ne.-1) then
200  !  get the size and data of the values array
201    call grib_get_real4_array(igrib,'values',zsec4,iret)
202    call grib_check(iret,gribFunction,gribErrorMsg)
203  endif
204
205  endif ! gribVer
206
207  if(ifield.eq.1) then
208
209  !get the required fields from section 2
210  !store compatible to gribex input
211  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
212       isec2(2),iret)
213  call grib_check(iret,gribFunction,gribErrorMsg)
214  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
215       isec2(3),iret)
216  call grib_check(iret,gribFunction,gribErrorMsg)
217  call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
218       xaux1in,iret)
219  call grib_check(iret,gribFunction,gribErrorMsg)
220  call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', &
221       xaux2in,iret)
222  call grib_check(iret,gribFunction,gribErrorMsg)
223  call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
224       yaux1in,iret)
225  call grib_check(iret,gribFunction,gribErrorMsg)
226  call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', &
227       yaux2in,iret)
228  call grib_check(iret,gribFunction,gribErrorMsg)
229
230  xaux1=xaux1in
231  xaux2=xaux2in
232  yaux1=yaux1in
233  yaux2=yaux2in
234
235  nxfield=isec2(2)
236  ny=isec2(3)
237  if((abs(xaux1).lt.eps).and.(xaux2.ge.359)) then ! NCEP DATA FROM 0 TO
238    xaux1=-179.0                             ! 359 DEG EAST ->
239    xaux2=-179.0+360.-360./real(nxfield)    ! TRANSFORMED TO -179
240  endif                                      ! TO 180 DEG EAST
241  if (xaux1.gt.180) xaux1=xaux1-360.0
242  if (xaux2.gt.180) xaux2=xaux2-360.0
243  if (xaux1.lt.-180) xaux1=xaux1+360.0
244  if (xaux2.lt.-180) xaux2=xaux2+360.0
245  if (xaux2.lt.xaux1) xaux2=xaux2+360.
246  xlon0=xaux1
247  ylat0=yaux1
248  dx=(xaux2-xaux1)/real(nxfield-1)
249  dy=(yaux2-yaux1)/real(ny-1)
250  dxconst=180./(dx*r_earth*pi)
251  dyconst=180./(dy*r_earth*pi)
252  !HSO end edits
253
254
255  ! Check whether fields are global
256  ! If they contain the poles, specify polar stereographic map
257  ! projections using the stlmbr- and stcm2p-calls
258  !***********************************************************
259
260    xauxa=abs(xaux2+dx-360.-xaux1)
261    if (xauxa.lt.0.001) then
262      nx=nxfield+1                 ! field is cyclic
263      xglobal=.true.
264      if (abs(nxshift).ge.nx) &
265           stop 'nxshift in file par_mod is too large'
266      xlon0=xlon0+real(nxshift)*dx
267    else
268      nx=nxfield
269      xglobal=.false.
270      if (nxshift.ne.0) &
271           stop 'nxshift (par_mod) must be zero for non-global domain'
272    endif
273    nxmin1=nx-1
274    nymin1=ny-1
275    if (xlon0.gt.180.) xlon0=xlon0-360.
276    xauxa=abs(yaux1+90.)
277    if (xglobal.and.xauxa.lt.0.001) then
278      sglobal=.true.               ! field contains south pole
279  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
280  ! map scale
281      sizesouth=6.*(switchsouth+90.)/dy
282      call stlmbr(southpolemap,-90.,0.)
283      call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, &
284           sizesouth,switchsouth,180.)
285      switchsouthg=(switchsouth-ylat0)/dy
286    else
287      sglobal=.false.
288      switchsouthg=999999.
289    endif
290    xauxa=abs(yaux2-90.)
291    if (xglobal.and.xauxa.lt.0.001) then
292      nglobal=.true.               ! field contains north pole
293  ! Enhance the map scale by factor 3 (*2=6) compared to north-south
294  ! map scale
295      sizenorth=6.*(90.-switchnorth)/dy
296      call stlmbr(northpolemap,90.,0.)
297      call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, &
298           sizenorth,switchnorth,180.)
299      switchnorthg=(switchnorth-ylat0)/dy
300    else
301      nglobal=.false.
302      switchnorthg=999999.
303    endif
304  endif ! ifield.eq.1
305
306  if (nxshift.lt.0) stop 'nxshift (par_mod) must not be negative'
307  if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large'
308
309  ! NCEP ISOBARIC LEVELS
310  !*********************
311
312  if((isec1(6).eq.33).and.(isec1(7).eq.100)) then ! check for U wind
313    iumax=iumax+1
314    pres(iumax)=real(isec1(8))*100.0
315  endif
316
317
318  i179=nint(179./dx)
319  if (dx.lt.0.7) then
320    i180=nint(180./dx)+1    ! 0.5 deg data
321  else
322    i180=nint(179./dx)+1    ! 1 deg data
323  endif
324  i181=i180+1
325
326
327  ! NCEP TERRAIN
328  !*************
329
330  if((isec1(6).eq.007).and.(isec1(7).eq.001)) then
331    do jy=0,ny-1
332      do ix=0,nxfield-1
333        help=zsec4(nxfield*(ny-jy-1)+ix+1)
334        if(ix.le.i180) then
335          oro(i179+ix,jy)=help
336          excessoro(i179+ix,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
337        else
338          oro(ix-i181,jy)=help
339          excessoro(ix-i181,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED
340        endif
341      end do
342    end do
343  endif
344
345  ! NCEP LAND SEA MASK
346  !*******************
347
348  if((isec1(6).eq.081).and.(isec1(7).eq.001)) then
349    do jy=0,ny-1
350      do ix=0,nxfield-1
351        help=zsec4(nxfield*(ny-jy-1)+ix+1)
352        if(ix.le.i180) then
353          lsm(i179+ix,jy)=help
354        else
355          lsm(ix-i181,jy)=help
356        endif
357      end do
358    end do
359  endif
360
361  call grib_release(igrib)
362
363  goto 10                      !! READ NEXT LEVEL OR PARAMETER
364  !
365  ! CLOSING OF INPUT DATA FILE
366  !
367
368  ! HSO
36930   continue
370  call grib_close_file(ifile)
371  ! HSO end edits
372
373  nuvz=iumax
374  nwz =iumax
375  nlev_ec=iumax
376
377  if (nx.gt.nxmax) then
378   write(*,*) 'FLEXPART error: Too many grid points in x direction.'
379    write(*,*) 'Reduce resolution of wind fields.'
380    write(*,*) 'Or change parameter settings in file par_mod.'
381    write(*,*) nx,nxmax
382    stop
383  endif
384
385  if (ny.gt.nymax) then
386   write(*,*) 'FLEXPART error: Too many grid points in y direction.'
387    write(*,*) 'Reduce resolution of wind fields.'
388    write(*,*) 'Or change parameter settings in file par_mod.'
389    write(*,*) ny,nymax
390    stop
391  endif
392
393  if (nuvz.gt.nuvzmax) then
394    write(*,*) 'FLEXPART error: Too many u,v grid points in z '// &
395         'direction.'
396    write(*,*) 'Reduce resolution of wind fields.'
397    write(*,*) 'Or change parameter settings in file par_mod.'
398    write(*,*) nuvz,nuvzmax
399    stop
400  endif
401
402  if (nwz.gt.nwzmax) then
403    write(*,*) 'FLEXPART error: Too many w grid points in z '// &
404         'direction.'
405    write(*,*) 'Reduce resolution of wind fields.'
406    write(*,*) 'Or change parameter settings in file par_mod.'
407    write(*,*) nwz,nwzmax
408    stop
409  endif
410
411  ! If desired, shift all grids by nxshift grid cells
412  !**************************************************
413
414  if (xglobal) then
415    call shift_field_0(oro,nxfield,ny)
416    call shift_field_0(lsm,nxfield,ny)
417    call shift_field_0(excessoro,nxfield,ny)
418  endif
419
420  ! Output of grid info
421  !********************
422
423  write(*,*)
424  write(*,*)
425  write(*,'(a,2i7)') '# of vertical levels in NCEP data: ', &
426       nuvz,nwz
427  write(*,*)
428  write(*,'(a)') 'Mother domain:'
429  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Longitude range: ', &
430       xlon0,' to ',xlon0+(nx-1)*dx,'   Grid distance: ',dx
431  write(*,'(a,f10.2,a1,f10.2,a,f10.2)') '  Latitude range: ', &
432       ylat0,' to ',ylat0+(ny-1)*dy,'   Grid distance: ',dy
433  write(*,*)
434
435
436  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
437  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
438
439  numskip=nlev_ec-nuvz  ! number of ecmwf model layers not used
440                        ! by trajectory model
441  do i=1,nwz
442    j=numskip+i
443    k=nlev_ec+1+numskip+i
444    akm_usort(nwz-i+1)=pres(nwz-i+1)
445    bkm(nwz-i+1)=0.0
446  end do
447
448  !******************************
449  ! change Sabine Eckhardt: akm should always be in descending order ... readwind adapted!
450  !******************************
451      do i=1,nwz
452         if (akm_usort(1).gt.akm_usort(2)) then
453            akm(i)=akm_usort(i)
454         else
455            akm(i)=akm_usort(nwz-i+1)
456         endif
457      end do
458
459!***************************             
460!JB sort akm which can be random
461      do i=1,nwzmax
462      akm2(i)=0.
463      enddo
464      do i=1,nwzmax
465      if (akm(i).gt.akm2(1)) akm2(1)=akm(i)
466      enddo
467      inc=2
468
46944    do i=1,nwzmax
470      if(akm(i).lt.akm2(inc-1).and.akm(i).gt.akm2(inc))akm2(inc)=akm(i)
471      enddo
472      inc=inc+1
473      if (inc.lt.nwzmax) goto 44
47445    continue
475      akm=akm2
476! end JB
477
478  !
479  ! CALCULATION OF AKZ, BKZ
480  ! AKZ,BKZ: model discretization parameters at the center of each model
481  !     layer
482  !
483  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
484  ! i.e. ground level
485  !*****************************************************************************
486
487  do i=1,nuvz
488     akz(i)=akm(i)
489     bkz(i)=bkm(i)
490  end do
491
492  ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled
493  ! upon the transformation to z levels. In order to save computer memory, this is
494  ! not done anymore in the standard version. However, this option can still be
495  ! switched on by replacing the following lines with those below, that are
496  ! currently commented out. For this, similar changes are necessary in
497  ! verttransform.f and verttranform_nests.f
498  !*****************************************************************************
499
500  nz=nuvz
501  if (nz.gt.nzmax) stop 'nzmax too small'
502  do i=1,nuvz
503    aknew(i)=akz(i)
504    bknew(i)=bkz(i)
505  end do
506
507  ! Switch on following lines to use doubled vertical resolution
508  !*************************************************************
509  !nz=nuvz+nwz-1
510  !if (nz.gt.nzmax) stop 'nzmax too small'
511  !do 100 i=1,nwz
512  !  aknew(2*(i-1)+1)=akm(i)
513  !00     bknew(2*(i-1)+1)=bkm(i)
514  !do 110 i=2,nuvz
515  !  aknew(2*(i-1))=akz(i)
516  !10     bknew(2*(i-1))=bkz(i)
517  ! End doubled vertical resolution
518
519
520  ! Determine the uppermost level for which the convection scheme shall be applied
521  ! by assuming that there is no convection above 50 hPa (for standard SLP)
522  !*****************************************************************************
523
524  do i=1,nuvz-2
525    pint=akz(i)+bkz(i)*101325.
526    if (pint.lt.5000.) goto 96
527  end do
52896   nconvlev=i
529  if (nconvlev.gt.nconvlevmax-1) then
530    nconvlev=nconvlevmax-1
531    write(*,*) 'Attention, convection only calculated up to ', &
532         akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa'
533  endif
534
535  return
536
537999   write(*,*)
538  write(*,*) ' ###########################################'// &
539       '###### '
540  write(*,*) '       TRAJECTORY MODEL SUBROUTINE GRIDCHECK:'
541  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn)
542  write(*,*) ' ###########################################'// &
543       '###### '
544  write(*,*)
545  write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND   !!!'
546  write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE...     !!!'
547  write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!'
548  write(*,'(a)') '!!! THE "X" KEY...                   !!!'
549  write(*,*)
550  read(*,'(a)') opt
551  if(opt.eq.'X') then
552    stop
553  else
554    goto 5
555  endif
556
557end subroutine gridcheck
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG