[16] | 1 | !*********************************************************************** |
---|
| 2 | !* Copyright 2012,2013 * |
---|
| 3 | !* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, * |
---|
| 4 | !* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,* |
---|
| 5 | !* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, * |
---|
| 6 | !* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, * |
---|
| 7 | !* * |
---|
| 8 | !* This file is part of FLEXPART WRF * |
---|
| 9 | !* * |
---|
| 10 | !* FLEXPART is free software: you can redistribute it and/or modify * |
---|
| 11 | !* it under the terms of the GNU General Public License as published by* |
---|
| 12 | !* the Free Software Foundation, either version 3 of the License, or * |
---|
| 13 | !* (at your option) any later version. * |
---|
| 14 | !* * |
---|
| 15 | !* FLEXPART is distributed in the hope that it will be useful, * |
---|
| 16 | !* but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
| 17 | !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
| 18 | !* GNU General Public License for more details. * |
---|
| 19 | !* * |
---|
| 20 | !* You should have received a copy of the GNU General Public License * |
---|
| 21 | !* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
| 22 | !*********************************************************************** |
---|
| 23 | |
---|
| 24 | subroutine boundcond_domainfill(itime,loutend) |
---|
| 25 | ! i i |
---|
| 26 | !******************************************************************************* |
---|
| 27 | ! * |
---|
| 28 | ! Note: This is the FLEXPART_WRF version of subr. boundcond_domainfill. * |
---|
| 29 | ! The computational grid is the WRF x-y grid rather than lat-lon. * |
---|
| 30 | ! * |
---|
| 31 | ! Particles are created by this subroutine continuously throughout the * |
---|
| 32 | ! simulation at the boundaries of the domain-filling box. * |
---|
| 33 | ! All particles carry the same amount of mass which alltogether comprises the * |
---|
| 34 | ! mass of air within the box, which remains (more or less) constant. * |
---|
| 35 | ! * |
---|
| 36 | ! Author: A. Stohl * |
---|
| 37 | ! * |
---|
| 38 | ! 16 October 2002 * |
---|
| 39 | ! * |
---|
| 40 | ! 26 Oct 2005, R. Easter - changes to calc. of boundarea * |
---|
| 41 | ! associated with WRF horizontal grid. * |
---|
| 42 | ! Also need to get true ylat for pv calcs. * |
---|
| 43 | ! 11 Nov 2005, R. Easter - fixed error involving xy to latlong * |
---|
| 44 | ! 2012, J. Brioude: coded in fortran 90 * |
---|
| 45 | !******************************************************************************* |
---|
| 46 | ! * |
---|
| 47 | ! Variables: * |
---|
| 48 | ! * |
---|
| 49 | ! nx_we(2) grid indices for western and eastern boundary of domain- * |
---|
| 50 | ! filling trajectory calculations * |
---|
| 51 | ! ny_sn(2) grid indices for southern and northern boundary of domain- * |
---|
| 52 | ! filling trajectory calculations * |
---|
| 53 | ! * |
---|
| 54 | !******************************************************************************* |
---|
| 55 | use point_mod |
---|
| 56 | use par_mod |
---|
| 57 | use com_mod |
---|
| 58 | |
---|
| 59 | implicit none |
---|
| 60 | |
---|
| 61 | real :: dz,dz1,dz2,ran1,dt1,dt2,dtt,xm,cosfact,accmasst |
---|
| 62 | integer :: itime,in,indz,indzp,i,loutend |
---|
| 63 | integer :: j,k,ix,jy,m,indzh,indexh,minpart,ipart,mmass |
---|
| 64 | integer :: numactiveparticles |
---|
| 65 | |
---|
| 66 | real :: windl(2),rhol(2),dumx,dumy,xlon,ylat |
---|
| 67 | real :: windhl(2),rhohl(2) |
---|
| 68 | real :: windx,rhox |
---|
| 69 | real :: deltaz,boundarea,fluxofmass |
---|
| 70 | |
---|
| 71 | integer :: ixm,ixp,jym,jyp,indzm,mm |
---|
| 72 | real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2),yh1(2) |
---|
| 73 | |
---|
| 74 | integer :: idummy = -11 |
---|
| 75 | |
---|
| 76 | |
---|
| 77 | |
---|
| 78 | ! If domain-filling is global, no boundary conditions are needed |
---|
| 79 | !*************************************************************** |
---|
| 80 | |
---|
| 81 | if (gdomainfill) return |
---|
| 82 | |
---|
| 83 | accmasst=0. |
---|
| 84 | numactiveparticles=0 |
---|
| 85 | |
---|
| 86 | ! Terminate trajectories that have left the domain, if domain-filling |
---|
| 87 | ! trajectory calculation domain is not global |
---|
| 88 | !******************************************************************** |
---|
| 89 | |
---|
| 90 | do i=1,numpart |
---|
| 91 | if (itra1(i).eq.itime) then |
---|
| 92 | if ((ytra1(i).gt.real(ny_sn(2))).or. & |
---|
| 93 | (ytra1(i).lt.real(ny_sn(1)))) itra1(i)=-999999999 |
---|
| 94 | if (((.not.xglobal).or.(nx_we(2).ne.(nx-2))).and. & |
---|
| 95 | ((xtra1(i).lt.real(nx_we(1))).or. & |
---|
| 96 | (xtra1(i).gt.real(nx_we(2))))) itra1(i)=-999999999 |
---|
| 97 | endif |
---|
| 98 | if (itra1(i).ne.-999999999) numactiveparticles= & |
---|
| 99 | + numactiveparticles+1 |
---|
| 100 | enddo |
---|
| 101 | |
---|
| 102 | |
---|
| 103 | ! Determine auxiliary variables for time interpolation |
---|
| 104 | !***************************************************** |
---|
| 105 | |
---|
| 106 | dt1=real(itime-memtime(1)) |
---|
| 107 | dt2=real(memtime(2)-itime) |
---|
| 108 | dtt=1./(dt1+dt2) |
---|
| 109 | |
---|
| 110 | ! Initialize auxiliary variable used to search for vacant storage space |
---|
| 111 | !********************************************************************** |
---|
| 112 | |
---|
| 113 | minpart=1 |
---|
| 114 | |
---|
| 115 | !*************************************** |
---|
| 116 | ! Western and eastern boundary condition |
---|
| 117 | !*************************************** |
---|
| 118 | |
---|
| 119 | ! Loop from south to north |
---|
| 120 | !************************* |
---|
| 121 | |
---|
| 122 | do jy=ny_sn(1),ny_sn(2) |
---|
| 123 | |
---|
| 124 | ! Loop over western (index 1) and eastern (index 2) boundary |
---|
| 125 | !*********************************************************** |
---|
| 126 | |
---|
| 127 | do k=1,2 |
---|
| 128 | |
---|
| 129 | ! for FLEXPART_WRF, x & y coords are in meters. |
---|
| 130 | ! In the "do 70" loop, ylat is only needed for for pv calcs, |
---|
| 131 | ! "if (ylat.lt.0.) pvpart=-1.*pvpart" |
---|
| 132 | ! Note: in the FLEXPART_ECMWF code, ylat was not defined |
---|
| 133 | ! in the "do 70" loop (a bug). |
---|
| 134 | dumx=real(nx_we(k)) |
---|
| 135 | dumy=real(jy) |
---|
| 136 | ! Are these dumx,dumy correct ??? |
---|
| 137 | call xyindex_to_ll_wrf( 0, dumx, dumy, xlon, ylat ) |
---|
| 138 | |
---|
| 139 | ! Loop over all release locations in a column |
---|
| 140 | !******************************************** |
---|
| 141 | |
---|
| 142 | do j=1,numcolumn_we(k,jy) |
---|
| 143 | |
---|
| 144 | ! Determine, for each release location, the area of the corresponding boundary |
---|
| 145 | !***************************************************************************** |
---|
| 146 | |
---|
| 147 | if (j.eq.1) then |
---|
| 148 | deltaz=(zcolumn_we(k,jy,2)+zcolumn_we(k,jy,1))/2. |
---|
| 149 | else if (j.eq.numcolumn_we(k,jy)) then |
---|
| 150 | ! deltaz=height(nz)-(zcolumn_we(k,jy,j-1)+ |
---|
| 151 | ! + zcolumn_we(k,jy,j))/2. |
---|
| 152 | ! In order to avoid taking a very high column for very many particles, |
---|
| 153 | ! use the deltaz from one particle below instead |
---|
| 154 | deltaz=(zcolumn_we(k,jy,j)-zcolumn_we(k,jy,j-2))/2. |
---|
| 155 | else |
---|
| 156 | deltaz=(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))/2. |
---|
| 157 | endif |
---|
| 158 | |
---|
| 159 | ! for FLEXPART_ECMWF, dy is in degrees-lat, and 111198.5 converts |
---|
| 160 | ! from degrees-latitude to m |
---|
| 161 | ! for FLEXPART_WRF, dy is in meters |
---|
| 162 | ! if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then |
---|
| 163 | ! boundarea=deltaz*111198.5/2.*dy |
---|
| 164 | ! else |
---|
| 165 | ! boundarea=deltaz*111198.5*dy |
---|
| 166 | ! endif |
---|
| 167 | if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then |
---|
| 168 | boundarea=deltaz/2.*dy |
---|
| 169 | else |
---|
| 170 | boundarea=deltaz*dy |
---|
| 171 | endif |
---|
| 172 | |
---|
| 173 | |
---|
| 174 | ! Interpolate the wind velocity and density to the release location |
---|
| 175 | !****************************************************************** |
---|
| 176 | |
---|
| 177 | ! Determine the model level below the release position |
---|
| 178 | !***************************************************** |
---|
| 179 | |
---|
| 180 | do i=2,nz |
---|
| 181 | if (height(i).gt.zcolumn_we(k,jy,j)) then |
---|
| 182 | indz=i-1 |
---|
| 183 | indzp=i |
---|
| 184 | goto 6 |
---|
| 185 | endif |
---|
| 186 | enddo |
---|
| 187 | 6 continue |
---|
| 188 | |
---|
| 189 | ! Vertical distance to the level below and above current position |
---|
| 190 | !**************************************************************** |
---|
| 191 | |
---|
| 192 | dz1=zcolumn_we(k,jy,j)-height(indz) |
---|
| 193 | dz2=height(indzp)-zcolumn_we(k,jy,j) |
---|
| 194 | dz=1./(dz1+dz2) |
---|
| 195 | |
---|
| 196 | ! Vertical and temporal interpolation |
---|
| 197 | !************************************ |
---|
| 198 | |
---|
| 199 | do m=1,2 |
---|
| 200 | indexh=memind(m) |
---|
| 201 | do in=1,2 |
---|
| 202 | indzh=indz+in-1 |
---|
| 203 | windl(in)=uu(nx_we(k),jy,indzh,indexh) |
---|
| 204 | rhol(in)=rho(nx_we(k),jy,indzh,indexh) |
---|
| 205 | enddo |
---|
| 206 | |
---|
| 207 | windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz |
---|
| 208 | rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz |
---|
| 209 | enddo |
---|
| 210 | |
---|
| 211 | windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt |
---|
| 212 | rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt |
---|
| 213 | |
---|
| 214 | ! Calculate mass flux |
---|
| 215 | !******************** |
---|
| 216 | |
---|
| 217 | fluxofmass=windx*rhox*boundarea*real(lsynctime) |
---|
| 218 | |
---|
| 219 | |
---|
| 220 | ! If the mass flux is directed into the domain, add it to previous mass fluxes; |
---|
| 221 | ! if it is out of the domain, set accumulated mass flux to zero |
---|
| 222 | !****************************************************************************** |
---|
| 223 | |
---|
| 224 | if (k.eq.1) then |
---|
| 225 | if (fluxofmass.ge.0.) then |
---|
| 226 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+fluxofmass |
---|
| 227 | else |
---|
| 228 | acc_mass_we(k,jy,j)=0. |
---|
| 229 | endif |
---|
| 230 | else |
---|
| 231 | if (fluxofmass.le.0.) then |
---|
| 232 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+abs(fluxofmass) |
---|
| 233 | else |
---|
| 234 | acc_mass_we(k,jy,j)=0. |
---|
| 235 | endif |
---|
| 236 | endif |
---|
| 237 | accmasst=accmasst+acc_mass_we(k,jy,j) |
---|
| 238 | |
---|
| 239 | ! If the accumulated mass exceeds half the mass that each particle shall carry, |
---|
| 240 | ! one (or more) particle(s) is (are) released and the accumulated mass is |
---|
| 241 | ! reduced by the mass of this (these) particle(s) |
---|
| 242 | !****************************************************************************** |
---|
| 243 | |
---|
| 244 | if (acc_mass_we(k,jy,j).ge.xmassperparticle/2.) then |
---|
| 245 | mmass=int((acc_mass_we(k,jy,j)+xmassperparticle/2.)/ & |
---|
| 246 | xmassperparticle) |
---|
| 247 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)- & |
---|
| 248 | real(mmass)*xmassperparticle |
---|
| 249 | else |
---|
| 250 | mmass=0 |
---|
| 251 | endif |
---|
| 252 | |
---|
| 253 | do m=1,mmass |
---|
| 254 | do ipart=minpart,maxpart |
---|
| 255 | |
---|
| 256 | ! If a vacant storage space is found, attribute everything to this array element |
---|
| 257 | !******************************************************************************* |
---|
| 258 | |
---|
| 259 | if (itra1(ipart).ne.itime) then |
---|
| 260 | |
---|
| 261 | ! Assign particle positions |
---|
| 262 | !************************** |
---|
| 263 | |
---|
| 264 | xtra1(ipart)=real(nx_we(k)) |
---|
| 265 | if (jy.eq.ny_sn(1)) then |
---|
| 266 | ytra1(ipart)=real(jy)+0.5*ran1(idummy) |
---|
| 267 | else if (jy.eq.ny_sn(2)) then |
---|
| 268 | ytra1(ipart)=real(jy)-0.5*ran1(idummy) |
---|
| 269 | else |
---|
| 270 | ytra1(ipart)=real(jy)+(ran1(idummy)-.5) |
---|
| 271 | endif |
---|
| 272 | if (j.eq.1) then |
---|
| 273 | ztra1(ipart)=zcolumn_we(k,jy,1)+(zcolumn_we(k,jy,2)- & |
---|
| 274 | zcolumn_we(k,jy,1))/4. |
---|
| 275 | else if (j.eq.numcolumn_we(k,jy)) then |
---|
| 276 | ztra1(ipart)=(2.*zcolumn_we(k,jy,j)+ & |
---|
| 277 | zcolumn_we(k,jy,j-1)+height(nz))/4. |
---|
| 278 | else |
---|
| 279 | ztra1(ipart)=zcolumn_we(k,jy,j-1)+ran1(idummy)* & |
---|
| 280 | (zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1)) |
---|
| 281 | endif |
---|
| 282 | |
---|
| 283 | ! Interpolate PV to the particle position |
---|
| 284 | !**************************************** |
---|
| 285 | ixm=int(xtra1(ipart)) |
---|
| 286 | jym=int(ytra1(ipart)) |
---|
| 287 | ixp=ixm+1 |
---|
| 288 | jyp=jym+1 |
---|
| 289 | ddx=xtra1(ipart)-real(ixm) |
---|
| 290 | ddy=ytra1(ipart)-real(jym) |
---|
| 291 | rddx=1.-ddx |
---|
| 292 | rddy=1.-ddy |
---|
| 293 | p1=rddx*rddy |
---|
| 294 | p2=ddx*rddy |
---|
| 295 | p3=rddx*ddy |
---|
| 296 | p4=ddx*ddy |
---|
| 297 | do i=2,nz |
---|
| 298 | if (height(i).gt.ztra1(ipart)) then |
---|
| 299 | indzm=i-1 |
---|
| 300 | indzp=i |
---|
| 301 | goto 26 |
---|
| 302 | endif |
---|
| 303 | enddo |
---|
| 304 | 26 continue |
---|
| 305 | dz1=ztra1(ipart)-height(indzm) |
---|
| 306 | dz2=height(indzp)-ztra1(ipart) |
---|
| 307 | dz=1./(dz1+dz2) |
---|
| 308 | do mm=1,2 |
---|
| 309 | indexh=memind(mm) |
---|
| 310 | do in=1,2 |
---|
| 311 | indzh=indzm+in-1 |
---|
| 312 | y1(in)=p1*pv(ixm,jym,indzh,indexh) & |
---|
| 313 | +p2*pv(ixp,jym,indzh,indexh) & |
---|
| 314 | +p3*pv(ixm,jyp,indzh,indexh) & |
---|
| 315 | +p4*pv(ixp,jyp,indzh,indexh) |
---|
| 316 | enddo |
---|
| 317 | yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz |
---|
| 318 | enddo |
---|
| 319 | pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt |
---|
| 320 | !JB |
---|
| 321 | ! ylat=ylat0+ytra1(ipart)*dy |
---|
| 322 | |
---|
| 323 | if (ylat.lt.0.) pvpart=-1.*pvpart |
---|
| 324 | |
---|
| 325 | |
---|
| 326 | ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere |
---|
| 327 | !************************************************************************************* |
---|
| 328 | |
---|
| 329 | if (((ztra1(ipart).gt.3000.).and. & |
---|
| 330 | (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then |
---|
| 331 | ! if (((ztra1(ipart).lt.8000.) |
---|
| 332 | ! + ).or.(mdomainfill.eq.1)) then |
---|
| 333 | nclass(ipart)=min(int(ran1(idummy)* & |
---|
| 334 | real(nclassunc))+1,nclassunc) |
---|
| 335 | numactiveparticles=numactiveparticles+1 |
---|
| 336 | numparticlecount=numparticlecount+1 |
---|
| 337 | npoint(ipart)=numparticlecount |
---|
| 338 | idt(ipart)=mintime |
---|
| 339 | itra1(ipart)=itime |
---|
| 340 | itramem(ipart)=itra1(ipart) |
---|
| 341 | itrasplit(ipart)=itra1(ipart)+ldirect*itsplit |
---|
| 342 | xmass1(ipart,1)=xmassperparticle |
---|
| 343 | if (mdomainfill.eq.2) xmass1(ipart,1)= & |
---|
| 344 | xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9 |
---|
| 345 | ! + xmass1(ipart,1)*60.*48./29./10.**9 |
---|
| 346 | else |
---|
| 347 | goto 71 |
---|
| 348 | endif |
---|
| 349 | |
---|
| 350 | |
---|
| 351 | ! Increase numpart, if necessary |
---|
| 352 | !******************************* |
---|
| 353 | |
---|
| 354 | numpart=max(numpart,ipart) |
---|
| 355 | goto 73 ! Storage space has been found, stop searching |
---|
| 356 | endif |
---|
| 357 | enddo |
---|
| 358 | if (ipart.gt.maxpart) & |
---|
| 359 | stop 'boundcond_domainfill.f: too many particles required' |
---|
| 360 | 73 minpart=ipart+1 |
---|
| 361 | 71 continue |
---|
| 362 | |
---|
| 363 | enddo |
---|
| 364 | |
---|
| 365 | enddo |
---|
| 366 | enddo |
---|
| 367 | enddo |
---|
| 368 | |
---|
| 369 | !***************************************** |
---|
| 370 | ! Southern and northern boundary condition |
---|
| 371 | !***************************************** |
---|
| 372 | |
---|
| 373 | ! Loop from west to east |
---|
| 374 | !*********************** |
---|
| 375 | |
---|
| 376 | do ix=nx_we(1),nx_we(2) |
---|
| 377 | |
---|
| 378 | ! Loop over southern (index 1) and northern (index 2) boundary |
---|
| 379 | !************************************************************* |
---|
| 380 | |
---|
| 381 | do k=1,2 |
---|
| 382 | |
---|
| 383 | ! for FLEXPART_WRF, x & y coords are in meters. |
---|
| 384 | ! ylat=ylat0+real(ny_sn(k))*dy |
---|
| 385 | ! cosfact=cos(ylat*pi180) |
---|
| 386 | ! In the "do 170" loop, ylat is only needed for for pv calcs, |
---|
| 387 | ! "if (ylat.lt.0.) pvpart=-1.*pvpart" |
---|
| 388 | dumx=real(ix) |
---|
| 389 | dumy=real(ny_sn(k)) |
---|
| 390 | ! Are these dumx,dumy correct ??? |
---|
| 391 | call xyindex_to_ll_wrf( 0, dumx, dumy, xlon, ylat ) |
---|
| 392 | |
---|
| 393 | ! Loop over all release locations in a column |
---|
| 394 | !******************************************** |
---|
| 395 | |
---|
| 396 | do j=1,numcolumn_sn(k,ix) |
---|
| 397 | |
---|
| 398 | ! Determine, for each release location, the area of the corresponding boundary |
---|
| 399 | !***************************************************************************** |
---|
| 400 | |
---|
| 401 | if (j.eq.1) then |
---|
| 402 | deltaz=(zcolumn_sn(k,ix,2)+zcolumn_sn(k,ix,1))/2. |
---|
| 403 | else if (j.eq.numcolumn_sn(k,ix)) then |
---|
| 404 | ! deltaz=height(nz)-(zcolumn_sn(k,ix,j-1)+ |
---|
| 405 | ! + zcolumn_sn(k,ix,j))/2. |
---|
| 406 | ! In order to avoid taking a very high column for very many particles, |
---|
| 407 | ! use the deltaz from one particle below instead |
---|
| 408 | deltaz=(zcolumn_sn(k,ix,j)-zcolumn_sn(k,ix,j-2))/2. |
---|
| 409 | else |
---|
| 410 | deltaz=(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))/2. |
---|
| 411 | endif |
---|
| 412 | |
---|
| 413 | ! for FLEXPART_ECMWF, dx is in degrees-long, and 111198.5*cosfact converts |
---|
| 414 | ! from degrees-longitude to m |
---|
| 415 | ! for FLEXPART_WRF, dx is in meters |
---|
| 416 | ! if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then |
---|
| 417 | ! boundarea=deltaz*111198.5/2.*cosfact*dx |
---|
| 418 | ! else |
---|
| 419 | ! boundarea=deltaz*111198.5*cosfact*dx |
---|
| 420 | ! endif |
---|
| 421 | if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then |
---|
| 422 | boundarea=deltaz/2.*dx |
---|
| 423 | else |
---|
| 424 | boundarea=deltaz*dx |
---|
| 425 | endif |
---|
| 426 | |
---|
| 427 | |
---|
| 428 | ! Interpolate the wind velocity and density to the release location |
---|
| 429 | !****************************************************************** |
---|
| 430 | |
---|
| 431 | ! Determine the model level below the release position |
---|
| 432 | !***************************************************** |
---|
| 433 | |
---|
| 434 | do i=2,nz |
---|
| 435 | if (height(i).gt.zcolumn_sn(k,ix,j)) then |
---|
| 436 | indz=i-1 |
---|
| 437 | indzp=i |
---|
| 438 | goto 16 |
---|
| 439 | endif |
---|
| 440 | enddo |
---|
| 441 | 16 continue |
---|
| 442 | |
---|
| 443 | ! Vertical distance to the level below and above current position |
---|
| 444 | !**************************************************************** |
---|
| 445 | |
---|
| 446 | dz1=zcolumn_sn(k,ix,j)-height(indz) |
---|
| 447 | dz2=height(indzp)-zcolumn_sn(k,ix,j) |
---|
| 448 | dz=1./(dz1+dz2) |
---|
| 449 | |
---|
| 450 | ! Vertical and temporal interpolation |
---|
| 451 | !************************************ |
---|
| 452 | |
---|
| 453 | do m=1,2 |
---|
| 454 | indexh=memind(m) |
---|
| 455 | do in=1,2 |
---|
| 456 | indzh=indz+in-1 |
---|
| 457 | windl(in)=vv(ix,ny_sn(k),indzh,indexh) |
---|
| 458 | rhol(in)=rho(ix,ny_sn(k),indzh,indexh) |
---|
| 459 | end do |
---|
| 460 | |
---|
| 461 | windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz |
---|
| 462 | rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz |
---|
| 463 | end do |
---|
| 464 | |
---|
| 465 | windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt |
---|
| 466 | rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt |
---|
| 467 | |
---|
| 468 | ! Calculate mass flux |
---|
| 469 | !******************** |
---|
| 470 | |
---|
| 471 | fluxofmass=windx*rhox*boundarea*real(lsynctime) |
---|
| 472 | |
---|
| 473 | ! If the mass flux is directed into the domain, add it to previous mass fluxes; |
---|
| 474 | ! if it is out of the domain, set accumulated mass flux to zero |
---|
| 475 | !****************************************************************************** |
---|
| 476 | |
---|
| 477 | if (k.eq.1) then |
---|
| 478 | if (fluxofmass.ge.0.) then |
---|
| 479 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+fluxofmass |
---|
| 480 | else |
---|
| 481 | acc_mass_sn(k,ix,j)=0. |
---|
| 482 | endif |
---|
| 483 | else |
---|
| 484 | if (fluxofmass.le.0.) then |
---|
| 485 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+abs(fluxofmass) |
---|
| 486 | else |
---|
| 487 | acc_mass_sn(k,ix,j)=0. |
---|
| 488 | endif |
---|
| 489 | endif |
---|
| 490 | accmasst=accmasst+acc_mass_sn(k,ix,j) |
---|
| 491 | |
---|
| 492 | ! If the accumulated mass exceeds half the mass that each particle shall carry, |
---|
| 493 | ! one (or more) particle(s) is (are) released and the accumulated mass is |
---|
| 494 | ! reduced by the mass of this (these) particle(s) |
---|
| 495 | !****************************************************************************** |
---|
| 496 | |
---|
| 497 | if (acc_mass_sn(k,ix,j).ge.xmassperparticle/2.) then |
---|
| 498 | mmass=int((acc_mass_sn(k,ix,j)+xmassperparticle/2.)/ & |
---|
| 499 | xmassperparticle) |
---|
| 500 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)- & |
---|
| 501 | real(mmass)*xmassperparticle |
---|
| 502 | else |
---|
| 503 | mmass=0 |
---|
| 504 | endif |
---|
| 505 | do m=1,mmass |
---|
| 506 | do ipart=minpart,maxpart |
---|
| 507 | |
---|
| 508 | ! If a vacant storage space is found, attribute everything to this array element |
---|
| 509 | !******************************************************************************* |
---|
| 510 | |
---|
| 511 | if (itra1(ipart).ne.itime) then |
---|
| 512 | |
---|
| 513 | ! Assign particle positions |
---|
| 514 | !************************** |
---|
| 515 | |
---|
| 516 | ytra1(ipart)=real(ny_sn(k)) |
---|
| 517 | if (ix.eq.nx_we(1)) then |
---|
| 518 | xtra1(ipart)=real(ix)+0.5*ran1(idummy) |
---|
| 519 | else if (ix.eq.nx_we(2)) then |
---|
| 520 | xtra1(ipart)=real(ix)-0.5*ran1(idummy) |
---|
| 521 | else |
---|
| 522 | xtra1(ipart)=real(ix)+(ran1(idummy)-.5) |
---|
| 523 | endif |
---|
| 524 | if (j.eq.1) then |
---|
| 525 | ztra1(ipart)=zcolumn_sn(k,ix,1)+(zcolumn_sn(k,ix,2)- & |
---|
| 526 | zcolumn_sn(k,ix,1))/4. |
---|
| 527 | else if (j.eq.numcolumn_sn(k,ix)) then |
---|
| 528 | ztra1(ipart)=(2.*zcolumn_sn(k,ix,j)+ & |
---|
| 529 | zcolumn_sn(k,ix,j-1)+height(nz))/4. |
---|
| 530 | else |
---|
| 531 | ztra1(ipart)=zcolumn_sn(k,ix,j-1)+ran1(idummy)* & |
---|
| 532 | (zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1)) |
---|
| 533 | endif |
---|
| 534 | |
---|
| 535 | |
---|
| 536 | ! Interpolate PV to the particle position |
---|
| 537 | !**************************************** |
---|
| 538 | ixm=int(xtra1(ipart)) |
---|
| 539 | jym=int(ytra1(ipart)) |
---|
| 540 | ixp=ixm+1 |
---|
| 541 | jyp=jym+1 |
---|
| 542 | ddx=xtra1(ipart)-real(ixm) |
---|
| 543 | ddy=ytra1(ipart)-real(jym) |
---|
| 544 | rddx=1.-ddx |
---|
| 545 | rddy=1.-ddy |
---|
| 546 | p1=rddx*rddy |
---|
| 547 | p2=ddx*rddy |
---|
| 548 | p3=rddx*ddy |
---|
| 549 | p4=ddx*ddy |
---|
| 550 | do i=2,nz |
---|
| 551 | if (height(i).gt.ztra1(ipart)) then |
---|
| 552 | indzm=i-1 |
---|
| 553 | indzp=i |
---|
| 554 | goto 126 |
---|
| 555 | endif |
---|
| 556 | enddo |
---|
| 557 | 126 continue |
---|
| 558 | dz1=ztra1(ipart)-height(indzm) |
---|
| 559 | dz2=height(indzp)-ztra1(ipart) |
---|
| 560 | dz=1./(dz1+dz2) |
---|
| 561 | do mm=1,2 |
---|
| 562 | indexh=memind(mm) |
---|
| 563 | do in=1,2 |
---|
| 564 | indzh=indzm+in-1 |
---|
| 565 | y1(in)=p1*pv(ixm,jym,indzh,indexh) & |
---|
| 566 | +p2*pv(ixp,jym,indzh,indexh) & |
---|
| 567 | +p3*pv(ixm,jyp,indzh,indexh) & |
---|
| 568 | +p4*pv(ixp,jyp,indzh,indexh) |
---|
| 569 | end do |
---|
| 570 | yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz |
---|
| 571 | end do |
---|
| 572 | pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt |
---|
| 573 | if (ylat.lt.0.) pvpart=-1.*pvpart |
---|
| 574 | |
---|
| 575 | |
---|
| 576 | ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere |
---|
| 577 | !************************************************************************************* |
---|
| 578 | |
---|
| 579 | if (((ztra1(ipart).gt.3000.).and. & |
---|
| 580 | (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then |
---|
| 581 | nclass(ipart)=min(int(ran1(idummy)* & |
---|
| 582 | real(nclassunc))+1,nclassunc) |
---|
| 583 | numactiveparticles=numactiveparticles+1 |
---|
| 584 | numparticlecount=numparticlecount+1 |
---|
| 585 | npoint(ipart)=numparticlecount |
---|
| 586 | idt(ipart)=mintime |
---|
| 587 | itra1(ipart)=itime |
---|
| 588 | itramem(ipart)=itra1(ipart) |
---|
| 589 | itrasplit(ipart)=itra1(ipart)+ldirect*itsplit |
---|
| 590 | xmass1(ipart,1)=xmassperparticle |
---|
| 591 | if (mdomainfill.eq.2) xmass1(ipart,1)= & |
---|
| 592 | xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9 |
---|
| 593 | else |
---|
| 594 | goto 171 |
---|
| 595 | endif |
---|
| 596 | |
---|
| 597 | |
---|
| 598 | |
---|
| 599 | ! Increase numpart, if necessary |
---|
| 600 | !******************************* |
---|
| 601 | numpart=max(numpart,ipart) |
---|
| 602 | goto 173 ! Storage space has been found, stop searching |
---|
| 603 | endif |
---|
| 604 | enddo |
---|
| 605 | if (ipart.gt.maxpart) & |
---|
| 606 | stop 'boundcond_domainfill.f: too many particles required' |
---|
| 607 | 173 minpart=ipart+1 |
---|
| 608 | 171 continue |
---|
| 609 | enddo |
---|
| 610 | |
---|
| 611 | enddo |
---|
| 612 | enddo |
---|
| 613 | enddo |
---|
| 614 | |
---|
| 615 | xm=0. |
---|
| 616 | do i=1,numpart |
---|
| 617 | if (itra1(i).eq.itime) xm=xm+xmass1(i,1) |
---|
| 618 | end do |
---|
| 619 | |
---|
| 620 | !write(*,*) itime,numactiveparticles,numparticlecount,numpart, |
---|
| 621 | ! +xm,accmasst,xm+accmasst |
---|
| 622 | |
---|
| 623 | |
---|
| 624 | ! If particles shall be dumped, then accumulated masses at the domain boundaries |
---|
| 625 | ! must be dumped, too, to be used for later runs |
---|
| 626 | !***************************************************************************** |
---|
| 627 | |
---|
| 628 | if ((ipout.gt.0).and.(itime.eq.loutend)) then |
---|
| 629 | open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', & |
---|
| 630 | form='unformatted') |
---|
| 631 | write(unitboundcond) numcolumn_we,numcolumn_sn, & |
---|
| 632 | zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn |
---|
| 633 | close(unitboundcond) |
---|
| 634 | endif |
---|
| 635 | |
---|
| 636 | end subroutine boundcond_domainfill |
---|
| 637 | |
---|