source: flexpart.git/src/gridcheck_gfs.f90 @ 6ecb30a

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 6ecb30a was 6ecb30a, checked in by Espen Sollum ATMOS <eso@…>, 7 years ago

Merged changes from CTBTO project

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