[8a65cb0] | 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 | |
---|
| 22 | subroutine init_domainfill |
---|
| 23 | ! |
---|
| 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! Initializes particles equally distributed over the first release location * |
---|
| 27 | ! specified in file RELEASES. This box is assumed to be the domain for doing * |
---|
| 28 | ! domain-filling trajectory calculations. * |
---|
| 29 | ! All particles carry the same amount of mass which alltogether comprises the* |
---|
| 30 | ! mass of air within the box. * |
---|
| 31 | ! * |
---|
| 32 | ! Author: A. Stohl * |
---|
| 33 | ! * |
---|
| 34 | ! 15 October 2002 * |
---|
| 35 | ! * |
---|
| 36 | ! CHANGES * |
---|
| 37 | ! 12/2014 eso: MPI version * |
---|
| 38 | !***************************************************************************** |
---|
| 39 | ! * |
---|
| 40 | ! Variables: * |
---|
| 41 | ! * |
---|
| 42 | ! numparticlecount consecutively counts the number of particles released * |
---|
| 43 | ! nx_we(2) grid indices for western and eastern boundary of domain- * |
---|
| 44 | ! filling trajectory calculations * |
---|
| 45 | ! ny_sn(2) grid indices for southern and northern boundary of domain- * |
---|
| 46 | ! filling trajectory calculations * |
---|
| 47 | ! * |
---|
| 48 | !***************************************************************************** |
---|
| 49 | |
---|
| 50 | use point_mod |
---|
| 51 | use par_mod |
---|
| 52 | use com_mod |
---|
| 53 | use random_mod, only: ran1 |
---|
| 54 | use mpi_mod |
---|
| 55 | |
---|
| 56 | ! :TODO: parallelize this |
---|
| 57 | |
---|
| 58 | implicit none |
---|
| 59 | |
---|
| 60 | integer :: j,ix,jy,kz,ncolumn,numparttot |
---|
| 61 | real :: gridarea(0:nymax-1),pp(nzmax),ylat,ylatp,ylatm,hzone |
---|
| 62 | real :: cosfactm,cosfactp,deltacol,dz1,dz2,dz,pnew,fractus |
---|
| 63 | real,parameter :: pih=pi/180. |
---|
| 64 | real :: colmass(0:nxmax-1,0:nymax-1),colmasstotal,zposition |
---|
| 65 | |
---|
| 66 | integer :: ixm,ixp,jym,jyp,indzm,indzp,in,indzh,i,jj |
---|
| 67 | real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2) |
---|
| 68 | |
---|
| 69 | integer :: idummy = -11 |
---|
| 70 | logical :: first_call=.true. |
---|
| 71 | |
---|
| 72 | ! Use different seed for each process |
---|
| 73 | if (first_call) then |
---|
| 74 | idummy=idummy+mp_seed |
---|
| 75 | first_call=.false. |
---|
| 76 | end if |
---|
| 77 | |
---|
| 78 | ! Determine the release region (only full grid cells), over which particles |
---|
| 79 | ! shall be initialized |
---|
| 80 | ! Use 2 fields for west/east and south/north boundary |
---|
| 81 | !************************************************************************** |
---|
| 82 | |
---|
| 83 | nx_we(1)=max(int(xpoint1(1)),0) |
---|
| 84 | nx_we(2)=min((int(xpoint2(1))+1),nxmin1) |
---|
| 85 | ny_sn(1)=max(int(ypoint1(1)),0) |
---|
| 86 | ny_sn(2)=min((int(ypoint2(1))+1),nymin1) |
---|
| 87 | |
---|
| 88 | ! For global simulations (both global wind data and global domain-filling), |
---|
| 89 | ! set a switch, such that no boundary conditions are used |
---|
| 90 | !************************************************************************** |
---|
| 91 | if (xglobal.and.sglobal.and.nglobal) then |
---|
| 92 | if ((nx_we(1).eq.0).and.(nx_we(2).eq.nxmin1).and. & |
---|
| 93 | (ny_sn(1).eq.0).and.(ny_sn(2).eq.nymin1)) then |
---|
| 94 | gdomainfill=.true. |
---|
| 95 | else |
---|
| 96 | gdomainfill=.false. |
---|
| 97 | endif |
---|
| 98 | endif |
---|
| 99 | |
---|
| 100 | ! Do not release particles twice (i.e., not at both in the leftmost and rightmost |
---|
| 101 | ! grid cell) for a global domain |
---|
| 102 | !***************************************************************************** |
---|
| 103 | if (xglobal) nx_we(2)=min(nx_we(2),nx-2) |
---|
| 104 | |
---|
| 105 | |
---|
| 106 | ! Calculate area of grid cell with formula M=2*pi*R*h*dx/360, |
---|
| 107 | ! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90 |
---|
| 108 | !************************************************************ |
---|
| 109 | |
---|
| 110 | do jy=ny_sn(1),ny_sn(2) ! loop about latitudes |
---|
| 111 | ylat=ylat0+real(jy)*dy |
---|
| 112 | ylatp=ylat+0.5*dy |
---|
| 113 | ylatm=ylat-0.5*dy |
---|
| 114 | if ((ylatm.lt.0).and.(ylatp.gt.0.)) then |
---|
| 115 | hzone=1./dyconst |
---|
| 116 | else |
---|
| 117 | cosfactp=cos(ylatp*pih)*r_earth |
---|
| 118 | cosfactm=cos(ylatm*pih)*r_earth |
---|
| 119 | if (cosfactp.lt.cosfactm) then |
---|
| 120 | hzone=sqrt(r_earth**2-cosfactp**2)- & |
---|
| 121 | sqrt(r_earth**2-cosfactm**2) |
---|
| 122 | else |
---|
| 123 | hzone=sqrt(r_earth**2-cosfactm**2)- & |
---|
| 124 | sqrt(r_earth**2-cosfactp**2) |
---|
| 125 | endif |
---|
| 126 | endif |
---|
| 127 | gridarea(jy)=2.*pi*r_earth*hzone*dx/360. |
---|
| 128 | end do |
---|
| 129 | |
---|
| 130 | ! Do the same for the south pole |
---|
| 131 | |
---|
| 132 | if (sglobal) then |
---|
| 133 | ylat=ylat0 |
---|
| 134 | ylatp=ylat+0.5*dy |
---|
| 135 | ylatm=ylat |
---|
| 136 | cosfactm=0. |
---|
| 137 | cosfactp=cos(ylatp*pih)*r_earth |
---|
| 138 | hzone=sqrt(r_earth**2-cosfactm**2)- & |
---|
| 139 | sqrt(r_earth**2-cosfactp**2) |
---|
| 140 | gridarea(0)=2.*pi*r_earth*hzone*dx/360. |
---|
| 141 | endif |
---|
| 142 | |
---|
| 143 | ! Do the same for the north pole |
---|
| 144 | |
---|
| 145 | if (nglobal) then |
---|
| 146 | ylat=ylat0+real(nymin1)*dy |
---|
| 147 | ylatp=ylat |
---|
| 148 | ylatm=ylat-0.5*dy |
---|
| 149 | cosfactp=0. |
---|
| 150 | cosfactm=cos(ylatm*pih)*r_earth |
---|
| 151 | hzone=sqrt(r_earth**2-cosfactp**2)- & |
---|
| 152 | sqrt(r_earth**2-cosfactm**2) |
---|
| 153 | gridarea(nymin1)=2.*pi*r_earth*hzone*dx/360. |
---|
| 154 | endif |
---|
| 155 | |
---|
| 156 | |
---|
| 157 | ! Calculate total mass of each grid column and of the whole atmosphere |
---|
| 158 | !********************************************************************* |
---|
| 159 | |
---|
| 160 | colmasstotal=0. |
---|
| 161 | do jy=ny_sn(1),ny_sn(2) ! loop about latitudes |
---|
| 162 | do ix=nx_we(1),nx_we(2) ! loop about longitudes |
---|
| 163 | pp(1)=rho(ix,jy,1,1)*r_air*tt(ix,jy,1,1) |
---|
| 164 | pp(nz)=rho(ix,jy,nz,1)*r_air*tt(ix,jy,nz,1) |
---|
| 165 | colmass(ix,jy)=(pp(1)-pp(nz))/ga*gridarea(jy) |
---|
| 166 | colmasstotal=colmasstotal+colmass(ix,jy) |
---|
| 167 | end do |
---|
| 168 | end do |
---|
| 169 | |
---|
| 170 | write(*,*) 'Atm. mass: ',colmasstotal |
---|
| 171 | |
---|
| 172 | |
---|
| 173 | if (ipin.eq.0) numpart=0 |
---|
| 174 | |
---|
| 175 | ! Determine the particle positions |
---|
| 176 | !********************************* |
---|
| 177 | |
---|
| 178 | numparttot=0 |
---|
| 179 | numcolumn=0 |
---|
| 180 | do jy=ny_sn(1),ny_sn(2) ! loop about latitudes |
---|
| 181 | ylat=ylat0+real(jy)*dy |
---|
| 182 | do ix=nx_we(1),nx_we(2) ! loop about longitudes |
---|
| 183 | ncolumn=nint(0.999*real(npart(1))*colmass(ix,jy)/ & |
---|
| 184 | colmasstotal) |
---|
| 185 | if (ncolumn.eq.0) goto 30 |
---|
| 186 | if (ncolumn.gt.numcolumn) numcolumn=ncolumn |
---|
| 187 | |
---|
| 188 | ! Calculate pressure at the altitudes of model surfaces, using the air density |
---|
| 189 | ! information, which is stored as a 3-d field |
---|
| 190 | !***************************************************************************** |
---|
| 191 | |
---|
| 192 | do kz=1,nz |
---|
| 193 | pp(kz)=rho(ix,jy,kz,1)*r_air*tt(ix,jy,kz,1) |
---|
| 194 | end do |
---|
| 195 | |
---|
| 196 | |
---|
| 197 | deltacol=(pp(1)-pp(nz))/real(ncolumn) |
---|
| 198 | pnew=pp(1)+deltacol/2. |
---|
| 199 | jj=0 |
---|
| 200 | do j=1,ncolumn |
---|
| 201 | jj=jj+1 |
---|
| 202 | |
---|
| 203 | |
---|
| 204 | ! For columns with many particles (i.e. around the equator), distribute |
---|
| 205 | ! the particles equally, for columns with few particles (i.e. around the |
---|
| 206 | ! poles), distribute the particles randomly |
---|
| 207 | !*********************************************************************** |
---|
| 208 | |
---|
| 209 | |
---|
| 210 | if (ncolumn.gt.20) then |
---|
| 211 | pnew=pnew-deltacol |
---|
| 212 | else |
---|
| 213 | pnew=pp(1)-ran1(idummy)*(pp(1)-pp(nz)) |
---|
| 214 | endif |
---|
| 215 | |
---|
| 216 | do kz=1,nz-1 |
---|
| 217 | if ((pp(kz).ge.pnew).and.(pp(kz+1).lt.pnew)) then |
---|
| 218 | dz1=pp(kz)-pnew |
---|
| 219 | dz2=pnew-pp(kz+1) |
---|
| 220 | dz=1./(dz1+dz2) |
---|
| 221 | |
---|
| 222 | ! Assign particle position |
---|
| 223 | !************************* |
---|
| 224 | ! Do the following steps only if particles are not read in from previous model run |
---|
| 225 | !***************************************************************************** |
---|
| 226 | if (ipin.eq.0) then |
---|
| 227 | xtra1(numpart+jj)=real(ix)-0.5+ran1(idummy) |
---|
| 228 | if (ix.eq.0) xtra1(numpart+jj)=ran1(idummy) |
---|
| 229 | if (ix.eq.nxmin1) xtra1(numpart+jj)= & |
---|
| 230 | real(nxmin1)-ran1(idummy) |
---|
| 231 | ytra1(numpart+jj)=real(jy)-0.5+ran1(idummy) |
---|
| 232 | ztra1(numpart+jj)=(height(kz)*dz2+height(kz+1)*dz1)*dz |
---|
| 233 | if (ztra1(numpart+jj).gt.height(nz)-0.5) & |
---|
| 234 | ztra1(numpart+jj)=height(nz)-0.5 |
---|
| 235 | |
---|
| 236 | |
---|
| 237 | ! Interpolate PV to the particle position |
---|
| 238 | !**************************************** |
---|
| 239 | ixm=int(xtra1(numpart+jj)) |
---|
| 240 | jym=int(ytra1(numpart+jj)) |
---|
| 241 | ixp=ixm+1 |
---|
| 242 | jyp=jym+1 |
---|
| 243 | ddx=xtra1(numpart+jj)-real(ixm) |
---|
| 244 | ddy=ytra1(numpart+jj)-real(jym) |
---|
| 245 | rddx=1.-ddx |
---|
| 246 | rddy=1.-ddy |
---|
| 247 | p1=rddx*rddy |
---|
| 248 | p2=ddx*rddy |
---|
| 249 | p3=rddx*ddy |
---|
| 250 | p4=ddx*ddy |
---|
| 251 | do i=2,nz |
---|
| 252 | if (height(i).gt.ztra1(numpart+jj)) then |
---|
| 253 | indzm=i-1 |
---|
| 254 | indzp=i |
---|
| 255 | goto 6 |
---|
| 256 | endif |
---|
| 257 | end do |
---|
| 258 | 6 continue |
---|
| 259 | dz1=ztra1(numpart+jj)-height(indzm) |
---|
| 260 | dz2=height(indzp)-ztra1(numpart+jj) |
---|
| 261 | dz=1./(dz1+dz2) |
---|
| 262 | do in=1,2 |
---|
| 263 | indzh=indzm+in-1 |
---|
| 264 | y1(in)=p1*pv(ixm,jym,indzh,1) & |
---|
| 265 | +p2*pv(ixp,jym,indzh,1) & |
---|
| 266 | +p3*pv(ixm,jyp,indzh,1) & |
---|
| 267 | +p4*pv(ixp,jyp,indzh,1) |
---|
| 268 | end do |
---|
| 269 | pvpart=(dz2*y1(1)+dz1*y1(2))*dz |
---|
| 270 | if (ylat.lt.0.) pvpart=-1.*pvpart |
---|
| 271 | |
---|
| 272 | |
---|
| 273 | ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere |
---|
| 274 | !***************************************************************************** |
---|
| 275 | |
---|
| 276 | if (((ztra1(numpart+jj).gt.3000.).and. & |
---|
| 277 | (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then |
---|
| 278 | |
---|
| 279 | ! Assign certain properties to the particle |
---|
| 280 | !****************************************** |
---|
| 281 | nclass(numpart+jj)=min(int(ran1(idummy)* & |
---|
| 282 | real(nclassunc))+1,nclassunc) |
---|
| 283 | numparticlecount=numparticlecount+1 |
---|
| 284 | npoint(numpart+jj)=numparticlecount |
---|
| 285 | idt(numpart+jj)=mintime |
---|
| 286 | itra1(numpart+jj)=0 |
---|
| 287 | itramem(numpart+jj)=0 |
---|
| 288 | itrasplit(numpart+jj)=itra1(numpart+jj)+ldirect* & |
---|
| 289 | itsplit |
---|
| 290 | xmass1(numpart+jj,1)=colmass(ix,jy)/real(ncolumn) |
---|
| 291 | if (mdomainfill.eq.2) xmass1(numpart+jj,1)= & |
---|
| 292 | xmass1(numpart+jj,1)*pvpart*48./29.*ozonescale/10.**9 |
---|
| 293 | else |
---|
| 294 | jj=jj-1 |
---|
| 295 | endif |
---|
| 296 | endif |
---|
| 297 | endif |
---|
| 298 | end do |
---|
| 299 | end do |
---|
| 300 | numparttot=numparttot+ncolumn |
---|
| 301 | if (ipin.eq.0) numpart=numpart+jj |
---|
| 302 | 30 continue |
---|
| 303 | end do |
---|
| 304 | end do |
---|
| 305 | |
---|
| 306 | |
---|
| 307 | ! Check whether numpart is really smaller than maxpart |
---|
| 308 | !***************************************************** |
---|
| 309 | |
---|
| 310 | if (numpart.gt.maxpart) then |
---|
| 311 | write(*,*) 'numpart too large: change source in init_atm_mass.f' |
---|
| 312 | write(*,*) 'numpart: ',numpart,' maxpart: ',maxpart |
---|
| 313 | endif |
---|
| 314 | |
---|
| 315 | |
---|
| 316 | xmassperparticle=colmasstotal/real(numparttot) |
---|
| 317 | |
---|
| 318 | |
---|
| 319 | ! Make sure that all particles are within domain |
---|
| 320 | !*********************************************** |
---|
| 321 | |
---|
| 322 | do j=1,numpart |
---|
| 323 | if ((xtra1(j).lt.0.).or.(xtra1(j).ge.real(nxmin1)).or. & |
---|
| 324 | (ytra1(j).lt.0.).or.(ytra1(j).ge.real(nymin1))) then |
---|
| 325 | itra1(j)=-999999999 |
---|
| 326 | endif |
---|
| 327 | end do |
---|
| 328 | |
---|
| 329 | |
---|
| 330 | |
---|
| 331 | |
---|
| 332 | ! For boundary conditions, we need fewer particle release heights per column, |
---|
| 333 | ! because otherwise it takes too long until enough mass has accumulated to |
---|
| 334 | ! release a particle at the boundary (would take dx/u seconds), leading to |
---|
| 335 | ! relatively large position errors of the order of one grid distance. |
---|
| 336 | ! It's better to release fewer particles per column, but to do so more often. |
---|
| 337 | ! Thus, use on the order of nz starting heights per column. |
---|
| 338 | ! We thus repeat the above to determine fewer starting heights, that are |
---|
| 339 | ! used furtheron in subroutine boundcond_domainfill.f. |
---|
| 340 | !**************************************************************************** |
---|
| 341 | |
---|
| 342 | fractus=real(numcolumn)/real(nz) |
---|
| 343 | write(*,*) 'Total number of particles at model start: ',numpart |
---|
| 344 | write(*,*) 'Maximum number of particles per column: ',numcolumn |
---|
| 345 | write(*,*) 'If ',fractus,' <1, better use more particles' |
---|
| 346 | fractus=sqrt(max(fractus,1.))/2. |
---|
| 347 | |
---|
| 348 | do jy=ny_sn(1),ny_sn(2) ! loop about latitudes |
---|
| 349 | do ix=nx_we(1),nx_we(2) ! loop about longitudes |
---|
| 350 | ncolumn=nint(0.999/fractus*real(npart(1))*colmass(ix,jy) & |
---|
| 351 | /colmasstotal) |
---|
| 352 | if (ncolumn.gt.maxcolumn) stop 'maxcolumn too small' |
---|
| 353 | if (ncolumn.eq.0) goto 80 |
---|
| 354 | |
---|
| 355 | |
---|
| 356 | ! Memorize how many particles per column shall be used for all boundaries |
---|
| 357 | ! This is further used in subroutine boundcond_domainfill.f |
---|
| 358 | ! Use 2 fields for west/east and south/north boundary |
---|
| 359 | !************************************************************************ |
---|
| 360 | |
---|
| 361 | if (ix.eq.nx_we(1)) numcolumn_we(1,jy)=ncolumn |
---|
| 362 | if (ix.eq.nx_we(2)) numcolumn_we(2,jy)=ncolumn |
---|
| 363 | if (jy.eq.ny_sn(1)) numcolumn_sn(1,ix)=ncolumn |
---|
| 364 | if (jy.eq.ny_sn(2)) numcolumn_sn(2,ix)=ncolumn |
---|
| 365 | |
---|
| 366 | ! Calculate pressure at the altitudes of model surfaces, using the air density |
---|
| 367 | ! information, which is stored as a 3-d field |
---|
| 368 | !***************************************************************************** |
---|
| 369 | |
---|
| 370 | do kz=1,nz |
---|
| 371 | pp(kz)=rho(ix,jy,kz,1)*r_air*tt(ix,jy,kz,1) |
---|
| 372 | end do |
---|
| 373 | |
---|
| 374 | ! Determine the reference starting altitudes |
---|
| 375 | !******************************************* |
---|
| 376 | |
---|
| 377 | deltacol=(pp(1)-pp(nz))/real(ncolumn) |
---|
| 378 | pnew=pp(1)+deltacol/2. |
---|
| 379 | do j=1,ncolumn |
---|
| 380 | pnew=pnew-deltacol |
---|
| 381 | do kz=1,nz-1 |
---|
| 382 | if ((pp(kz).ge.pnew).and.(pp(kz+1).lt.pnew)) then |
---|
| 383 | dz1=pp(kz)-pnew |
---|
| 384 | dz2=pnew-pp(kz+1) |
---|
| 385 | dz=1./(dz1+dz2) |
---|
| 386 | zposition=(height(kz)*dz2+height(kz+1)*dz1)*dz |
---|
| 387 | if (zposition.gt.height(nz)-0.5) zposition=height(nz)-0.5 |
---|
| 388 | |
---|
| 389 | ! Memorize vertical positions where particles are introduced |
---|
| 390 | ! This is further used in subroutine boundcond_domainfill.f |
---|
| 391 | !*********************************************************** |
---|
| 392 | |
---|
| 393 | if (ix.eq.nx_we(1)) zcolumn_we(1,jy,j)=zposition |
---|
| 394 | if (ix.eq.nx_we(2)) zcolumn_we(2,jy,j)=zposition |
---|
| 395 | if (jy.eq.ny_sn(1)) zcolumn_sn(1,ix,j)=zposition |
---|
| 396 | if (jy.eq.ny_sn(2)) zcolumn_sn(2,ix,j)=zposition |
---|
| 397 | |
---|
| 398 | ! Initialize mass that has accumulated at boundary to zero |
---|
| 399 | !********************************************************* |
---|
| 400 | |
---|
| 401 | acc_mass_we(1,jy,j)=0. |
---|
| 402 | acc_mass_we(2,jy,j)=0. |
---|
| 403 | acc_mass_sn(1,jy,j)=0. |
---|
| 404 | acc_mass_sn(2,jy,j)=0. |
---|
| 405 | endif |
---|
| 406 | end do |
---|
| 407 | end do |
---|
| 408 | 80 continue |
---|
| 409 | end do |
---|
| 410 | end do |
---|
| 411 | |
---|
| 412 | ! If particles shall be read in to continue an existing run, |
---|
| 413 | ! then the accumulated masses at the domain boundaries must be read in, too. |
---|
| 414 | ! This overrides any previous calculations. |
---|
| 415 | !*************************************************************************** |
---|
| 416 | |
---|
| 417 | if (ipin.eq.1) then |
---|
| 418 | open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', & |
---|
| 419 | form='unformatted') |
---|
| 420 | read(unitboundcond) numcolumn_we,numcolumn_sn, & |
---|
| 421 | zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn |
---|
| 422 | close(unitboundcond) |
---|
| 423 | endif |
---|
| 424 | |
---|
| 425 | |
---|
| 426 | |
---|
| 427 | |
---|
| 428 | end subroutine init_domainfill |
---|