source: flexpart.git/src/gridcheck_nests.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: 16.3 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine gridcheck_nests
5
6  !*****************************************************************************
7  !                                                                            *
8  !     This routine checks the grid specification for the nested model        *
9  !     domains. It is similar to subroutine gridcheck, which checks the       *
10  !     mother domain.                                                         *
11  !                                                                            *
12  !     Authors: A. Stohl, G. Wotawa                                           *
13  !                                                                            *
14  !     8 February 1999                                                        *
15  !                                                                            *
16  !*****************************************************************************
17  !  CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with ECMWF grib_api    *
18  !  CHANGE: 03/12/2008, Harald Sodemann, change to f90 grib_api               *
19  !*****************************************************************************
20
21  use grib_api
22  use par_mod
23  use com_mod
24
25  implicit none
26
27  !HSO  parameters for grib_api
28  integer :: ifile
29  integer :: iret
30  integer :: igrib
31  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
32  integer :: parID !added by mc for making it consistent with new gridcheck.f90
33  integer :: gotGrib
34  !HSO  end
35  integer :: i,j,k,l,ifn,ifield,iumax,iwmax,numskip,nlev_ecn
36  integer :: nuvzn,nwzn
37  real :: akmn(nwzmax),bkmn(nwzmax),akzn(nuvzmax),bkzn(nuvzmax)
38  real(kind=4) :: xaux1,xaux2,yaux1,yaux2
39  real(kind=8) :: xaux1in,xaux2in,yaux1in,yaux2in
40  real :: conversion_factor !added by mc to make it consistent with new gridchek.f90
41
42  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
43
44  ! dimension of isec2 at least (22+n), where n is the number of parallels or
45  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
46
47  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
48  ! coordinate parameters
49
50  integer :: isec1(56),isec2(22+nxmaxn+nymaxn)
51  real(kind=4) :: zsec2(60+2*nuvzmax),zsec4(jpunp)
52
53  !HSO  grib api error messages
54  character(len=24) :: gribErrorMsg = 'Error reading grib file'
55  character(len=20) :: gribFunction = 'gridcheck_nests'
56
57  xresoln(0)=1.       ! resolution enhancement for mother grid
58  yresoln(0)=1.       ! resolution enhancement for mother grid
59
60  ! Loop about all nesting levels
61  !******************************
62
63  do l=1,numbnests
64
65    iumax=0
66    iwmax=0
67
68    if(ideltas.gt.0) then
69      ifn=1
70    else
71      ifn=numbwf
72    endif
73  !
74  ! OPENING OF DATA FILE (GRIB CODE)
75  !
76  ifile=0
77  igrib=0
78  iret=0
79
80  write(*,*) 'Reading: '//path(numpath+2*(l-1)+1) &
81         (1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,ifn))
82
835   call grib_open_file(ifile,path(numpath+2*(l-1)+1) &
84         (1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,ifn)),'r',iret)
85  if (iret.ne.GRIB_SUCCESS) then
86    goto 999   ! ERROR DETECTED
87  endif
88  !turn on support for multi fields messages
89  !call grib_multi_support_on
90
91  gotGrib=0
92  ifield=0
9310   ifield=ifield+1
94
95  !
96  ! GET NEXT FIELDS
97  !
98  call grib_new_from_file(ifile,igrib,iret)
99  if (iret.eq.GRIB_END_OF_FILE)  then
100    goto 30    ! EOF DETECTED
101  elseif (iret.ne.GRIB_SUCCESS) then
102    goto 999   ! ERROR DETECTED
103  endif
104
105  !first see if we read GRIB1 or GRIB2
106  call grib_get_int(igrib,'editionNumber',gribVer,iret)
107  call grib_check(iret,gribFunction,gribErrorMsg)
108
109  if (gribVer.eq.1) then ! GRIB Edition 1
110
111  !print*,'GRiB Edition 1'
112  !read the grib2 identifiers
113  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
114  call grib_check(iret,gribFunction,gribErrorMsg)
115  call grib_get_int(igrib,'level',isec1(8),iret)
116  call grib_check(iret,gribFunction,gribErrorMsg)
117
118  !change code for etadot to code for omega
119  if (isec1(6).eq.77) then
120    isec1(6)=135
121  endif
122
123  !print*,isec1(6),isec1(8)
124
125  else
126
127  !print*,'GRiB Edition 2'
128  !read the grib2 identifiers
129  call grib_get_int(igrib,'discipline',discipl,iret)
130  call grib_check(iret,gribFunction,gribErrorMsg)
131  call grib_get_int(igrib,'parameterCategory',parCat,iret)
132  call grib_check(iret,gribFunction,gribErrorMsg)
133  call grib_get_int(igrib,'parameterNumber',parNum,iret)
134  call grib_check(iret,gribFunction,gribErrorMsg)
135  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
136  call grib_check(iret,gribFunction,gribErrorMsg)
137  call grib_get_int(igrib,'level',valSurf,iret)
138  call grib_check(iret,gribFunction,gribErrorMsg)
139  call grib_get_int(igrib,'paramId',parId,iret) !added by mc to make it consisitent with new grid_check.f90
140  call grib_check(iret,gribFunction,gribErrorMsg) !added by mc to make it consisitent with new  grid_check.f90
141
142  !print*,discipl,parCat,parNum,typSurf,valSurf
143
144  !convert to grib1 identifiers
145  isec1(6)=-1
146  isec1(7)=-1
147  isec1(8)=-1
148  isec1(8)=valSurf     ! level
149  if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
150    isec1(6)=130         ! indicatorOfParameter
151  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
152    isec1(6)=131         ! indicatorOfParameter
153  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
154    isec1(6)=132         ! indicatorOfParameter
155  elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
156    isec1(6)=133         ! indicatorOfParameter
157  elseif ((parCat.eq.1).and.(parNum.eq.83).and.(typSurf.eq.105)) then ! clwc
158    isec1(6)=246         ! indicatorOfParameter
159  elseif ((parCat.eq.1).and.(parNum.eq.84).and.(typSurf.eq.105)) then ! ciwc
160    isec1(6)=247         ! indicatorOfParameter
161!ZHG end
162! ESO qc(=clwc+ciwc)
163  elseif ((parCat.eq.201).and.(parNum.eq.31).and.(typSurf.eq.105)) then ! qc
164    isec1(6)=201031      ! indicatorOfParameter
165  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
166    isec1(6)=134         ! indicatorOfParameter
167  elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
168    isec1(6)=135         ! indicatorOfParameter
169  elseif ((parCat.eq.128).and.(parNum.eq.77)) then ! W, actually eta dot !added bymc to make it consistent with new gridcheck.f90
170    isec1(6)=135         ! indicatorOfParameter    !
171  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
172    isec1(6)=151         ! indicatorOfParameter
173  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
174    isec1(6)=165         ! indicatorOfParameter
175  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
176    isec1(6)=166         ! indicatorOfParameter
177  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
178    isec1(6)=167         ! indicatorOfParameter
179  elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
180    isec1(6)=168         ! indicatorOfParameter
181  elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
182    isec1(6)=141         ! indicatorOfParameter
183  elseif ((parCat.eq.6).and.(parNum.eq.1) .or. parId .eq. 164) then ! CC !added by mc to make it consistent with new gridchek.f90
184    isec1(6)=164         ! indicatorOfParameter
185 elseif ((parCat.eq.1).and.(parNum.eq.9) .or. parId .eq. 142) then ! LSP !added by mc to make it consistent with new gridchek.f90
186    isec1(6)=142         ! indicatorOfParameter
187  elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
188    isec1(6)=143         ! indicatorOfParameter
189  elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
190    isec1(6)=146         ! indicatorOfParameter
191  elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
192    isec1(6)=176         ! indicatorOfParameter
193  elseif ((parCat.eq.2).and.(parNum.eq.17) .or. parId .eq. 180) then ! EWSS !added by mc to make it consistent with new gridchek.f90
194    isec1(6)=180         ! indicatorOfParameter
195  elseif ((parCat.eq.2).and.(parNum.eq.18) .or. parId .eq. 181) then ! NSSS !added by mc to make it consistent with new gridchek.f90
196    isec1(6)=181         ! indicatorOfParameter
197  elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
198    isec1(6)=129         ! indicatorOfParameter
199  elseif ((parCat.eq.3).and.(parNum.eq.7) .or. parId .eq. 160) then ! SDO !added by mc to make it consistent with new gridchek.f90
200    isec1(6)=160         ! indicatorOfParameter
201  elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
202       (typSurf.eq.1)) then ! LSM
203    isec1(6)=172         ! indicatorOfParameter
204  else
205    print*,'***ERROR: undefined GRiB2 message found!',discipl, &
206         parCat,parNum,typSurf
207  endif
208  if(parId .ne. isec1(6) .and. parId .ne. 77) then !added by mc to make it consistent with new gridchek.f90
209    write(*,*) 'parId',parId, 'isec1(6)',isec1(6)
210!    stop
211  endif
212
213  endif
214
215  !get the size and data of the values array
216  if (isec1(6).ne.-1) then
217    call grib_get_real4_array(igrib,'values',zsec4,iret)
218    call grib_check(iret,gribFunction,gribErrorMsg)
219  endif
220
221  !HSO  get the required fields from section 2 in a gribex compatible manner
222  if (ifield.eq.1) then
223    call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
224         isec2(2),iret)
225    call grib_check(iret,gribFunction,gribErrorMsg)
226    call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
227         isec2(3),iret)
228    call grib_check(iret,gribFunction,gribErrorMsg)
229    call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
230         isec2(12),iret)
231    call grib_check(iret,gribFunction,gribErrorMsg)
232  !HSO    get the size and data of the vertical coordinate array
233    call grib_get_real4_array(igrib,'pv',zsec2,iret)
234    call grib_check(iret,gribFunction,gribErrorMsg)
235
236    nxn(l)=isec2(2)
237    nyn(l)=isec2(3)
238    nlev_ecn=isec2(12)/2-1
239  endif ! ifield
240
241  if (nxn(l).gt.nxmaxn) then
242  write(*,*) 'FLEXPART error: Too many grid points in x direction.'
243  write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)'
244  write(*,*) 'for nesting level ',l
245  write(*,*) 'Or change parameter settings in file par_mod.'
246  write(*,*) nxn(l),nxmaxn
247  stop
248  endif
249
250  if (nyn(l).gt.nymaxn) then
251  write(*,*) 'FLEXPART error: Too many grid points in y direction.'
252  write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)'
253  write(*,*) 'for nesting level ',l
254  write(*,*) 'Or change parameter settings in file par_mod.'
255  write(*,*) nyn(l),nymaxn
256  stop
257  endif
258
259  !HSO  get the second part of the grid dimensions only from GRiB1 messages
260  if (isec1(6) .eq. 167 .and. (gotGrib.eq.0)) then !added by mc to make it consistent with new gridchek.f90 note that gotGrid must be changed in gotGrib!!
261    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', & !comment by mc: note that this was in the (if (ifield.eq.1) ..end above in gridchek.f90 see line 257
262         xaux1in,iret)
263    call grib_check(iret,gribFunction,gribErrorMsg)
264    call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', &
265         xaux2in,iret)
266    call grib_check(iret,gribFunction,gribErrorMsg)
267    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
268         yaux1in,iret)
269    call grib_check(iret,gribFunction,gribErrorMsg)
270    call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', &
271         yaux2in,iret)
272    call grib_check(iret,gribFunction,gribErrorMsg)
273    xaux1=xaux1in
274    xaux2=xaux2in
275    yaux1=yaux1in
276    yaux2=yaux2in
277    if(xaux1.gt.180.) xaux1=xaux1-360.0
278    if(xaux2.gt.180.) xaux2=xaux2-360.0
279    if(xaux1.lt.-180.) xaux1=xaux1+360.0
280    if(xaux2.lt.-180.) xaux2=xaux2+360.0
281    if (xaux2.lt.xaux1) xaux2=xaux2+360.0
282    xlon0n(l)=xaux1
283    ylat0n(l)=yaux1
284    dxn(l)=(xaux2-xaux1)/real(nxn(l)-1)
285    dyn(l)=(yaux2-yaux1)/real(nyn(l)-1)
286    gotGrib=1 !commetn by mc note tahthere gotGRIB is used instead of gotGrid!!!
287  endif ! ifield.eq.1
288
289  k=isec1(8)
290  if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
291  if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
292
293  if(isec1(6).eq.129) then
294    do j=0,nyn(l)-1
295      do i=0,nxn(l)-1
296        oron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
297      end do
298    end do
299  endif
300  if(isec1(6).eq.172) then
301    do j=0,nyn(l)-1
302      do i=0,nxn(l)-1
303        lsmn(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
304      end do
305    end do
306  endif
307  if(isec1(6).eq.160) then
308    do j=0,nyn(l)-1
309      do i=0,nxn(l)-1
310        excessoron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
311      end do
312    end do
313  endif
314
315  call grib_release(igrib)
316  goto 10                 !! READ NEXT LEVEL OR PARAMETER
317  !
318  ! CLOSING OF INPUT DATA FILE
319  !
320
32130   call grib_close_file(ifile)
322
323  !error message if no fields found with correct first longitude in it
324  if (gotGrib.eq.0) then
325    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
326         'messages'
327    stop
328  endif
329
330  nuvzn=iumax
331  nwzn=iwmax + 1       ! +1 needed if input contains only lower levels
332  if(nuvzn.eq.nlev_ec) nwzn=nlev_ecn+1
333
334  if ((nuvzn.gt.nuvzmax).or.(nwzn.gt.nwzmax)) then
335  write(*,*) 'FLEXPART error: Nested wind fields have too many'// &
336       'vertical levels.'
337  write(*,*) 'Problem was encountered for nesting level ',l
338  stop
339  endif
340
341
342  ! Output of grid info
343  !********************
344
345  write(*,'(a,i2,a)') ' Nested domain ',l,':'
346  write(*,'(a,f10.5,a,f10.5,a,f10.5)') '  Longitude range: ', &
347       xlon0n(l),' to ',xlon0n(l)+(nxn(l)-1)*dxn(l), &
348       '   Grid distance: ',dxn(l)
349  write(*,'(a,f10.5,a,f10.5,a,f10.5)') '  Latitude range : ', &
350       ylat0n(l),' to ',ylat0n(l)+(nyn(l)-1)*dyn(l), &
351       '   Grid distance: ',dyn(l)
352  write(*,*)
353
354  ! Determine, how much the resolutions in the nests are enhanced as
355  ! compared to the mother grid
356  !*****************************************************************
357
358  xresoln(l)=dx/dxn(l)
359  yresoln(l)=dy/dyn(l)
360
361  ! Determine the mother grid coordinates of the corner points of the
362  ! nested grids
363  ! Convert first to geographical coordinates, then to grid coordinates
364  !********************************************************************
365
366  xaux1=xlon0n(l)
367  xaux2=xlon0n(l)+real(nxn(l)-1)*dxn(l)
368  yaux1=ylat0n(l)
369  yaux2=ylat0n(l)+real(nyn(l)-1)*dyn(l)
370
371  xln(l)=(xaux1-xlon0)/dx
372  xrn(l)=(xaux2-xlon0)/dx
373  yln(l)=(yaux1-ylat0)/dy
374  yrn(l)=(yaux2-ylat0)/dy
375
376
377  if ((xln(l).lt.0.).or.(yln(l).lt.0.).or. &
378       (xrn(l).gt.real(nxmin1)).or.(yrn(l).gt.real(nymin1))) then
379    write(*,*) 'Nested domain does not fit into mother domain'
380    write(*,*) 'For global mother domain fields, you can shift'
381    write(*,*) 'shift the mother domain into x-direction'
382    write(*,*) 'by setting nxshift (file par_mod) to a'
383    write(*,*) 'positive value. Execution is terminated.'
384    stop
385  endif
386
387
388  ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL
389  ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM
390
391  numskip=nlev_ecn-nuvzn ! number of ecmwf model layers not used by FLEXPART
392  do i=1,nwzn
393    j=numskip+i
394    k=nlev_ecn+1+numskip+i
395    akmn(nwzn-i+1)=zsec2(j)
396    bkmn(nwzn-i+1)=zsec2(k)
397  end do
398
399  !
400  ! CALCULATION OF AKZ, BKZ
401  ! AKZ,BKZ: model discretization parameters at the center of each model
402  !     layer
403  !
404  ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0,
405  ! i.e. ground level
406  !*****************************************************************************
407
408    akzn(1)=0.
409    bkzn(1)=1.0
410    do i=1,nuvzn
411      akzn(i+1)=0.5*(akmn(i+1)+akmn(i))
412      bkzn(i+1)=0.5*(bkmn(i+1)+bkmn(i))
413    end do
414    nuvzn=nuvzn+1
415
416  ! Check, whether the heights of the model levels of the nested
417  ! wind fields are consistent with those of the mother domain.
418  ! If not, terminate model run.
419  !*************************************************************
420
421    do i=1,nuvz
422      if ((akzn(i).ne.akz(i)).or.(bkzn(i).ne.bkz(i))) then
423  write(*,*) 'FLEXPART error: The wind fields of nesting level',l
424  write(*,*) 'are not consistent with the mother domain:'
425  write(*,*) 'Differences in vertical levels detected.'
426        stop
427      endif
428    end do
429
430    do i=1,nwz
431      if ((akmn(i).ne.akm(i)).or.(bkmn(i).ne.bkm(i))) then
432  write(*,*) 'FLEXPART error: The wind fields of nesting level',l
433  write(*,*) 'are not consistent with the mother domain:'
434  write(*,*) 'Differences in vertical levels detected.'
435        stop
436      endif
437    end do
438
439  end do
440
441  return
442
443999   write(*,*)
444  write(*,*) ' ###########################################'// &
445       '###### '
446  write(*,*) '                FLEXPART SUBROUTINE GRIDCHECK:'
447  write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfnamen(l,ifn)
448  write(*,*) ' FOR NESTING LEVEL ',k
449  write(*,*) ' ###########################################'// &
450       '###### '
451  stop
452
453end subroutine gridcheck_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG