[e200b7a] | 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 | |
---|
[6ecb30a] | 22 | subroutine gridcheck_gfs |
---|
[e200b7a] | 23 | |
---|
| 24 | !********************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! FLEXPART MODEL SUBROUTINE GRIDCHECK * |
---|
| 27 | ! * |
---|
| 28 | !********************************************************************** |
---|
| 29 | ! * |
---|
| 30 | ! AUTHOR: G. WOTAWA * |
---|
| 31 | ! DATE: 1997-08-06 * |
---|
| 32 | ! LAST UPDATE: 1997-10-10 * |
---|
| 33 | ! * |
---|
| 34 | ! Update: 1999-02-08, global fields allowed, A. Stohl* |
---|
| 35 | ! CHANGE: 17/11/2005, Caroline Forster, GFS data * |
---|
| 36 | ! CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with * |
---|
| 37 | ! ECMWF grib_api * |
---|
| 38 | ! CHANGE: 03/12/2008, Harald Sodemann, update to f90 with * |
---|
| 39 | ! ECMWF grib_api * |
---|
| 40 | ! * |
---|
[6ecb30a] | 41 | ! Unified ECMWF and GFS builds * |
---|
| 42 | ! Marian Harustak, 12.5.2017 * |
---|
| 43 | ! - Renamed routine from gridcheck to gridcheck_gfs * |
---|
| 44 | ! * |
---|
[e200b7a] | 45 | !********************************************************************** |
---|
| 46 | ! * |
---|
| 47 | ! DESCRIPTION: * |
---|
| 48 | ! * |
---|
| 49 | ! THIS SUBROUTINE DETERMINES THE GRID SPECIFICATIONS (LOWER LEFT * |
---|
| 50 | ! LONGITUDE, LOWER LEFT LATITUDE, NUMBER OF GRID POINTS, GRID DIST- * |
---|
| 51 | ! ANCE AND VERTICAL DISCRETIZATION OF THE ECMWF MODEL) FROM THE * |
---|
| 52 | ! GRIB HEADER OF THE FIRST INPUT FILE. THE CONSISTANCY (NO CHANGES * |
---|
| 53 | ! WITHIN ONE FLEXPART RUN) IS CHECKED IN THE ROUTINE "READWIND" AT * |
---|
| 54 | ! ANY CALL. * |
---|
| 55 | ! * |
---|
| 56 | ! XLON0 geographical longitude of lower left gridpoint * |
---|
| 57 | ! YLAT0 geographical latitude of lower left gridpoint * |
---|
| 58 | ! NX number of grid points x-direction * |
---|
| 59 | ! NY number of grid points y-direction * |
---|
| 60 | ! DX grid distance x-direction * |
---|
| 61 | ! DY grid distance y-direction * |
---|
| 62 | ! NUVZ number of grid points for horizontal wind * |
---|
| 63 | ! components in z direction * |
---|
| 64 | ! NWZ number of grid points for vertical wind * |
---|
| 65 | ! component in z direction * |
---|
| 66 | ! sizesouth, sizenorth give the map scale (i.e. number of virtual grid* |
---|
| 67 | ! points of the polar stereographic grid): * |
---|
| 68 | ! used to check the CFL criterion * |
---|
| 69 | ! UVHEIGHT(1)- heights of gridpoints where u and v are * |
---|
| 70 | ! UVHEIGHT(NUVZ) given * |
---|
| 71 | ! WHEIGHT(1)- heights of gridpoints where w is given * |
---|
| 72 | ! WHEIGHT(NWZ) * |
---|
| 73 | ! * |
---|
| 74 | !********************************************************************** |
---|
| 75 | |
---|
| 76 | use grib_api |
---|
| 77 | use par_mod |
---|
| 78 | use com_mod |
---|
| 79 | use conv_mod |
---|
| 80 | use cmapf_mod, only: stlmbr,stcm2p |
---|
| 81 | |
---|
| 82 | implicit none |
---|
| 83 | |
---|
| 84 | !HSO parameters for grib_api |
---|
| 85 | integer :: ifile |
---|
| 86 | integer :: iret |
---|
| 87 | integer :: igrib |
---|
| 88 | real(kind=4) :: xaux1,xaux2,yaux1,yaux2 |
---|
| 89 | real(kind=8) :: xaux1in,xaux2in,yaux1in,yaux2in |
---|
| 90 | integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl |
---|
| 91 | !HSO end |
---|
| 92 | integer :: ix,jy,i,ifn,ifield,j,k,iumax,iwmax,numskip |
---|
| 93 | real :: sizesouth,sizenorth,xauxa,pint |
---|
| 94 | real :: akm_usort(nwzmax) |
---|
| 95 | real,parameter :: eps=0.0001 |
---|
| 96 | |
---|
| 97 | ! NCEP GFS |
---|
| 98 | real :: pres(nwzmax), help |
---|
| 99 | |
---|
| 100 | integer :: i179,i180,i181 |
---|
| 101 | |
---|
| 102 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
| 103 | |
---|
| 104 | integer :: isec1(8),isec2(3) |
---|
| 105 | real(kind=4) :: zsec4(jpunp) |
---|
| 106 | character(len=1) :: opt |
---|
| 107 | |
---|
| 108 | !HSO grib api error messages |
---|
| 109 | character(len=24) :: gribErrorMsg = 'Error reading grib file' |
---|
| 110 | character(len=20) :: gribFunction = 'gridcheckwind_gfs' |
---|
| 111 | ! |
---|
| 112 | if (numbnests.ge.1) then |
---|
| 113 | write(*,*) ' ###########################################' |
---|
| 114 | write(*,*) ' FLEXPART ERROR SUBROUTINE GRIDCHECK:' |
---|
| 115 | write(*,*) ' NO NESTED WINDFIELDAS ALLOWED FOR GFS! ' |
---|
| 116 | write(*,*) ' ###########################################' |
---|
| 117 | stop |
---|
| 118 | endif |
---|
| 119 | |
---|
| 120 | iumax=0 |
---|
| 121 | iwmax=0 |
---|
| 122 | |
---|
| 123 | if(ideltas.gt.0) then |
---|
| 124 | ifn=1 |
---|
| 125 | else |
---|
| 126 | ifn=numbwf |
---|
| 127 | endif |
---|
| 128 | ! |
---|
| 129 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
| 130 | ! |
---|
| 131 | 5 call grib_open_file(ifile,path(3)(1:length(3)) & |
---|
| 132 | //trim(wfname(ifn)),'r',iret) |
---|
| 133 | if (iret.ne.GRIB_SUCCESS) then |
---|
| 134 | goto 999 ! ERROR DETECTED |
---|
| 135 | endif |
---|
| 136 | !turn on support for multi fields messages |
---|
| 137 | call grib_multi_support_on |
---|
| 138 | |
---|
| 139 | ifield=0 |
---|
| 140 | 10 ifield=ifield+1 |
---|
| 141 | ! |
---|
| 142 | ! GET NEXT FIELDS |
---|
| 143 | ! |
---|
| 144 | call grib_new_from_file(ifile,igrib,iret) |
---|
| 145 | if (iret.eq.GRIB_END_OF_FILE ) then |
---|
| 146 | goto 30 ! EOF DETECTED |
---|
| 147 | elseif (iret.ne.GRIB_SUCCESS) then |
---|
| 148 | goto 999 ! ERROR DETECTED |
---|
| 149 | endif |
---|
| 150 | |
---|
| 151 | !first see if we read GRIB1 or GRIB2 |
---|
| 152 | call grib_get_int(igrib,'editionNumber',gribVer,iret) |
---|
| 153 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 154 | |
---|
| 155 | if (gribVer.eq.1) then ! GRIB Edition 1 |
---|
| 156 | |
---|
| 157 | !read the grib1 identifiers |
---|
| 158 | call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret) |
---|
| 159 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 160 | call grib_get_int(igrib,'indicatorOfTypeOfLevel',isec1(7),iret) |
---|
| 161 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 162 | call grib_get_int(igrib,'level',isec1(8),iret) |
---|
| 163 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 164 | |
---|
| 165 | !get the size and data of the values array |
---|
| 166 | call grib_get_real4_array(igrib,'values',zsec4,iret) |
---|
| 167 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 168 | |
---|
| 169 | else ! GRIB Edition 2 |
---|
| 170 | |
---|
| 171 | !read the grib2 identifiers |
---|
| 172 | call grib_get_int(igrib,'discipline',discipl,iret) |
---|
| 173 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 174 | call grib_get_int(igrib,'parameterCategory',parCat,iret) |
---|
| 175 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 176 | call grib_get_int(igrib,'parameterNumber',parNum,iret) |
---|
| 177 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 178 | call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret) |
---|
| 179 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 180 | call grib_get_int(igrib,'scaledValueOfFirstFixedSurface', & |
---|
| 181 | valSurf,iret) |
---|
| 182 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 183 | |
---|
| 184 | !convert to grib1 identifiers |
---|
| 185 | isec1(6)=-1 |
---|
| 186 | isec1(7)=-1 |
---|
| 187 | isec1(8)=-1 |
---|
| 188 | if ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.100)) then ! U |
---|
| 189 | isec1(6)=33 ! indicatorOfParameter |
---|
| 190 | isec1(7)=100 ! indicatorOfTypeOfLevel |
---|
| 191 | isec1(8)=valSurf/100 ! level, convert to hPa |
---|
| 192 | elseif ((parCat.eq.3).and.(parNum.eq.5).and.(typSurf.eq.1)) then ! TOPO |
---|
| 193 | isec1(6)=7 ! indicatorOfParameter |
---|
| 194 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
| 195 | isec1(8)=0 |
---|
| 196 | elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.1) & |
---|
| 197 | .and.(discipl.eq.2)) then ! LSM |
---|
| 198 | isec1(6)=81 ! indicatorOfParameter |
---|
| 199 | isec1(7)=1 ! indicatorOfTypeOfLevel |
---|
| 200 | isec1(8)=0 |
---|
| 201 | endif |
---|
| 202 | |
---|
| 203 | if (isec1(6).ne.-1) then |
---|
| 204 | ! get the size and data of the values array |
---|
| 205 | call grib_get_real4_array(igrib,'values',zsec4,iret) |
---|
| 206 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 207 | endif |
---|
| 208 | |
---|
| 209 | endif ! gribVer |
---|
| 210 | |
---|
| 211 | if(ifield.eq.1) then |
---|
| 212 | |
---|
| 213 | !get the required fields from section 2 |
---|
| 214 | !store compatible to gribex input |
---|
| 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_real8(igrib,'longitudeOfFirstGridPointInDegrees', & |
---|
| 222 | xaux1in,iret) |
---|
| 223 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 224 | call grib_get_real8(igrib,'longitudeOfLastGridPointInDegrees', & |
---|
| 225 | xaux2in,iret) |
---|
| 226 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 227 | call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', & |
---|
| 228 | yaux1in,iret) |
---|
| 229 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 230 | call grib_get_real8(igrib,'latitudeOfFirstGridPointInDegrees', & |
---|
| 231 | yaux2in,iret) |
---|
| 232 | call grib_check(iret,gribFunction,gribErrorMsg) |
---|
| 233 | |
---|
[d8eed02] | 234 | ! Fix for flexpart.eu ticket #48 |
---|
| 235 | if (xaux2in.lt.0) xaux2in = 359.0 |
---|
| 236 | |
---|
[e200b7a] | 237 | xaux1=xaux1in |
---|
| 238 | xaux2=xaux2in |
---|
| 239 | yaux1=yaux1in |
---|
| 240 | yaux2=yaux2in |
---|
| 241 | |
---|
| 242 | nxfield=isec2(2) |
---|
| 243 | ny=isec2(3) |
---|
| 244 | if((abs(xaux1).lt.eps).and.(xaux2.ge.359)) then ! NCEP DATA FROM 0 TO |
---|
| 245 | xaux1=-179.0 ! 359 DEG EAST -> |
---|
| 246 | xaux2=-179.0+360.-360./real(nxfield) ! TRANSFORMED TO -179 |
---|
| 247 | endif ! TO 180 DEG EAST |
---|
| 248 | if (xaux1.gt.180) xaux1=xaux1-360.0 |
---|
| 249 | if (xaux2.gt.180) xaux2=xaux2-360.0 |
---|
| 250 | if (xaux1.lt.-180) xaux1=xaux1+360.0 |
---|
| 251 | if (xaux2.lt.-180) xaux2=xaux2+360.0 |
---|
| 252 | if (xaux2.lt.xaux1) xaux2=xaux2+360. |
---|
| 253 | xlon0=xaux1 |
---|
| 254 | ylat0=yaux1 |
---|
| 255 | dx=(xaux2-xaux1)/real(nxfield-1) |
---|
| 256 | dy=(yaux2-yaux1)/real(ny-1) |
---|
| 257 | dxconst=180./(dx*r_earth*pi) |
---|
| 258 | dyconst=180./(dy*r_earth*pi) |
---|
| 259 | !HSO end edits |
---|
| 260 | |
---|
| 261 | |
---|
| 262 | ! Check whether fields are global |
---|
| 263 | ! If they contain the poles, specify polar stereographic map |
---|
| 264 | ! projections using the stlmbr- and stcm2p-calls |
---|
| 265 | !*********************************************************** |
---|
| 266 | |
---|
| 267 | xauxa=abs(xaux2+dx-360.-xaux1) |
---|
| 268 | if (xauxa.lt.0.001) then |
---|
| 269 | nx=nxfield+1 ! field is cyclic |
---|
| 270 | xglobal=.true. |
---|
| 271 | if (abs(nxshift).ge.nx) & |
---|
| 272 | stop 'nxshift in file par_mod is too large' |
---|
| 273 | xlon0=xlon0+real(nxshift)*dx |
---|
| 274 | else |
---|
| 275 | nx=nxfield |
---|
| 276 | xglobal=.false. |
---|
| 277 | if (nxshift.ne.0) & |
---|
| 278 | stop 'nxshift (par_mod) must be zero for non-global domain' |
---|
| 279 | endif |
---|
| 280 | nxmin1=nx-1 |
---|
| 281 | nymin1=ny-1 |
---|
| 282 | if (xlon0.gt.180.) xlon0=xlon0-360. |
---|
| 283 | xauxa=abs(yaux1+90.) |
---|
| 284 | if (xglobal.and.xauxa.lt.0.001) then |
---|
| 285 | sglobal=.true. ! field contains south pole |
---|
| 286 | ! Enhance the map scale by factor 3 (*2=6) compared to north-south |
---|
| 287 | ! map scale |
---|
| 288 | sizesouth=6.*(switchsouth+90.)/dy |
---|
| 289 | call stlmbr(southpolemap,-90.,0.) |
---|
| 290 | call stcm2p(southpolemap,0.,0.,switchsouth,0.,sizesouth, & |
---|
| 291 | sizesouth,switchsouth,180.) |
---|
| 292 | switchsouthg=(switchsouth-ylat0)/dy |
---|
| 293 | else |
---|
| 294 | sglobal=.false. |
---|
| 295 | switchsouthg=999999. |
---|
| 296 | endif |
---|
| 297 | xauxa=abs(yaux2-90.) |
---|
| 298 | if (xglobal.and.xauxa.lt.0.001) then |
---|
| 299 | nglobal=.true. ! field contains north pole |
---|
| 300 | ! Enhance the map scale by factor 3 (*2=6) compared to north-south |
---|
| 301 | ! map scale |
---|
| 302 | sizenorth=6.*(90.-switchnorth)/dy |
---|
| 303 | call stlmbr(northpolemap,90.,0.) |
---|
| 304 | call stcm2p(northpolemap,0.,0.,switchnorth,0.,sizenorth, & |
---|
| 305 | sizenorth,switchnorth,180.) |
---|
| 306 | switchnorthg=(switchnorth-ylat0)/dy |
---|
| 307 | else |
---|
| 308 | nglobal=.false. |
---|
| 309 | switchnorthg=999999. |
---|
| 310 | endif |
---|
| 311 | endif ! ifield.eq.1 |
---|
| 312 | |
---|
| 313 | if (nxshift.lt.0) stop 'nxshift (par_mod) must not be negative' |
---|
| 314 | if (nxshift.ge.nxfield) stop 'nxshift (par_mod) too large' |
---|
| 315 | |
---|
| 316 | ! NCEP ISOBARIC LEVELS |
---|
| 317 | !********************* |
---|
| 318 | |
---|
| 319 | if((isec1(6).eq.33).and.(isec1(7).eq.100)) then ! check for U wind |
---|
| 320 | iumax=iumax+1 |
---|
| 321 | pres(iumax)=real(isec1(8))*100.0 |
---|
| 322 | endif |
---|
| 323 | |
---|
| 324 | |
---|
| 325 | i179=nint(179./dx) |
---|
| 326 | if (dx.lt.0.7) then |
---|
| 327 | i180=nint(180./dx)+1 ! 0.5 deg data |
---|
| 328 | else |
---|
| 329 | i180=nint(179./dx)+1 ! 1 deg data |
---|
| 330 | endif |
---|
| 331 | i181=i180+1 |
---|
| 332 | |
---|
| 333 | |
---|
| 334 | ! NCEP TERRAIN |
---|
| 335 | !************* |
---|
| 336 | |
---|
| 337 | if((isec1(6).eq.007).and.(isec1(7).eq.001)) then |
---|
| 338 | do jy=0,ny-1 |
---|
| 339 | do ix=0,nxfield-1 |
---|
| 340 | help=zsec4(nxfield*(ny-jy-1)+ix+1) |
---|
| 341 | if(ix.le.i180) then |
---|
| 342 | oro(i179+ix,jy)=help |
---|
| 343 | excessoro(i179+ix,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
| 344 | else |
---|
| 345 | oro(ix-i181,jy)=help |
---|
| 346 | excessoro(ix-i181,jy)=0.0 ! ISOBARIC SURFACES: SUBGRID TERRAIN DISREGARDED |
---|
| 347 | endif |
---|
| 348 | end do |
---|
| 349 | end do |
---|
| 350 | endif |
---|
| 351 | |
---|
| 352 | ! NCEP LAND SEA MASK |
---|
| 353 | !******************* |
---|
| 354 | |
---|
| 355 | if((isec1(6).eq.081).and.(isec1(7).eq.001)) then |
---|
| 356 | do jy=0,ny-1 |
---|
| 357 | do ix=0,nxfield-1 |
---|
| 358 | help=zsec4(nxfield*(ny-jy-1)+ix+1) |
---|
| 359 | if(ix.le.i180) then |
---|
| 360 | lsm(i179+ix,jy)=help |
---|
| 361 | else |
---|
| 362 | lsm(ix-i181,jy)=help |
---|
| 363 | endif |
---|
| 364 | end do |
---|
| 365 | end do |
---|
| 366 | endif |
---|
| 367 | |
---|
| 368 | call grib_release(igrib) |
---|
| 369 | |
---|
| 370 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
| 371 | ! |
---|
| 372 | ! CLOSING OF INPUT DATA FILE |
---|
| 373 | ! |
---|
| 374 | |
---|
| 375 | ! HSO |
---|
| 376 | 30 continue |
---|
| 377 | call grib_close_file(ifile) |
---|
| 378 | ! HSO end edits |
---|
| 379 | |
---|
| 380 | nuvz=iumax |
---|
| 381 | nwz =iumax |
---|
| 382 | nlev_ec=iumax |
---|
| 383 | |
---|
| 384 | if (nx.gt.nxmax) then |
---|
| 385 | write(*,*) 'FLEXPART error: Too many grid points in x direction.' |
---|
| 386 | write(*,*) 'Reduce resolution of wind fields.' |
---|
| 387 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
| 388 | write(*,*) nx,nxmax |
---|
| 389 | stop |
---|
| 390 | endif |
---|
| 391 | |
---|
| 392 | if (ny.gt.nymax) then |
---|
| 393 | write(*,*) 'FLEXPART error: Too many grid points in y direction.' |
---|
| 394 | write(*,*) 'Reduce resolution of wind fields.' |
---|
| 395 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
| 396 | write(*,*) ny,nymax |
---|
| 397 | stop |
---|
| 398 | endif |
---|
| 399 | |
---|
| 400 | if (nuvz.gt.nuvzmax) then |
---|
| 401 | write(*,*) 'FLEXPART error: Too many u,v grid points in z '// & |
---|
| 402 | 'direction.' |
---|
| 403 | write(*,*) 'Reduce resolution of wind fields.' |
---|
| 404 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
| 405 | write(*,*) nuvz,nuvzmax |
---|
| 406 | stop |
---|
| 407 | endif |
---|
| 408 | |
---|
| 409 | if (nwz.gt.nwzmax) then |
---|
| 410 | write(*,*) 'FLEXPART error: Too many w grid points in z '// & |
---|
| 411 | 'direction.' |
---|
| 412 | write(*,*) 'Reduce resolution of wind fields.' |
---|
| 413 | write(*,*) 'Or change parameter settings in file par_mod.' |
---|
| 414 | write(*,*) nwz,nwzmax |
---|
| 415 | stop |
---|
| 416 | endif |
---|
| 417 | |
---|
| 418 | ! If desired, shift all grids by nxshift grid cells |
---|
| 419 | !************************************************** |
---|
| 420 | |
---|
| 421 | if (xglobal) then |
---|
| 422 | call shift_field_0(oro,nxfield,ny) |
---|
| 423 | call shift_field_0(lsm,nxfield,ny) |
---|
| 424 | call shift_field_0(excessoro,nxfield,ny) |
---|
| 425 | endif |
---|
| 426 | |
---|
| 427 | ! Output of grid info |
---|
| 428 | !******************** |
---|
| 429 | |
---|
[4c64400] | 430 | if (lroot) then |
---|
| 431 | write(*,*) |
---|
| 432 | write(*,*) |
---|
| 433 | write(*,'(a,2i7)') 'Vertical levels in NCEP data: ', & |
---|
| 434 | nuvz,nwz |
---|
| 435 | write(*,*) |
---|
| 436 | write(*,'(a)') 'Mother domain:' |
---|
| 437 | write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Longitude range: ', & |
---|
| 438 | xlon0,' to ',xlon0+(nx-1)*dx,' Grid distance: ',dx |
---|
| 439 | write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Latitude range : ', & |
---|
| 440 | ylat0,' to ',ylat0+(ny-1)*dy,' Grid distance: ',dy |
---|
| 441 | write(*,*) |
---|
| 442 | end if |
---|
[e200b7a] | 443 | |
---|
| 444 | ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL |
---|
| 445 | ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM |
---|
| 446 | |
---|
| 447 | numskip=nlev_ec-nuvz ! number of ecmwf model layers not used |
---|
| 448 | ! by trajectory model |
---|
| 449 | do i=1,nwz |
---|
| 450 | j=numskip+i |
---|
| 451 | k=nlev_ec+1+numskip+i |
---|
| 452 | akm_usort(nwz-i+1)=pres(nwz-i+1) |
---|
| 453 | bkm(nwz-i+1)=0.0 |
---|
| 454 | end do |
---|
| 455 | |
---|
| 456 | !****************************** |
---|
| 457 | ! change Sabine Eckhardt: akm should always be in descending order ... readwind adapted! |
---|
| 458 | !****************************** |
---|
| 459 | do i=1,nwz |
---|
| 460 | if (akm_usort(1).gt.akm_usort(2)) then |
---|
| 461 | akm(i)=akm_usort(i) |
---|
| 462 | else |
---|
| 463 | akm(i)=akm_usort(nwz-i+1) |
---|
| 464 | endif |
---|
| 465 | end do |
---|
| 466 | |
---|
| 467 | ! |
---|
| 468 | ! CALCULATION OF AKZ, BKZ |
---|
| 469 | ! AKZ,BKZ: model discretization parameters at the center of each model |
---|
| 470 | ! layer |
---|
| 471 | ! |
---|
| 472 | ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0, |
---|
| 473 | ! i.e. ground level |
---|
| 474 | !***************************************************************************** |
---|
| 475 | |
---|
| 476 | do i=1,nuvz |
---|
| 477 | akz(i)=akm(i) |
---|
| 478 | bkz(i)=bkm(i) |
---|
| 479 | end do |
---|
| 480 | |
---|
| 481 | ! NOTE: In FLEXPART versions up to 4.0, the number of model levels was doubled |
---|
| 482 | ! upon the transformation to z levels. In order to save computer memory, this is |
---|
| 483 | ! not done anymore in the standard version. However, this option can still be |
---|
| 484 | ! switched on by replacing the following lines with those below, that are |
---|
| 485 | ! currently commented out. For this, similar changes are necessary in |
---|
| 486 | ! verttransform.f and verttranform_nests.f |
---|
| 487 | !***************************************************************************** |
---|
| 488 | |
---|
| 489 | nz=nuvz |
---|
| 490 | if (nz.gt.nzmax) stop 'nzmax too small' |
---|
| 491 | do i=1,nuvz |
---|
| 492 | aknew(i)=akz(i) |
---|
| 493 | bknew(i)=bkz(i) |
---|
| 494 | end do |
---|
| 495 | |
---|
| 496 | ! Switch on following lines to use doubled vertical resolution |
---|
| 497 | !************************************************************* |
---|
| 498 | !nz=nuvz+nwz-1 |
---|
| 499 | !if (nz.gt.nzmax) stop 'nzmax too small' |
---|
| 500 | !do 100 i=1,nwz |
---|
| 501 | ! aknew(2*(i-1)+1)=akm(i) |
---|
| 502 | !00 bknew(2*(i-1)+1)=bkm(i) |
---|
| 503 | !do 110 i=2,nuvz |
---|
| 504 | ! aknew(2*(i-1))=akz(i) |
---|
| 505 | !10 bknew(2*(i-1))=bkz(i) |
---|
| 506 | ! End doubled vertical resolution |
---|
| 507 | |
---|
| 508 | |
---|
| 509 | ! Determine the uppermost level for which the convection scheme shall be applied |
---|
| 510 | ! by assuming that there is no convection above 50 hPa (for standard SLP) |
---|
| 511 | !***************************************************************************** |
---|
| 512 | |
---|
| 513 | do i=1,nuvz-2 |
---|
| 514 | pint=akz(i)+bkz(i)*101325. |
---|
| 515 | if (pint.lt.5000.) goto 96 |
---|
| 516 | end do |
---|
| 517 | 96 nconvlev=i |
---|
| 518 | if (nconvlev.gt.nconvlevmax-1) then |
---|
| 519 | nconvlev=nconvlevmax-1 |
---|
| 520 | write(*,*) 'Attention, convection only calculated up to ', & |
---|
| 521 | akz(nconvlev)+bkz(nconvlev)*1013.25,' hPa' |
---|
| 522 | endif |
---|
| 523 | |
---|
| 524 | return |
---|
| 525 | |
---|
| 526 | 999 write(*,*) |
---|
| 527 | write(*,*) ' ###########################################'// & |
---|
| 528 | '###### ' |
---|
| 529 | write(*,*) ' TRAJECTORY MODEL SUBROUTINE GRIDCHECK:' |
---|
| 530 | write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfname(ifn) |
---|
| 531 | write(*,*) ' ###########################################'// & |
---|
| 532 | '###### ' |
---|
| 533 | write(*,*) |
---|
| 534 | write(*,'(a)') '!!! PLEASE INSERT A NEW CD-ROM AND !!!' |
---|
| 535 | write(*,'(a)') '!!! PRESS ANY KEY TO CONTINUE... !!!' |
---|
| 536 | write(*,'(a)') '!!! ...OR TERMINATE FLEXPART PRESSING!!!' |
---|
| 537 | write(*,'(a)') '!!! THE "X" KEY... !!!' |
---|
| 538 | write(*,*) |
---|
| 539 | read(*,'(a)') opt |
---|
| 540 | if(opt.eq.'X') then |
---|
| 541 | stop |
---|
| 542 | else |
---|
| 543 | goto 5 |
---|
| 544 | endif |
---|
| 545 | |
---|
[6ecb30a] | 546 | end subroutine gridcheck_gfs |
---|