[e200b7a] | 1 | subroutine readwind(indj,n,uuh,vvh,wwh) |
---|
| 2 | |
---|
| 3 | !********************************************************************** |
---|
| 4 | ! * |
---|
| 5 | ! TRAJECTORY MODEL SUBROUTINE READWIND * |
---|
| 6 | ! * |
---|
| 7 | !********************************************************************** |
---|
| 8 | ! * |
---|
| 9 | ! AUTHOR: G. WOTAWA * |
---|
| 10 | ! DATE: 1997-08-05 * |
---|
| 11 | ! LAST UPDATE: 2000-10-17, Andreas Stohl * |
---|
| 12 | ! * |
---|
| 13 | !********************************************************************** |
---|
| 14 | ! Changes, Bernd C. Krueger, Feb. 2001: |
---|
| 15 | ! Variables tth and qvh (on eta coordinates) in common block |
---|
| 16 | !********************************************************************** |
---|
| 17 | ! * |
---|
| 18 | ! DESCRIPTION: * |
---|
| 19 | ! * |
---|
| 20 | ! READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE * |
---|
| 21 | ! INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE * |
---|
| 22 | ! * |
---|
| 23 | ! INPUT: * |
---|
| 24 | ! indj indicates number of the wind field to be read in * |
---|
| 25 | ! n temporal index for meteorological fields (1 to 3)* |
---|
| 26 | ! * |
---|
| 27 | ! IMPORTANT VARIABLES FROM COMMON BLOCK: * |
---|
| 28 | ! * |
---|
| 29 | ! wfname File name of data to be read in * |
---|
| 30 | ! nx,ny,nuvz,nwz expected field dimensions * |
---|
| 31 | ! nlev_ec number of vertical levels ecmwf model * |
---|
| 32 | ! uu,vv,ww wind fields * |
---|
| 33 | ! tt,qv temperature and specific humidity * |
---|
| 34 | ! ps surface pressure * |
---|
| 35 | ! * |
---|
| 36 | !********************************************************************** |
---|
| 37 | |
---|
| 38 | use par_mod |
---|
| 39 | use com_mod |
---|
| 40 | |
---|
| 41 | implicit none |
---|
| 42 | |
---|
| 43 | real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 44 | real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
| 45 | real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
| 46 | integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,lunit |
---|
| 47 | |
---|
| 48 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
| 49 | |
---|
| 50 | ! dimension of isec2 at least (22+n), where n is the number of parallels or |
---|
| 51 | ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid |
---|
| 52 | |
---|
| 53 | ! dimension of zsec2 at least (10+nn), where nn is the number of vertical |
---|
| 54 | ! coordinate parameters |
---|
| 55 | |
---|
| 56 | integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2) |
---|
| 57 | integer :: isec4(64),inbuff(jpack),ilen,ierr,iword |
---|
| 58 | !integer iswap |
---|
| 59 | real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp) |
---|
| 60 | real :: xaux,yaux,xaux0,yaux0 |
---|
| 61 | real,parameter :: eps=1.e-4 |
---|
| 62 | real :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1) |
---|
| 63 | real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1 |
---|
| 64 | |
---|
| 65 | character(len=1) :: yoper = 'D' |
---|
| 66 | logical :: hflswitch,strswitch |
---|
| 67 | |
---|
| 68 | hflswitch=.false. |
---|
| 69 | strswitch=.false. |
---|
| 70 | levdiff2=nlev_ec-nwz+1 |
---|
| 71 | iumax=0 |
---|
| 72 | iwmax=0 |
---|
| 73 | |
---|
| 74 | ! |
---|
| 75 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
| 76 | ! |
---|
| 77 | 5 call pbopen(lunit,path(3)(1:length(3))//wfname(indj),'r',ierr) |
---|
| 78 | if(ierr.lt.0) goto 999 |
---|
| 79 | |
---|
| 80 | ifield=0 |
---|
| 81 | 10 ifield=ifield+1 |
---|
| 82 | ! |
---|
| 83 | ! GET NEXT FIELDS |
---|
| 84 | ! |
---|
| 85 | call pbgrib(lunit,inbuff,jpack,ilen,ierr) |
---|
| 86 | if(ierr.eq.-1) goto 50 ! EOF DETECTED |
---|
| 87 | if(ierr.lt.-1) goto 888 ! ERROR DETECTED |
---|
| 88 | |
---|
| 89 | ierr=1 |
---|
| 90 | |
---|
| 91 | ! Check whether we are on a little endian or on a big endian computer |
---|
| 92 | !******************************************************************** |
---|
| 93 | |
---|
| 94 | !if (inbuff(1).eq.1112101447) then ! little endian, swap bytes |
---|
| 95 | ! iswap=1+ilen/4 |
---|
| 96 | ! call swap32(inbuff,iswap) |
---|
| 97 | !else if (inbuff(1).ne.1196575042) then ! big endian |
---|
| 98 | ! stop 'subroutine gridcheck: corrupt GRIB data' |
---|
| 99 | !endif |
---|
| 100 | |
---|
| 101 | |
---|
| 102 | call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, & |
---|
| 103 | zsec4,jpunp,inbuff,jpack,iword,yoper,ierr) |
---|
| 104 | if (ierr.ne.0) goto 888 ! ERROR DETECTED |
---|
| 105 | |
---|
| 106 | if(ifield.eq.1) then |
---|
| 107 | |
---|
| 108 | ! CHECK GRID SPECIFICATIONS |
---|
| 109 | |
---|
| 110 | if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT' |
---|
| 111 | if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT' |
---|
| 112 | if(isec2(12)/2-1.ne.nlev_ec) & |
---|
| 113 | stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT' |
---|
| 114 | xaux=real(isec2(5))/1000.+real(nxshift)*dx |
---|
| 115 | yaux=real(isec2(7))/1000. |
---|
| 116 | xaux0=xlon0 |
---|
| 117 | yaux0=ylat0 |
---|
| 118 | if(xaux.lt.0.) xaux=xaux+360. |
---|
| 119 | if(yaux.lt.0.) yaux=yaux+360. |
---|
| 120 | if(xaux0.lt.0.) xaux0=xaux0+360. |
---|
| 121 | if(yaux0.lt.0.) yaux0=yaux0+360. |
---|
| 122 | if(abs(xaux-xaux0).gt.eps) & |
---|
| 123 | stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT' |
---|
| 124 | if(abs(yaux-yaux0).gt.eps) & |
---|
| 125 | stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT' |
---|
| 126 | endif |
---|
| 127 | |
---|
| 128 | do j=0,nymin1 |
---|
| 129 | do i=0,nxfield-1 |
---|
| 130 | k=isec1(8) |
---|
| 131 | if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE |
---|
| 132 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 133 | if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY |
---|
| 134 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 135 | if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY |
---|
| 136 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 137 | if(isec1(6).eq.133) then !! SPEC. HUMIDITY |
---|
| 138 | qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 139 | if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) & |
---|
| 140 | qvh(i,j,nlev_ec-k+2,n) = 0. |
---|
| 141 | ! this is necessary because the gridded data may contain |
---|
| 142 | ! spurious negative values |
---|
| 143 | endif |
---|
| 144 | if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS. |
---|
| 145 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 146 | |
---|
| 147 | if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY |
---|
| 148 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 149 | if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH |
---|
| 150 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 151 | if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS. |
---|
| 152 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 153 | if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER |
---|
| 154 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 155 | if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY |
---|
| 156 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 157 | if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY |
---|
| 158 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 159 | if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE |
---|
| 160 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 161 | if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT |
---|
| 162 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 163 | if(isec1(6).eq.142) then !! LARGE SCALE PREC. |
---|
| 164 | lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 165 | if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0. |
---|
| 166 | endif |
---|
| 167 | if(isec1(6).eq.143) then !! CONVECTIVE PREC. |
---|
| 168 | convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 169 | if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0. |
---|
| 170 | endif |
---|
| 171 | if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX |
---|
| 172 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 173 | if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) & |
---|
| 174 | hflswitch=.true. ! Heat flux available |
---|
| 175 | if(isec1(6).eq.176) then !! SOLAR RADIATION |
---|
| 176 | ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 177 | if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0. |
---|
| 178 | endif |
---|
| 179 | if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS |
---|
| 180 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 181 | if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS |
---|
| 182 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 183 | if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. & |
---|
| 184 | (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true. ! stress available |
---|
| 185 | if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY |
---|
| 186 | zsec4(nxfield*(ny-j-1)+i+1)/ga |
---|
| 187 | if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY |
---|
| 188 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 189 | if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK |
---|
| 190 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
| 191 | if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1) |
---|
| 192 | if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1) |
---|
| 193 | end do |
---|
| 194 | end do |
---|
| 195 | |
---|
| 196 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
| 197 | ! |
---|
| 198 | ! CLOSING OF INPUT DATA FILE |
---|
| 199 | ! |
---|
| 200 | 50 call pbclose(lunit,ierr) !! FINNISHED READING / CLOSING GRIB FILE |
---|
| 201 | |
---|
| 202 | if(levdiff2.eq.0) then |
---|
| 203 | iwmax=nlev_ec+1 |
---|
| 204 | do i=0,nxmin1 |
---|
| 205 | do j=0,nymin1 |
---|
| 206 | wwh(i,j,nlev_ec+1)=0. |
---|
| 207 | end do |
---|
| 208 | end do |
---|
| 209 | endif |
---|
| 210 | |
---|
| 211 | |
---|
| 212 | ! For global fields, assign the leftmost data column also to the rightmost |
---|
| 213 | ! data column; if required, shift whole grid by nxshift grid points |
---|
| 214 | !************************************************************************* |
---|
| 215 | |
---|
| 216 | if (xglobal) then |
---|
| 217 | call shift_field_0(ewss,nxfield,ny) |
---|
| 218 | call shift_field_0(nsss,nxfield,ny) |
---|
| 219 | call shift_field_0(oro,nxfield,ny) |
---|
| 220 | call shift_field_0(excessoro,nxfield,ny) |
---|
| 221 | call shift_field_0(lsm,nxfield,ny) |
---|
| 222 | call shift_field(ps,nxfield,ny,1,1,2,n) |
---|
| 223 | call shift_field(sd,nxfield,ny,1,1,2,n) |
---|
| 224 | call shift_field(msl,nxfield,ny,1,1,2,n) |
---|
| 225 | call shift_field(tcc,nxfield,ny,1,1,2,n) |
---|
| 226 | call shift_field(u10,nxfield,ny,1,1,2,n) |
---|
| 227 | call shift_field(v10,nxfield,ny,1,1,2,n) |
---|
| 228 | call shift_field(tt2,nxfield,ny,1,1,2,n) |
---|
| 229 | call shift_field(td2,nxfield,ny,1,1,2,n) |
---|
| 230 | call shift_field(lsprec,nxfield,ny,1,1,2,n) |
---|
| 231 | call shift_field(convprec,nxfield,ny,1,1,2,n) |
---|
| 232 | call shift_field(sshf,nxfield,ny,1,1,2,n) |
---|
| 233 | call shift_field(ssr,nxfield,ny,1,1,2,n) |
---|
| 234 | call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
| 235 | call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
| 236 | call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
| 237 | call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
| 238 | call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1) |
---|
| 239 | endif |
---|
| 240 | |
---|
| 241 | do i=0,nxmin1 |
---|
| 242 | do j=0,nymin1 |
---|
| 243 | surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2) |
---|
| 244 | end do |
---|
| 245 | end do |
---|
| 246 | |
---|
| 247 | if ((.not.hflswitch).or.(.not.strswitch)) then |
---|
| 248 | write(*,*) 'WARNING: No flux data contained in GRIB file ', & |
---|
| 249 | wfname(indj) |
---|
| 250 | |
---|
| 251 | ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD |
---|
| 252 | ! As ECMWF has increased the model resolution, such that now the first model |
---|
| 253 | ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level |
---|
| 254 | ! (3rd model level in FLEXPART) for the profile method |
---|
| 255 | !*************************************************************************** |
---|
| 256 | |
---|
| 257 | do i=0,nxmin1 |
---|
| 258 | do j=0,nymin1 |
---|
| 259 | plev1=akz(3)+bkz(3)*ps(i,j,1,n) |
---|
| 260 | pmean=0.5*(ps(i,j,1,n)+plev1) |
---|
| 261 | tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n)) |
---|
| 262 | fu=-r_air*tv/ga/pmean |
---|
| 263 | hlev1=fu*(plev1-ps(i,j,1,n)) ! HEIGTH OF FIRST MODEL LAYER |
---|
| 264 | ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2) |
---|
| 265 | fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2) |
---|
| 266 | call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, & |
---|
| 267 | tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, & |
---|
| 268 | surfstr(i,j,1,n),sshf(i,j,1,n)) |
---|
| 269 | if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200. |
---|
| 270 | if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400. |
---|
| 271 | end do |
---|
| 272 | end do |
---|
| 273 | endif |
---|
| 274 | |
---|
| 275 | |
---|
| 276 | ! Assign 10 m wind to model level at eta=1.0 to have one additional model |
---|
| 277 | ! level at the ground |
---|
| 278 | ! Specific humidity is taken the same as at one level above |
---|
| 279 | ! Temperature is taken as 2 m temperature |
---|
| 280 | !************************************************************************** |
---|
| 281 | |
---|
| 282 | do i=0,nxmin1 |
---|
| 283 | do j=0,nymin1 |
---|
| 284 | uuh(i,j,1)=u10(i,j,1,n) |
---|
| 285 | vvh(i,j,1)=v10(i,j,1,n) |
---|
| 286 | qvh(i,j,1,n)=qvh(i,j,2,n) |
---|
| 287 | tth(i,j,1,n)=tt2(i,j,1,n) |
---|
| 288 | end do |
---|
| 289 | end do |
---|
| 290 | |
---|
| 291 | if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT' |
---|
| 292 | if(iwmax.ne.nwz) stop 'READWIND: NWZ NOT CONSISTENT' |
---|
| 293 | |
---|
| 294 | return |
---|
| 295 | 888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
| 296 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
| 297 | write(*,*) ' #### IS NOT GRIB FORMAT !!! #### ' |
---|
| 298 | stop 'Execution terminated' |
---|
| 299 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
| 300 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
| 301 | write(*,*) ' #### CANNOT BE OPENED !!! #### ' |
---|
| 302 | stop 'Execution terminated' |
---|
| 303 | |
---|
| 304 | end subroutine readwind |
---|