1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420
|
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
! *
! This file is part of FLEXPART. *
! *
! FLEXPART is free software: you can redistribute it and/or modify *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or *
! (at your option) any later version. *
! *
! FLEXPART is distributed in the hope that it will be useful, *
! but WITHOUT ANY WARRANTY; without even the implied warranty of *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
! GNU General Public License for more details. *
! *
! You should have received a copy of the GNU General Public License *
! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
!**********************************************************************
subroutine readwind_nests(indj,n,uuhn,vvhn,wwhn)
! i i o o o
!*****************************************************************************
! *
! This routine reads the wind fields for the nested model domains. *
! It is similar to subroutine readwind, which reads the mother domain. *
! *
! Authors: A. Stohl, G. Wotawa *
! *
! 8 February 1999 *
! *
! Last update: 17 October 2000, A. Stohl *
! *
!*****************************************************************************
! Changes, Bernd C. Krueger, Feb. 2001: *
! Variables tthn and qvhn (on eta coordinates) in common block *
! CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with ECMWF grib_api *
! CHANGE: 03/12/2008, Harald Sodemann, update to f90 with ECMWF grib_api *
!*****************************************************************************
use grib_api
use par_mod
use com_mod
implicit none
!HSO parameters for grib_api
integer :: ifile
integer :: iret
integer :: igrib
integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
integer :: gotGrid
!HSO end
real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,l
! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
! dimension of isec2 at least (22+n), where n is the number of parallels or
! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
! dimension of zsec2 at least (10+nn), where nn is the number of vertical
! coordinate parameters
integer :: isec1(56),isec2(22+nxmaxn+nymaxn)
real(kind=4) :: zsec4(jpunp)
real(kind=8) :: xauxin,yauxin
real(kind=4) :: xaux,yaux,xaux0,yaux0
real :: ewss(0:nxmaxn-1,0:nymaxn-1),nsss(0:nxmaxn-1,0:nymaxn-1)
real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
logical :: hflswitch,strswitch
!HSO grib api error messages
character(len=24) :: gribErrorMsg = 'Error reading grib file'
character(len=20) :: gribFunction = 'readwind_nests'
do l=1,numbnests
hflswitch=.false.
strswitch=.false.
levdiff2=nlev_ec-nwz+1
iumax=0
iwmax=0
ifile=0
igrib=0
iret=0
!
! OPENING OF DATA FILE (GRIB CODE)
!
5 call grib_open_file(ifile,path(numpath+2*(l-1)+1) &
(1:length(numpath+2*(l-1)+1))//trim(wfnamen(l,indj)),'r')
if (iret.ne.GRIB_SUCCESS) then
goto 888 ! ERROR DETECTED
endif
!turn on support for multi fields messages */
!call grib_multi_support_on
gotGrid=0
ifield=0
10 ifield=ifield+1
!
! GET NEXT FIELDS
!
call grib_new_from_file(ifile,igrib,iret)
if (iret.eq.GRIB_END_OF_FILE) then
goto 50 ! EOF DETECTED
elseif (iret.ne.GRIB_SUCCESS) then
goto 888 ! ERROR DETECTED
endif
!first see if we read GRIB1 or GRIB2
call grib_get_int(igrib,'editionNumber',gribVer,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
if (gribVer.eq.1) then ! GRIB Edition 1
!print*,'GRiB Edition 1'
!read the grib2 identifiers
call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'level',isec1(8),iret)
call grib_check(iret,gribFunction,gribErrorMsg)
!change code for etadot to code for omega
if (isec1(6).eq.77) then
isec1(6)=135
endif
else
!print*,'GRiB Edition 2'
!read the grib2 identifiers
call grib_get_int(igrib,'discipline',discipl,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'parameterCategory',parCat,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'parameterNumber',parNum,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'level',valSurf,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
!print*,discipl,parCat,parNum,typSurf,valSurf
!convert to grib1 identifiers
isec1(6)=-1
isec1(7)=-1
isec1(8)=-1
isec1(8)=valSurf ! level
if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
isec1(6)=130 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
isec1(6)=131 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
isec1(6)=132 ! indicatorOfParameter
elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
isec1(6)=133 ! indicatorOfParameter
elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
isec1(6)=134 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
isec1(6)=135 ! indicatorOfParameter
elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
isec1(6)=151 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
isec1(6)=165 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
isec1(6)=166 ! indicatorOfParameter
elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
isec1(6)=167 ! indicatorOfParameter
elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
isec1(6)=168 ! indicatorOfParameter
elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
isec1(6)=141 ! indicatorOfParameter
elseif ((parCat.eq.6).and.(parNum.eq.1)) then ! CC
isec1(6)=164 ! indicatorOfParameter
elseif ((parCat.eq.1).and.(parNum.eq.9)) then ! LSP
isec1(6)=142 ! indicatorOfParameter
elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
isec1(6)=143 ! indicatorOfParameter
elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
isec1(6)=146 ! indicatorOfParameter
elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
isec1(6)=176 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.17)) then ! EWSS
isec1(6)=180 ! indicatorOfParameter
elseif ((parCat.eq.2).and.(parNum.eq.18)) then ! NSSS
isec1(6)=181 ! indicatorOfParameter
elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
isec1(6)=129 ! indicatorOfParameter
elseif ((parCat.eq.3).and.(parNum.eq.7)) then ! SDO
isec1(6)=160 ! indicatorOfParameter
elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
(typSurf.eq.1)) then ! LSM
isec1(6)=172 ! indicatorOfParameter
else
print*,'***ERROR: undefined GRiB2 message found!',discipl, &
parCat,parNum,typSurf
endif
endif
!HSO get the size and data of the values array
if (isec1(6).ne.-1) then
call grib_get_real4_array(igrib,'values',zsec4,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
endif
!HSO get the required fields from section 2 in a gribex compatible manner
if(ifield.eq.1) then
call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
isec2(2),iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
isec2(3),iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
isec2(12))
call grib_check(iret,gribFunction,gribErrorMsg)
! CHECK GRID SPECIFICATIONS
if(isec2(2).ne.nxn(l)) stop &
'READWIND: NX NOT CONSISTENT FOR A NESTING LEVEL'
if(isec2(3).ne.nyn(l)) stop &
'READWIND: NY NOT CONSISTENT FOR A NESTING LEVEL'
if(isec2(12)/2-1.ne.nlev_ec) stop 'READWIND: VERTICAL DISCRET&
&IZATION NOT CONSISTENT FOR A NESTING LEVEL'
endif ! ifield
!HSO get the second part of the grid dimensions only from GRiB1 messages
if ((gribVer.eq.1).and.(gotGrid.eq.0)) then
call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
xauxin,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
yauxin,iret)
call grib_check(iret,gribFunction,gribErrorMsg)
xaux=xauxin
yaux=yauxin
xaux0=xlon0n(l)
yaux0=ylat0n(l)
if(xaux.lt.0.) xaux=xaux+360.
if(yaux.lt.0.) yaux=yaux+360.
if(xaux0.lt.0.) xaux0=xaux0+360.
if(yaux0.lt.0.) yaux0=yaux0+360.
if(xaux.ne.xaux0) &
stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT FOR A NES&
&TING LEVEL'
if(yaux.ne.yaux0) &
stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT FOR A NEST&
&ING LEVEL'
gotGrid=1
endif
do j=0,nyn(l)-1
do i=0,nxn(l)-1
k=isec1(8)
if(isec1(6).eq.130) tthn(i,j,nlev_ec-k+2,n,l)= &!! TEMPERATURE
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.131) uuhn(i,j,nlev_ec-k+2,l)= &!! U VELOCITY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.132) vvhn(i,j,nlev_ec-k+2,l)= &!! V VELOCITY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.133) then !! SPEC. HUMIDITY
qvhn(i,j,nlev_ec-k+2,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if (qvhn(i,j,nlev_ec-k+2,n,l) .lt. 0.) &
qvhn(i,j,nlev_ec-k+2,n,l) = 0.
! this is necessary because the gridded data may contain
! spurious negative values
endif
if(isec1(6).eq.134) psn(i,j,1,n,l)= &!! SURF. PRESS.
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.135) wwhn(i,j,nlev_ec-k+1,l)= &!! W VELOCITY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.141) sdn(i,j,1,n,l)= &!! SNOW DEPTH
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.151) msln(i,j,1,n,l)= &!! SEA LEVEL PRESS.
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.164) tccn(i,j,1,n,l)= &!! CLOUD COVER
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.165) u10n(i,j,1,n,l)= &!! 10 M U VELOCITY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.166) v10n(i,j,1,n,l)= &!! 10 M V VELOCITY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.167) tt2n(i,j,1,n,l)= &!! 2 M TEMPERATURE
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.168) td2n(i,j,1,n,l)= &!! 2 M DEW POINT
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.142) then !! LARGE SCALE PREC.
lsprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if (lsprecn(i,j,1,n,l).lt.0.) lsprecn(i,j,1,n,l)=0.
endif
if(isec1(6).eq.143) then !! CONVECTIVE PREC.
convprecn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if (convprecn(i,j,1,n,l).lt.0.) convprecn(i,j,1,n,l)=0.
endif
if(isec1(6).eq.146) sshfn(i,j,1,n,l)= &!! SENS. HEAT FLUX
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if((isec1(6).eq.146).and. &
(zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) hflswitch=.true. ! Heat flux available
if(isec1(6).eq.176) then !! SOLAR RADIATION
ssrn(i,j,1,n,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if (ssrn(i,j,1,n,l).lt.0.) ssrn(i,j,1,n,l)=0.
endif
if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
(zsec4(nxn(l)*(nyn(l)-j-1)+i+1).ne.0.)) strswitch=.true. ! stress available
if(isec1(6).eq.129) oron(i,j,l)= &!! ECMWF OROGRAPHY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga
if(isec1(6).eq.160) excessoron(i,j,l)= &!! STANDARD DEVIATION OF OROGRAPHY
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.172) lsmn(i,j,l)= &!! ECMWF LAND SEA MASK
zsec4(nxn(l)*(nyn(l)-j-1)+i+1)
if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
end do
end do
call grib_release(igrib)
goto 10 !! READ NEXT LEVEL OR PARAMETER
!
! CLOSING OF INPUT DATA FILE
!
50 call grib_close_file(ifile)
!error message if no fields found with correct first longitude in it
if (gotGrid.eq.0) then
print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
'messages'
stop
endif
if(levdiff2.eq.0) then
iwmax=nlev_ec+1
do i=0,nxn(l)-1
do j=0,nyn(l)-1
wwhn(i,j,nlev_ec+1,l)=0.
end do
end do
endif
do i=0,nxn(l)-1
do j=0,nyn(l)-1
surfstrn(i,j,1,n,l)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
end do
end do
if ((.not.hflswitch).or.(.not.strswitch)) then
write(*,*) 'WARNING: No flux data contained in GRIB file ', &
wfnamen(l,indj)
! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
! As ECMWF has increased the model resolution, such that now the first model
! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
! (3rd model level in FLEXPART) for the profile method
!***************************************************************************
do i=0,nxn(l)-1
do j=0,nyn(l)-1
plev1=akz(3)+bkz(3)*psn(i,j,1,n,l)
pmean=0.5*(psn(i,j,1,n,l)+plev1)
tv=tthn(i,j,3,n,l)*(1.+0.61*qvhn(i,j,3,n,l))
fu=-r_air*tv/ga/pmean
hlev1=fu*(plev1-psn(i,j,1,n,l)) ! HEIGTH OF FIRST MODEL LAYER
ff10m= sqrt(u10n(i,j,1,n,l)**2+v10n(i,j,1,n,l)**2)
fflev1=sqrt(uuhn(i,j,3,l)**2+vvhn(i,j,3,l)**2)
call pbl_profile(psn(i,j,1,n,l),td2n(i,j,1,n,l),hlev1, &
tt2n(i,j,1,n,l),tthn(i,j,3,n,l),ff10m,fflev1, &
surfstrn(i,j,1,n,l),sshfn(i,j,1,n,l))
if(sshfn(i,j,1,n,l).gt.200.) sshfn(i,j,1,n,l)=200.
if(sshfn(i,j,1,n,l).lt.-400.) sshfn(i,j,1,n,l)=-400.
end do
end do
endif
! Assign 10 m wind to model level at eta=1.0 to have one additional model
! level at the ground
! Specific humidity is taken the same as at one level above
! Temperature is taken as 2 m temperature
!**************************************************************************
do i=0,nxn(l)-1
do j=0,nyn(l)-1
uuhn(i,j,1,l)=u10n(i,j,1,n,l)
vvhn(i,j,1,l)=v10n(i,j,1,n,l)
qvhn(i,j,1,n,l)=qvhn(i,j,2,n,l)
tthn(i,j,1,n,l)=tt2n(i,j,1,n,l)
end do
end do
if(iumax.ne.nuvz-1) stop &
'READWIND: NUVZ NOT CONSISTENT FOR A NESTING LEVEL'
if(iwmax.ne.nwz) stop &
'READWIND: NWZ NOT CONSISTENT FOR A NESTING LEVEL'
end do
return
888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
write(*,*) ' #### ',wfnamen(l,indj),' FOR NESTING LEVEL #### '
write(*,*) ' #### ',l,' IS NOT GRIB FORMAT !!! #### '
stop 'Execution terminated'
999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### '
write(*,*) ' #### ',wfnamen(l,indj),' #### '
write(*,*) ' #### CANNOT BE OPENED FOR NESTING LEVEL ',l,'####'
end subroutine readwind_nests
|