source: branches/flexpart91_hasod/src_parallel/gridcheck_nests.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

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