source: flexpart.git/src/gridcheck_ecmwf.f90 @ 3ea93bb

GFS_025dev
Last change on this file since 3ea93bb was 3ea93bb, checked in by Ignacio Pisso <ip@…>, 4 years ago

Merge branch 'release-10' into dev

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