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