source: flexpart.git/src/gridcheck_ecmwf.f90

bugfixes+enhancements
Last change on this file was dba4221, checked in by Pirmin Kaufmann <pirmin.kaufmann@…>, 18 months ago

Bugfixes:

  • options/SPECIES/SPECIES_009: corrected wrong number format, replaced comma by decimal point
  • options/SPECIES/SPECIES_028: corrected wrong number format, moved sign of exponent to after the E
  • options/SPECIES/specoverview.f90: added namelist parameters that appear in SPECIES files but were missing here
  • src/FLEXPART.f90: replaced compiler-specific command line argument routines by standard Fortran intrinsic routines
  • src/FLEXPART_MPI.f90: ditto
  • src/gridcheck_ecmwf.f90: corrected handling of vertical levels when input files do not contain uppermost layers
  • src/gridcheck_nests.f90: ditto
  • src/readwind_ecmwf.f90: corrected handling of vertical levels when input files do not contain uppermost layers
  • readwind_ecmwf_mpi.f90: ditto

Code enhancements:

  • options/OUTGRID: added comments describing contents
  • options/SPECIES/SPECIES_*: aligned comments
  • options/SPECIES/specoverview.f90: removed commented lines, rectified lines indenting
  • src/FLEXPART.f90: rectified lines indenting, updated date in version string
  • src/FLEXPART_MPI.f90: ditto, and realigned code with src/FLEXPART.f90
  • src/gridcheck_*.f90: added code to write out name of file before it is opened (helps a lot when an input file causes troubles)
  • src/par_mod.f90: added comment explaining relevance of nuvzmax for GRIB input
  • src/readreleases.f90: write out warning if too few particles are used to randomize release
  • src/readspecies.f90: write out name of SPECIES file before it is read
  • src/readwind_*.f90: write out name of input file before opening it
  • src/writeheader_txt.f90: removed wrong comment
  • Property mode set to 100644
File size: 20.8 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  !
105  write(*,*)   'Reading: '//path(3)(1:length(3)) &
106       //trim(wfname(ifn))
1075 call grib_open_file(ifile,path(3)(1:length(3)) &
108       //trim(wfname(ifn)),'r',iret)
109  if (iret.ne.GRIB_SUCCESS) then
110    goto 999   ! ERROR DETECTED
111  endif
112  !turn on support for multi fields messages
113  !call grib_multi_support_on
114
115  gotGrid=0
116  ifield=0
11710 ifield=ifield+1
118
119  !
120  ! GET NEXT FIELDS
121  !
122  call grib_new_from_file(ifile,igrib,iret)
123  if (iret.eq.GRIB_END_OF_FILE )  then
124    goto 30    ! EOF DETECTED
125  elseif (iret.ne.GRIB_SUCCESS) then
126    goto 999   ! ERROR DETECTED
127  endif
128
129  !first see if we read GRIB1 or GRIB2
130  call grib_get_int(igrib,'editionNumber',gribVer,iret)
131  call grib_check(iret,gribFunction,gribErrorMsg)
132
133  if (gribVer.eq.1) then ! GRIB Edition 1
134
135    !print*,'GRiB Edition 1'
136    !read the grib2 identifiers
137    call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
138    call grib_check(iret,gribFunction,gribErrorMsg)
139    call grib_get_int(igrib,'level',isec1(8),iret)
140    call grib_check(iret,gribFunction,gribErrorMsg)
141
142    !change code for etadot to code for omega
143    if (isec1(6).eq.77) then
144      isec1(6)=135
145    endif
146
147    !print*,isec1(6),isec1(8)
148
149  else
150
151    !print*,'GRiB Edition 2'
152    !read the grib2 identifiers
153    call grib_get_int(igrib,'discipline',discipl,iret)
154    call grib_check(iret,gribFunction,gribErrorMsg)
155    call grib_get_int(igrib,'parameterCategory',parCat,iret)
156    call grib_check(iret,gribFunction,gribErrorMsg)
157    call grib_get_int(igrib,'parameterNumber',parNum,iret)
158    call grib_check(iret,gribFunction,gribErrorMsg)
159    call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
160    call grib_check(iret,gribFunction,gribErrorMsg)
161    call grib_get_int(igrib,'level',valSurf,iret)
162    call grib_check(iret,gribFunction,gribErrorMsg)
163    call grib_get_int(igrib,'paramId',parId,iret)
164    call grib_check(iret,gribFunction,gribErrorMsg)
165
166    !print*,discipl,parCat,parNum,typSurf,valSurf
167
168    !convert to grib1 identifiers
169    isec1(6)=-1
170    isec1(7)=-1
171    isec1(8)=-1
172    isec1(8)=valSurf     ! level
173    if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
174      isec1(6)=130         ! indicatorOfParameter
175    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
176      isec1(6)=131         ! indicatorOfParameter
177    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
178      isec1(6)=132         ! indicatorOfParameter
179    elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
180      isec1(6)=133         ! indicatorOfParameter
181      !ZHG FOR CLOUDS FROM GRIB
182    elseif ((parCat.eq.1).and.(parNum.eq.83).and.(typSurf.eq.105)) then ! clwc
183      isec1(6)=246         ! indicatorOfParameter
184    elseif ((parCat.eq.1).and.(parNum.eq.84).and.(typSurf.eq.105)) then ! ciwc
185      isec1(6)=247         ! indicatorOfParameter
186      !ZHG end
187      ! ESO qc(=clwc+ciwc)
188    elseif ((parCat.eq.201).and.(parNum.eq.31).and.(typSurf.eq.105)) then ! qc
189      isec1(6)=201031      ! indicatorOfParameter
190    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
191      isec1(6)=134         ! indicatorOfParameter
192    elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
193      isec1(6)=135         ! indicatorOfParameter
194    elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually eta dot
195      isec1(6)=135         ! indicatorOfParameter
196    elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
197      isec1(6)=151         ! indicatorOfParameter
198    elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
199      isec1(6)=165         ! indicatorOfParameter
200    elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
201      isec1(6)=166         ! indicatorOfParameter
202    elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
203      isec1(6)=167         ! indicatorOfParameter
204    elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
205      isec1(6)=168         ! indicatorOfParameter
206    elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
207      isec1(6)=141         ! indicatorOfParameter
208    elseif ((parCat.eq.6).and.(parNum.eq.1) .or. parId .eq. 164) then ! CC
209      isec1(6)=164         ! indicatorOfParameter
210    elseif ((parCat.eq.1).and.(parNum.eq.9) .or. parId .eq. 142) then ! LSP
211      isec1(6)=142         ! indicatorOfParameter
212    elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
213      isec1(6)=143         ! indicatorOfParameter
214    elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
215      isec1(6)=146         ! indicatorOfParameter
216    elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
217      isec1(6)=176         ! indicatorOfParameter
218    elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS
219      isec1(6)=180         ! indicatorOfParameter
220    elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS
221      isec1(6)=181         ! indicatorOfParameter
222    elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
223      isec1(6)=129         ! indicatorOfParameter
224    elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO
225      isec1(6)=160         ! indicatorOfParameter
226    elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
227         (typSurf.eq.1)) then ! LSM
228      isec1(6)=172         ! indicatorOfParameter
229    else
230      print*,'***ERROR: undefined GRiB2 message found!',discipl, &
231           parCat,parNum,typSurf
232    endif
233    if(parId .ne. isec1(6) .and. parId .ne. 77) then
234      write(*,*) 'parId',parId, 'isec1(6)',isec1(6)
235      !    stop
236    endif
237
238  endif
239
240  CALL grib_get_int(igrib,'numberOfPointsAlongAParallel', &
241       isec2(2),iret)
242  ! nx=isec2(2)
243  ! WRITE(*,*) nx,nxmax
244  IF (isec2(2).GT.nxmax) THEN
245    WRITE(*,*) 'FLEXPART error: Too many grid points in x direction.'
246    WRITE(*,*) 'Reduce resolution of wind fields.'
247    WRITE(*,*) 'Or change parameter settings in file ecmwf_mod.'
248    WRITE(*,*) nx,nxmax
249!    STOP
250  ENDIF
251
252  !get the size and data of the values array
253  if (isec1(6).ne.-1) then
254    call grib_get_real4_array(igrib,'values',zsec4,iret)
255    call grib_check(iret,gribFunction,gribErrorMsg)
256  endif
257
258  if (ifield.eq.1) then
259
260    !HSO  get the required fields from section 2 in a gribex compatible manner
261    call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
262         isec2(2),iret)
263    call grib_check(iret,gribFunction,gribErrorMsg)
264    call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
265         isec2(3),iret)
266    call grib_check(iret,gribFunction,gribErrorMsg)
267    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
268         xaux1in,iret)
269    call grib_check(iret,gribFunction,gribErrorMsg)
270    call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
271         isec2(12),iret)
272    call grib_check(iret,gribFunction,gribErrorMsg)
273
274    !  get the size and data of the vertical coordinate array
275    call grib_get_real4_array(igrib,'pv',zsec2,iret)
276    call grib_check(iret,gribFunction,gribErrorMsg)
277
278    nxfield=isec2(2)
279    ny=isec2(3)
280    nlev_ec=isec2(12)/2-1
281  endif
282
283  !HSO  get the second part of the grid dimensions only from GRiB1 messages
284  if (isec1(6) .eq. 167 .and. (gotGrid.eq.0)) then
285    call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', &
286         xaux2in,iret)
287    call grib_check(iret,gribFunction,gribErrorMsg)
288    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
289         yaux1in,iret)
290    call grib_check(iret,gribFunction,gribErrorMsg)
291    call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', &
292         yaux2in,iret)
293    call grib_check(iret,gribFunction,gribErrorMsg)
294    xaux1=xaux1in
295    xaux2=xaux2in
296    yaux1=yaux1in
297    yaux2=yaux2in
298    if (xaux1.gt.180.) xaux1=xaux1-360.0
299    if (xaux2.gt.180.) xaux2=xaux2-360.0
300    if (xaux1.lt.-180.) xaux1=xaux1+360.0
301    if (xaux2.lt.-180.) xaux2=xaux2+360.0
302    if (xaux2.lt.xaux1) xaux2=xaux2+360.0
303    xlon0=xaux1
304    ylat0=yaux1
305    dx=(xaux2-xaux1)/real(nxfield-1)
306    dy=(yaux2-yaux1)/real(ny-1)
307    dxconst=180./(dx*r_earth*pi)
308    dyconst=180./(dy*r_earth*pi)
309    gotGrid=1
310    ! Check whether fields are global
311    ! If they contain the poles, specify polar stereographic map
312    ! projections using the stlmbr- and stcm2p-calls
313    !***********************************************************
314
315    xauxa=abs(xaux2+dx-360.-xaux1)
316    if (xauxa.lt.0.001) then
317      nx=nxfield+1                 ! field is cyclic
318      xglobal=.true.
319      if (abs(nxshift).ge.nx) &
320           stop 'nxshift in file par_mod is too large'
321      xlon0=xlon0+real(nxshift)*dx
322    else
323      nx=nxfield
324      xglobal=.false.
325      if (nxshift.ne.0) &
326           stop 'nxshift (par_mod) must be zero for non-global domain'
327    endif
328    nxmin1=nx-1
329    nymin1=ny-1
330    if (xlon0.gt.180.) xlon0=xlon0-360.
331    xauxa=abs(yaux1+90.)
332    if (xglobal.and.xauxa.lt.0.001) then
333      sglobal=.true.               ! field contains south pole
334      ! Enhance the map scale by factor 3 (*2=6) compared to north-south
335      ! map scale
336      sizesouth=6.*(switchsouth+90.)/dy
337      call stlmbr(southpolemap,-90.,0.)
338      call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, &
339           sizesouth,switchsouth,180.)
340      switchsouthg=(switchsouth-ylat0)/dy
341    else
342      sglobal=.false.
343      switchsouthg=999999.
344    endif
345    xauxa=abs(yaux2-90.)
346    if (xglobal.and.xauxa.lt.0.001) then
347      nglobal=.true.               ! field contains north pole
348      ! Enhance the map scale by factor 3 (*2=6) compared to north-south
349      ! map scale
350      sizenorth=6.*(90.-switchnorth)/dy
351      call stlmbr(northpolemap,90.,0.)
352      call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, &
353           sizenorth,switchnorth,180.)
354      switchnorthg=(switchnorth-ylat0)/dy
355    else
356      nglobal=.false.
357      switchnorthg=999999.
358    endif
359    if (nxshift.lt.0) &
360         stop 'nxshift (par_mod) must not be negative'
361    if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large'
362  endif ! gotGrid
363
364  if (nx.gt.nxmax) then
365    write(*,*) 'FLEXPART error: Too many grid points in x direction.'
366    write(*,*) 'Reduce resolution of wind fields.'
367    write(*,*) 'Or change parameter settings in file par_mod.'
368    write(*,*) nx,nxmax
369    stop
370  endif
371
372  if (ny.gt.nymax) then
373    write(*,*) 'FLEXPART error: Too many grid points in y direction.'
374    write(*,*) 'Reduce resolution of wind fields.'
375    write(*,*) 'Or change parameter settings in file par_mod.'
376    write(*,*) ny,nymax
377    stop
378  endif
379
380  k=isec1(8)
381  if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
382  if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
383
384  if(isec1(6).eq.129) then
385    do jy=0,ny-1
386      do ix=0,nxfield-1
387        oro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)/ga
388      end do
389    end do
390  endif
391  if(isec1(6).eq.172) then
392    do jy=0,ny-1
393      do ix=0,nxfield-1
394        lsm(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
395      end do
396    end do
397  endif
398  if(isec1(6).eq.160) then
399    do jy=0,ny-1
400      do ix=0,nxfield-1
401        excessoro(ix,jy)=zsec4(nxfield*(ny-jy-1)+ix+1)
402      end do
403    end do
404  endif
405
406  call grib_release(igrib)
407  goto 10                      !! READ NEXT LEVEL OR PARAMETER
408  !
409  ! CLOSING OF INPUT DATA FILE
410  !
411
41230 call grib_close_file(ifile)
413
414  !error message if no fields found with correct first longitude in it
415  if (gotGrid.eq.0) then
416    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
417         'messages'
418    stop
419  endif
420
421  nuvz=iumax
422  nwz =iwmax + 1       ! +1 needed if input contains only lower levels
423  if(nuvz.eq.nlev_ec) nwz=nlev_ec+1
424
425  if (nuvz+1.gt.nuvzmax) then
426    write(*,*) 'FLEXPART error: Too many u,v grid points in z '// &
427         'direction.'
428    write(*,*) 'Reduce resolution of wind fields.'
429    write(*,*) 'Or change parameter settings in file par_mod.'
430    write(*,*) nuvz+1,nuvzmax
431    stop
432  endif
433
434  if (nwz.gt.nwzmax) then
435    write(*,*) 'FLEXPART error: Too many w grid points in z '// &
436         'direction.'
437    write(*,*) 'Reduce resolution of wind fields.'
438    write(*,*) 'Or change parameter settings in file par_mod.'
439    write(*,*) nwz,nwzmax
440    stop
441  endif
442
443  ! If desired, shift all grids by nxshift grid cells
444  !**************************************************
445
446  if (xglobal) then
447    call shift_field_0(oro,nxfield,ny)
448    call shift_field_0(lsm,nxfield,ny)
449    call shift_field_0(excessoro,nxfield,ny)
450  endif
451
452  ! Output of grid info
453  !********************
454
455  if (lroot) then
456    write(*,'(a,2i7)') ' Vertical levels in ECMWF data: ', &
457         nuvz+1,nwz
458    write(*,*)
459    write(*,'(a)') ' Mother domain:'
460    write(*,'(a,f10.5,a,f10.5,a,f10.5)') '  Longitude range: ', &
461         xlon0,' to ',xlon0+(nx-1)*dx,'   Grid distance: ',dx
462    write(*,'(a,f10.5,a,f10.5,a,f10.5)') '  Latitude range : ', &
463         ylat0,' to ',ylat0+(ny-1)*dy,'   Grid distance: ',dy
464    write(*,*)
465  end if
466
467  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
468  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
469
470  numskip=nlev_ec-nuvz  ! number of ecmwf model layers not used
471  ! by trajectory model
472  !do 8940 i=1,244
473  !   write (*,*) 'zsec2:',i,ifield,zsec2(i),numskip
474  !940  continue
475  !   stop
476  ! SEC SEC SEC
477  ! for unknown reason zsec 1 to 10 is filled in this version
478  ! compared to the old one
479  ! SEC SEC SE
480  do i=1,nwz
481    j=numskip+i
482    k=nlev_ec+1+numskip+i
483    akm(nwz-i+1)=zsec2(j)
484    !   write (*,*) 'ifield:',ifield,k,j,zsec2(10+j)
485    bkm(nwz-i+1)=zsec2(k)
486  end do
487
488  !
489  ! CALCULATION OF AKZ, BKZ
490  ! AKZ,BKZ: model discretization parameters at the center of each model
491  !     layer
492  !
493  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
494  ! i.e. ground level
495  !*****************************************************************************
496
497  akz(1)=0.
498  bkz(1)=1.0
499  do i=1,nuvz
500    akz(i+1)=0.5*(akm(i+1)+akm(i))
501    bkz(i+1)=0.5*(bkm(i+1)+bkm(i))
502  end do
503  nuvz=nuvz+1
504
505  ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled
506  ! upon the transformation to z levels. In order to save computer memory, this is
507  ! not done anymore in the standard version. However, this option can still be
508  ! switched on by replacing the following lines with those below, that are
509  ! currently commented out. For this, similar changes are necessary in
510  ! verttransform.f and verttranform_nests.f
511  !*****************************************************************************
512
513  nz=nuvz
514  if (nz.gt.nzmax) stop 'nzmax too small'
515  do i=1,nuvz
516    aknew(i)=akz(i)
517    bknew(i)=bkz(i)
518  end do
519
520  ! Switch on following lines to use doubled vertical resolution
521  !*************************************************************
522  !nz=nuvz+nwz-1
523  !if (nz.gt.nzmax) stop 'nzmax too small'
524  !do 100 i=1,nwz
525  !  aknew(2*(i-1)+1)=akm(i)
526  !00     bknew(2*(i-1)+1)=bkm(i)
527  !do 110 i=2,nuvz
528  !  aknew(2*(i-1))=akz(i)
529  !10     bknew(2*(i-1))=bkz(i)
530  ! End doubled vertical resolution
531
532
533  ! Determine the uppermost level for which the convection scheme shall be applied
534  ! by assuming that there is no convection above 50 hPa (for standard SLP)
535  !*****************************************************************************
536
537  do i=1,nuvz-2
538    pint=akz(i)+bkz(i)*101325.
539    if (pint.lt.5000.) goto 96
540  end do
54196 nconvlev=i
542  if (nconvlev.gt.nconvlevmax-1) then
543    nconvlev=nconvlevmax-1
544    write(*,*) 'Attention, convection only calculated up to ', &
545         akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa'
546  endif
547
548  return
549
550999 write(*,*)
551  write(*,*) ' ###########################################'// &
552       '###### '
553  write(*,*) '       TRAJECTORY MODEL SUBROUTINE GRIDCHECK:'
554  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn)
555  write(*,*) ' ###########################################'// &
556       '###### '
557  write(*,*)
558  write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND   !!!'
559  write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE...     !!!'
560  write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!'
561  write(*,'(a)') '!!! THE "X" KEY...                   !!!'
562  write(*,*)
563  read(*,'(a)') opt
564  if(opt.eq.'X') then
565    stop
566  else
567    goto 5
568  endif
569
570end subroutine gridcheck_ecmwf
571
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG