source: flexpart.git/src/gridcheck.f90 @ b0434e1

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

Initial code to handle cloud water in nested wind fields (it is not completely implemented in this commit)

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