source: flexpart.git/src/gridcheck_ecmwf.f90 @ 02095e3

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

Adding files from CTBTO project wo_17

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