[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 | |
---|
| 22 | subroutine timemanager |
---|
| 23 | |
---|
| 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! Handles the computation of trajectories, i.e. determines which * |
---|
| 27 | ! trajectories have to be computed at what time. * |
---|
| 28 | ! Manages dry+wet deposition routines, radioactive decay and the computation * |
---|
| 29 | ! of concentrations. * |
---|
| 30 | ! * |
---|
| 31 | ! Author: A. Stohl * |
---|
| 32 | ! * |
---|
| 33 | ! 20 May 1996 * |
---|
| 34 | ! * |
---|
| 35 | !***************************************************************************** |
---|
| 36 | ! Changes, Bernd C. Krueger, Feb. 2001: * |
---|
| 37 | ! Call of convmix when new windfield is read * |
---|
| 38 | !------------------------------------ * |
---|
| 39 | ! Changes Petra Seibert, Sept 2002 * |
---|
| 40 | ! fix wet scavenging problem * |
---|
| 41 | ! Code may not be correct for decay of deposition! * |
---|
| 42 | ! Changes Petra Seibert, Nov 2002 * |
---|
| 43 | ! call convection BEFORE new fields are read in BWD mode * |
---|
| 44 | ! Changes Caroline Forster, Feb 2005 * |
---|
[8a65cb0] | 45 | ! new interface between flexpart and convection scheme * |
---|
| 46 | ! Emanuel's latest subroutine convect43c.f is used * |
---|
| 47 | ! Changes Stefan Henne, Harald Sodemann, 2013-2014 * |
---|
| 48 | ! added netcdf output code * |
---|
| 49 | ! Changes Espen Sollum 2014 * |
---|
| 50 | ! For compatibility with MPI version, * |
---|
| 51 | ! variables uap,ucp,uzp,us,vs,ws,cbt now in module com_mod * |
---|
[e200b7a] | 52 | !***************************************************************************** |
---|
| 53 | ! * |
---|
| 54 | ! Variables: * |
---|
| 55 | ! DEP .true. if either wet or dry deposition is switched on * |
---|
| 56 | ! decay(maxspec) [1/s] decay constant for radioactive decay * |
---|
| 57 | ! DRYDEP .true. if dry deposition is switched on * |
---|
| 58 | ! ideltas [s] modelling period * |
---|
| 59 | ! itime [s] actual temporal position of calculation * |
---|
| 60 | ! ldeltat [s] time since computation of radioact. decay of depositions* |
---|
| 61 | ! loutaver [s] averaging period for concentration calculations * |
---|
| 62 | ! loutend [s] end of averaging for concentration calculations * |
---|
| 63 | ! loutnext [s] next time at which output fields shall be centered * |
---|
| 64 | ! loutsample [s] sampling interval for averaging of concentrations * |
---|
| 65 | ! loutstart [s] start of averaging for concentration calculations * |
---|
| 66 | ! loutstep [s] time interval for which concentrations shall be * |
---|
| 67 | ! calculated * |
---|
| 68 | ! npoint(maxpart) index, which starting point the trajectory has * |
---|
| 69 | ! starting positions of trajectories * |
---|
| 70 | ! nstop serves as indicator for fate of particles * |
---|
| 71 | ! in the particle loop * |
---|
| 72 | ! nstop1 serves as indicator for wind fields (see getfields) * |
---|
| 73 | ! outnum number of samples for each concentration calculation * |
---|
| 74 | ! outnum number of samples for each concentration calculation * |
---|
| 75 | ! prob probability of absorption at ground due to dry * |
---|
| 76 | ! deposition * |
---|
| 77 | ! WETDEP .true. if wet deposition is switched on * |
---|
| 78 | ! weight weight for each concentration sample (1/2 or 1) * |
---|
| 79 | ! uap(maxpart),ucp(maxpart),uzp(maxpart) = random velocities due to * |
---|
| 80 | ! turbulence * |
---|
| 81 | ! us(maxpart),vs(maxpart),ws(maxpart) = random velocities due to inter- * |
---|
| 82 | ! polation * |
---|
| 83 | ! xtra1(maxpart), ytra1(maxpart), ztra1(maxpart) = * |
---|
| 84 | ! spatial positions of trajectories * |
---|
| 85 | ! * |
---|
| 86 | ! Constants: * |
---|
| 87 | ! maxpart maximum number of trajectories * |
---|
| 88 | ! * |
---|
| 89 | !***************************************************************************** |
---|
| 90 | |
---|
| 91 | use unc_mod |
---|
| 92 | use point_mod |
---|
| 93 | use xmass_mod |
---|
| 94 | use flux_mod |
---|
| 95 | use outg_mod |
---|
| 96 | use oh_mod |
---|
| 97 | use par_mod |
---|
| 98 | use com_mod |
---|
[8a65cb0] | 99 | use netcdf_output_mod, only: concoutput_netcdf,concoutput_nest_netcdf,& |
---|
| 100 | &concoutput_surf_netcdf,concoutput_surf_nest_netcdf |
---|
[e200b7a] | 101 | |
---|
| 102 | implicit none |
---|
| 103 | |
---|
| 104 | integer :: j,ks,kp,l,n,itime,nstop,nstop1 |
---|
| 105 | ! integer :: ksp |
---|
| 106 | integer :: loutnext,loutstart,loutend |
---|
| 107 | integer :: ix,jy,ldeltat,itage,nage |
---|
[8a65cb0] | 108 | integer :: i_nan=0,ii_nan,total_nan_intl=0 !added by mc to check instability in CBL scheme |
---|
| 109 | real :: outnum,weight,prob(maxspec),decfact |
---|
| 110 | ! real :: uap(maxpart),ucp(maxpart),uzp(maxpart) |
---|
| 111 | ! real :: us(maxpart),vs(maxpart),ws(maxpart) |
---|
| 112 | ! integer(kind=2) :: cbt(maxpart) |
---|
[e200b7a] | 113 | real :: drydeposit(maxspec),gridtotalunc,wetgridtotalunc |
---|
| 114 | real :: drygridtotalunc,xold,yold,zold,xmassfract |
---|
| 115 | !double precision xm(maxspec,maxpointspec_act), |
---|
| 116 | ! + xm_depw(maxspec,maxpointspec_act), |
---|
| 117 | ! + xm_depd(maxspec,maxpointspec_act) |
---|
| 118 | |
---|
| 119 | |
---|
| 120 | !open(88,file='TEST.dat') |
---|
| 121 | |
---|
| 122 | ! First output for time 0 |
---|
| 123 | !************************ |
---|
| 124 | |
---|
| 125 | loutnext=loutstep/2 |
---|
| 126 | outnum=0. |
---|
| 127 | loutstart=loutnext-loutaver/2 |
---|
| 128 | loutend=loutnext+loutaver/2 |
---|
| 129 | |
---|
| 130 | ! open(127,file=path(2)(1:length(2))//'depostat.dat' |
---|
| 131 | ! + ,form='unformatted') |
---|
| 132 | !write (*,*) 'writing deposition statistics depostat.dat!' |
---|
| 133 | |
---|
| 134 | !********************************************************************** |
---|
| 135 | ! Loop over the whole modelling period in time steps of mintime seconds |
---|
| 136 | !********************************************************************** |
---|
| 137 | |
---|
[f13406c] | 138 | |
---|
| 139 | if (verbosity.gt.0) then |
---|
| 140 | write (*,*) 'timemanager> starting simulation' |
---|
| 141 | if (verbosity.gt.1) then |
---|
| 142 | CALL SYSTEM_CLOCK(count_clock) |
---|
| 143 | WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate) |
---|
| 144 | endif |
---|
| 145 | endif |
---|
| 146 | |
---|
[e200b7a] | 147 | do itime=0,ideltas,lsynctime |
---|
| 148 | |
---|
| 149 | ! Computation of wet deposition, OH reaction and mass transfer |
---|
| 150 | ! between two species every lsynctime seconds |
---|
| 151 | ! maybe wet depo frequency can be relaxed later but better be on safe side |
---|
| 152 | ! wetdepo must be called BEFORE new fields are read in but should not |
---|
| 153 | ! be called in the very beginning before any fields are loaded, or |
---|
| 154 | ! before particles are in the system |
---|
| 155 | ! Code may not be correct for decay of deposition |
---|
| 156 | ! changed by Petra Seibert 9/02 |
---|
| 157 | !******************************************************************** |
---|
| 158 | |
---|
[f13406c] | 159 | if (WETDEP .and. itime .ne. 0 .and. numpart .gt. 0) then |
---|
| 160 | if (verbosity.gt.0) then |
---|
| 161 | write (*,*) 'timemanager> call wetdepo' |
---|
| 162 | endif |
---|
[e200b7a] | 163 | call wetdepo(itime,lsynctime,loutnext) |
---|
[f13406c] | 164 | endif |
---|
[e200b7a] | 165 | |
---|
| 166 | if (OHREA .and. itime .ne. 0 .and. numpart .gt. 0) & |
---|
| 167 | call ohreaction(itime,lsynctime,loutnext) |
---|
| 168 | |
---|
| 169 | if (ASSSPEC .and. itime .ne. 0 .and. numpart .gt. 0) then |
---|
| 170 | stop 'associated species not yet implemented!' |
---|
| 171 | ! call transferspec(itime,lsynctime,loutnext) |
---|
| 172 | endif |
---|
| 173 | |
---|
| 174 | ! compute convection for backward runs |
---|
| 175 | !************************************* |
---|
| 176 | |
---|
[f13406c] | 177 | if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) then |
---|
| 178 | if (verbosity.gt.0) then |
---|
| 179 | write (*,*) 'timemanager> call convmix -- backward' |
---|
| 180 | endif |
---|
| 181 | call convmix(itime) |
---|
| 182 | if (verbosity.gt.1) then |
---|
| 183 | !CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) |
---|
| 184 | CALL SYSTEM_CLOCK(count_clock) |
---|
| 185 | WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate) |
---|
| 186 | endif |
---|
| 187 | endif |
---|
[e200b7a] | 188 | |
---|
| 189 | ! Get necessary wind fields if not available |
---|
| 190 | !******************************************* |
---|
[f13406c] | 191 | if (verbosity.gt.0) then |
---|
| 192 | write (*,*) 'timemanager> call getfields' |
---|
| 193 | endif |
---|
[e200b7a] | 194 | call getfields(itime,nstop1) |
---|
[f13406c] | 195 | if (verbosity.gt.1) then |
---|
| 196 | CALL SYSTEM_CLOCK(count_clock) |
---|
| 197 | WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate) |
---|
| 198 | endif |
---|
[e200b7a] | 199 | if (nstop1.gt.1) stop 'NO METEO FIELDS AVAILABLE' |
---|
[f13406c] | 200 | |
---|
[8a65cb0] | 201 | ! Get hourly OH fields if not available |
---|
| 202 | !**************************************************** |
---|
| 203 | if (OHREA) then |
---|
| 204 | if (verbosity.gt.0) then |
---|
| 205 | write (*,*) 'timemanager> call gethourlyOH' |
---|
| 206 | endif |
---|
| 207 | call gethourlyOH(itime) |
---|
| 208 | if (verbosity.gt.1) then |
---|
| 209 | CALL SYSTEM_CLOCK(count_clock) |
---|
| 210 | WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate) |
---|
| 211 | endif |
---|
| 212 | endif |
---|
| 213 | |
---|
[e200b7a] | 214 | ! Release particles |
---|
| 215 | !****************** |
---|
| 216 | |
---|
[8a65cb0] | 217 | if (verbosity.gt.0) then |
---|
| 218 | write (*,*) 'timemanager> Release particles' |
---|
| 219 | endif |
---|
[f13406c] | 220 | |
---|
[e200b7a] | 221 | if (mdomainfill.ge.1) then |
---|
| 222 | if (itime.eq.0) then |
---|
[f13406c] | 223 | if (verbosity.gt.0) then |
---|
| 224 | write (*,*) 'timemanager> call init_domainfill' |
---|
| 225 | endif |
---|
[e200b7a] | 226 | call init_domainfill |
---|
| 227 | else |
---|
[f13406c] | 228 | if (verbosity.gt.0) then |
---|
| 229 | write (*,*) 'timemanager> call boundcond_domainfill' |
---|
| 230 | endif |
---|
[e200b7a] | 231 | call boundcond_domainfill(itime,loutend) |
---|
| 232 | endif |
---|
| 233 | else |
---|
[f13406c] | 234 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 235 | print*,'call releaseparticles' |
---|
[f13406c] | 236 | endif |
---|
[e200b7a] | 237 | call releaseparticles(itime) |
---|
[f13406c] | 238 | if (verbosity.gt.1) then |
---|
| 239 | CALL SYSTEM_CLOCK(count_clock) |
---|
| 240 | WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate) |
---|
| 241 | endif |
---|
[e200b7a] | 242 | endif |
---|
| 243 | |
---|
| 244 | |
---|
| 245 | ! Compute convective mixing for forward runs |
---|
| 246 | ! for backward runs it is done before next windfield is read in |
---|
| 247 | !************************************************************** |
---|
| 248 | |
---|
[f13406c] | 249 | if ((ldirect.eq.1).and.(lconvection.eq.1)) then |
---|
| 250 | if (verbosity.gt.0) then |
---|
| 251 | write (*,*) 'timemanager> call convmix -- forward' |
---|
| 252 | endif |
---|
| 253 | call convmix(itime) |
---|
| 254 | endif |
---|
[e200b7a] | 255 | |
---|
| 256 | ! If middle of averaging period of output fields is reached, accumulated |
---|
| 257 | ! deposited mass radioactively decays |
---|
| 258 | !*********************************************************************** |
---|
| 259 | |
---|
| 260 | if (DEP.and.(itime.eq.loutnext).and.(ldirect.gt.0)) then |
---|
| 261 | do ks=1,nspec |
---|
| 262 | do kp=1,maxpointspec_act |
---|
| 263 | if (decay(ks).gt.0.) then |
---|
| 264 | do nage=1,nageclass |
---|
| 265 | do l=1,nclassunc |
---|
| 266 | ! Mother output grid |
---|
| 267 | do jy=0,numygrid-1 |
---|
| 268 | do ix=0,numxgrid-1 |
---|
| 269 | wetgridunc(ix,jy,ks,kp,l,nage)= & |
---|
| 270 | wetgridunc(ix,jy,ks,kp,l,nage)* & |
---|
| 271 | exp(-1.*outstep*decay(ks)) |
---|
| 272 | drygridunc(ix,jy,ks,kp,l,nage)= & |
---|
| 273 | drygridunc(ix,jy,ks,kp,l,nage)* & |
---|
| 274 | exp(-1.*outstep*decay(ks)) |
---|
| 275 | end do |
---|
| 276 | end do |
---|
| 277 | ! Nested output grid |
---|
| 278 | if (nested_output.eq.1) then |
---|
| 279 | do jy=0,numygridn-1 |
---|
| 280 | do ix=0,numxgridn-1 |
---|
| 281 | wetgriduncn(ix,jy,ks,kp,l,nage)= & |
---|
| 282 | wetgriduncn(ix,jy,ks,kp,l,nage)* & |
---|
| 283 | exp(-1.*outstep*decay(ks)) |
---|
| 284 | drygriduncn(ix,jy,ks,kp,l,nage)= & |
---|
| 285 | drygriduncn(ix,jy,ks,kp,l,nage)* & |
---|
| 286 | exp(-1.*outstep*decay(ks)) |
---|
| 287 | end do |
---|
| 288 | end do |
---|
| 289 | endif |
---|
| 290 | end do |
---|
| 291 | end do |
---|
| 292 | endif |
---|
| 293 | end do |
---|
| 294 | end do |
---|
| 295 | endif |
---|
| 296 | |
---|
| 297 | !!! CHANGE: These lines may be switched on to check the conservation |
---|
| 298 | !!! of mass within FLEXPART |
---|
| 299 | ! if (itime.eq.loutnext) then |
---|
| 300 | ! do 247 ksp=1, nspec |
---|
| 301 | ! do 247 kp=1, maxpointspec_act |
---|
| 302 | !47 xm(ksp,kp)=0. |
---|
| 303 | |
---|
| 304 | ! do 249 ksp=1, nspec |
---|
| 305 | ! do 249 j=1,numpart |
---|
| 306 | ! if (ioutputforeachrelease.eq.1) then |
---|
| 307 | ! kp=npoint(j) |
---|
| 308 | ! else |
---|
| 309 | ! kp=1 |
---|
| 310 | ! endif |
---|
| 311 | ! if (itra1(j).eq.itime) then |
---|
| 312 | ! xm(ksp,kp)=xm(ksp,kp)+xmass1(j,ksp) |
---|
| 313 | ! write(*,*) 'xmass: ',xmass1(j,ksp),j,ksp,nspec |
---|
| 314 | ! endif |
---|
| 315 | !49 continue |
---|
| 316 | ! do 248 ksp=1,nspec |
---|
| 317 | ! do 248 kp=1,maxpointspec_act |
---|
| 318 | ! xm_depw(ksp,kp)=0. |
---|
| 319 | ! xm_depd(ksp,kp)=0. |
---|
| 320 | ! do 248 nage=1,nageclass |
---|
| 321 | ! do 248 ix=0,numxgrid-1 |
---|
| 322 | ! do 248 jy=0,numygrid-1 |
---|
| 323 | ! do 248 l=1,nclassunc |
---|
| 324 | ! xm_depw(ksp,kp)=xm_depw(ksp,kp) |
---|
| 325 | ! + +wetgridunc(ix,jy,ksp,kp,l,nage) |
---|
| 326 | !48 xm_depd(ksp,kp)=xm_depd(ksp,kp) |
---|
| 327 | ! + +drygridunc(ix,jy,ksp,kp,l,nage) |
---|
| 328 | ! do 246 ksp=1,nspec |
---|
| 329 | !46 write(88,'(2i10,3e12.3)') |
---|
| 330 | ! + itime,ksp,(xm(ksp,kp),kp=1,maxpointspec_act), |
---|
| 331 | ! + (xm_depw(ksp,kp),kp=1,maxpointspec_act), |
---|
| 332 | ! + (xm_depd(ksp,kp),kp=1,maxpointspec_act) |
---|
| 333 | ! endif |
---|
| 334 | !!! CHANGE |
---|
| 335 | |
---|
| 336 | |
---|
| 337 | |
---|
| 338 | ! Check whether concentrations are to be calculated |
---|
| 339 | !************************************************** |
---|
| 340 | |
---|
| 341 | if ((ldirect*itime.ge.ldirect*loutstart).and. & |
---|
| 342 | (ldirect*itime.le.ldirect*loutend)) then ! add to grid |
---|
| 343 | if (mod(itime-loutstart,loutsample).eq.0) then |
---|
| 344 | |
---|
| 345 | ! If we are exactly at the start or end of the concentration averaging interval, |
---|
| 346 | ! give only half the weight to this sample |
---|
| 347 | !***************************************************************************** |
---|
| 348 | |
---|
| 349 | if ((itime.eq.loutstart).or.(itime.eq.loutend)) then |
---|
| 350 | weight=0.5 |
---|
| 351 | else |
---|
| 352 | weight=1.0 |
---|
| 353 | endif |
---|
| 354 | outnum=outnum+weight |
---|
| 355 | call conccalc(itime,weight) |
---|
| 356 | endif |
---|
| 357 | |
---|
| 358 | |
---|
| 359 | if ((mquasilag.eq.1).and.(itime.eq.(loutstart+loutend)/2)) & |
---|
| 360 | call partoutput_short(itime) ! dump particle positions in extremely compressed format |
---|
| 361 | |
---|
| 362 | |
---|
| 363 | ! Output and reinitialization of grid |
---|
| 364 | ! If necessary, first sample of new grid is also taken |
---|
| 365 | !***************************************************** |
---|
| 366 | |
---|
| 367 | if ((itime.eq.loutend).and.(outnum.gt.0.)) then |
---|
| 368 | if ((iout.le.3.).or.(iout.eq.5)) then |
---|
[f13406c] | 369 | if (surf_only.ne.1) then |
---|
[8a65cb0] | 370 | if (lnetcdfout.eq.1) then |
---|
| 371 | call concoutput_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc) |
---|
| 372 | else |
---|
| 373 | call concoutput(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc) |
---|
| 374 | endif |
---|
[f13406c] | 375 | else |
---|
[8a65cb0] | 376 | if (verbosity.eq.1) then |
---|
| 377 | print*,'call concoutput_surf ' |
---|
| 378 | call system_clock(count_clock) |
---|
| 379 | write(*,*) 'system clock',count_clock - count_clock0 |
---|
| 380 | endif |
---|
| 381 | if (lnetcdfout.eq.1) then |
---|
| 382 | call concoutput_surf_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc) |
---|
| 383 | else |
---|
| 384 | call concoutput_surf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc) |
---|
| 385 | if (verbosity.eq.1) then |
---|
| 386 | print*,'called concoutput_surf ' |
---|
| 387 | call system_clock(count_clock) |
---|
| 388 | write(*,*) 'system clock',count_clock - count_clock0 |
---|
| 389 | endif |
---|
| 390 | endif |
---|
[f13406c] | 391 | endif |
---|
| 392 | |
---|
[8a65cb0] | 393 | if (nested_output .eq. 1) then |
---|
| 394 | if (lnetcdfout.eq.0) then |
---|
| 395 | if (surf_only.ne.1) then |
---|
| 396 | call concoutput_nest(itime,outnum) |
---|
| 397 | else |
---|
| 398 | call concoutput_surf_nest(itime,outnum) |
---|
| 399 | endif |
---|
| 400 | else |
---|
| 401 | if (surf_only.ne.1) then |
---|
| 402 | call concoutput_nest_netcdf(itime,outnum) |
---|
| 403 | else |
---|
| 404 | call concoutput_surf_nest_netcdf(itime,outnum) |
---|
| 405 | endif |
---|
| 406 | endif |
---|
| 407 | endif |
---|
[e200b7a] | 408 | outnum=0. |
---|
| 409 | endif |
---|
| 410 | if ((iout.eq.4).or.(iout.eq.5)) call plumetraj(itime) |
---|
| 411 | if (iflux.eq.1) call fluxoutput(itime) |
---|
[8a65cb0] | 412 | write(*,45) itime,numpart,gridtotalunc,wetgridtotalunc,drygridtotalunc |
---|
| 413 | !write(*,46) float(itime)/3600,itime,numpart |
---|
[b4d29ce] | 414 | 45 format(i9,' SECONDS SIMULATED: ',i8, ' PARTICLES: Uncertainty: ',3f7.3) |
---|
| 415 | 46 format(' Simulated ',f7.1,' hours (',i9,' s), ',i8, ' particles') |
---|
[e200b7a] | 416 | if (ipout.ge.1) call partoutput(itime) ! dump particle positions |
---|
| 417 | loutnext=loutnext+loutstep |
---|
| 418 | loutstart=loutnext-loutaver/2 |
---|
| 419 | loutend=loutnext+loutaver/2 |
---|
| 420 | if (itime.eq.loutstart) then |
---|
| 421 | weight=0.5 |
---|
| 422 | outnum=outnum+weight |
---|
| 423 | call conccalc(itime,weight) |
---|
| 424 | endif |
---|
| 425 | |
---|
| 426 | |
---|
| 427 | ! Check, whether particles are to be split: |
---|
| 428 | ! If so, create new particles and attribute all information from the old |
---|
| 429 | ! particles also to the new ones; old and new particles both get half the |
---|
| 430 | ! mass of the old ones |
---|
| 431 | !************************************************************************ |
---|
| 432 | |
---|
| 433 | if (ldirect*itime.ge.ldirect*itsplit) then |
---|
| 434 | n=numpart |
---|
| 435 | do j=1,numpart |
---|
| 436 | if (ldirect*itime.ge.ldirect*itrasplit(j)) then |
---|
| 437 | if (n.lt.maxpart) then |
---|
| 438 | n=n+1 |
---|
| 439 | itrasplit(j)=2*(itrasplit(j)-itramem(j))+itramem(j) |
---|
| 440 | itrasplit(n)=itrasplit(j) |
---|
| 441 | itramem(n)=itramem(j) |
---|
| 442 | itra1(n)=itra1(j) |
---|
| 443 | idt(n)=idt(j) |
---|
| 444 | npoint(n)=npoint(j) |
---|
| 445 | nclass(n)=nclass(j) |
---|
| 446 | xtra1(n)=xtra1(j) |
---|
| 447 | ytra1(n)=ytra1(j) |
---|
| 448 | ztra1(n)=ztra1(j) |
---|
| 449 | uap(n)=uap(j) |
---|
| 450 | ucp(n)=ucp(j) |
---|
| 451 | uzp(n)=uzp(j) |
---|
| 452 | us(n)=us(j) |
---|
| 453 | vs(n)=vs(j) |
---|
| 454 | ws(n)=ws(j) |
---|
| 455 | cbt(n)=cbt(j) |
---|
| 456 | do ks=1,nspec |
---|
| 457 | xmass1(j,ks)=xmass1(j,ks)/2. |
---|
| 458 | xmass1(n,ks)=xmass1(j,ks) |
---|
| 459 | end do |
---|
| 460 | endif |
---|
| 461 | endif |
---|
| 462 | end do |
---|
| 463 | numpart=n |
---|
| 464 | endif |
---|
| 465 | endif |
---|
| 466 | endif |
---|
| 467 | |
---|
| 468 | |
---|
| 469 | if (itime.eq.ideltas) exit ! almost finished |
---|
| 470 | |
---|
| 471 | ! Compute interval since radioactive decay of deposited mass was computed |
---|
| 472 | !************************************************************************ |
---|
| 473 | |
---|
| 474 | if (itime.lt.loutnext) then |
---|
| 475 | ldeltat=itime-(loutnext-loutstep) |
---|
| 476 | else ! first half of next interval |
---|
| 477 | ldeltat=itime-loutnext |
---|
| 478 | endif |
---|
| 479 | |
---|
| 480 | |
---|
| 481 | ! Loop over all particles |
---|
| 482 | !************************ |
---|
[8a65cb0] | 483 | ! Various variables for testing reason of CBL scheme, by mc |
---|
| 484 | well_mixed_vector=0. !erase vector to test well mixed condition: modified by mc |
---|
| 485 | well_mixed_norm=0. !erase normalization to test well mixed condition: modified by mc |
---|
| 486 | avg_ol=0. |
---|
| 487 | avg_wst=0. |
---|
| 488 | avg_h=0. |
---|
| 489 | avg_air_dens=0. !erase vector to obtain air density at particle positions: modified by mc |
---|
| 490 | !----------------------------------------------------------------------------- |
---|
[e200b7a] | 491 | do j=1,numpart |
---|
| 492 | |
---|
| 493 | |
---|
| 494 | ! If integration step is due, do it |
---|
| 495 | !********************************** |
---|
| 496 | |
---|
| 497 | if (itra1(j).eq.itime) then |
---|
| 498 | |
---|
| 499 | if (ioutputforeachrelease.eq.1) then |
---|
| 500 | kp=npoint(j) |
---|
| 501 | else |
---|
| 502 | kp=1 |
---|
| 503 | endif |
---|
| 504 | ! Determine age class of the particle |
---|
| 505 | itage=abs(itra1(j)-itramem(j)) |
---|
| 506 | do nage=1,nageclass |
---|
| 507 | if (itage.lt.lage(nage)) exit |
---|
| 508 | end do |
---|
| 509 | |
---|
| 510 | ! Initialize newly released particle |
---|
| 511 | !*********************************** |
---|
| 512 | |
---|
| 513 | if ((itramem(j).eq.itime).or.(itime.eq.0)) & |
---|
| 514 | call initialize(itime,idt(j),uap(j),ucp(j),uzp(j), & |
---|
| 515 | us(j),vs(j),ws(j),xtra1(j),ytra1(j),ztra1(j),cbt(j)) |
---|
| 516 | |
---|
| 517 | ! Memorize particle positions |
---|
| 518 | !**************************** |
---|
| 519 | |
---|
| 520 | xold=xtra1(j) |
---|
| 521 | yold=ytra1(j) |
---|
| 522 | zold=ztra1(j) |
---|
| 523 | |
---|
| 524 | ! Integrate Lagevin equation for lsynctime seconds |
---|
| 525 | !************************************************* |
---|
| 526 | |
---|
| 527 | call advance(itime,npoint(j),idt(j),uap(j),ucp(j),uzp(j), & |
---|
| 528 | us(j),vs(j),ws(j),nstop,xtra1(j),ytra1(j),ztra1(j),prob, & |
---|
| 529 | cbt(j)) |
---|
| 530 | |
---|
| 531 | ! Calculate the gross fluxes across layer interfaces |
---|
| 532 | !*************************************************** |
---|
| 533 | |
---|
| 534 | if (iflux.eq.1) call calcfluxes(nage,j,xold,yold,zold) |
---|
| 535 | |
---|
| 536 | |
---|
| 537 | ! Determine, when next time step is due |
---|
| 538 | ! If trajectory is terminated, mark it |
---|
| 539 | !************************************** |
---|
| 540 | |
---|
| 541 | if (nstop.gt.1) then |
---|
| 542 | if (linit_cond.ge.1) call initial_cond_calc(itime,j) |
---|
| 543 | itra1(j)=-999999999 |
---|
| 544 | else |
---|
| 545 | itra1(j)=itime+lsynctime |
---|
| 546 | |
---|
| 547 | |
---|
| 548 | ! Dry deposition and radioactive decay for each species |
---|
| 549 | ! Also check maximum (of all species) of initial mass remaining on the particle; |
---|
| 550 | ! if it is below a threshold value, terminate particle |
---|
| 551 | !***************************************************************************** |
---|
| 552 | |
---|
| 553 | xmassfract=0. |
---|
| 554 | do ks=1,nspec |
---|
| 555 | if (decay(ks).gt.0.) then ! radioactive decay |
---|
| 556 | decfact=exp(-real(abs(lsynctime))*decay(ks)) |
---|
| 557 | else |
---|
| 558 | decfact=1. |
---|
| 559 | endif |
---|
| 560 | |
---|
| 561 | if (DRYDEPSPEC(ks)) then ! dry deposition |
---|
| 562 | drydeposit(ks)=xmass1(j,ks)*prob(ks)*decfact |
---|
| 563 | xmass1(j,ks)=xmass1(j,ks)*(1.-prob(ks))*decfact |
---|
| 564 | if (decay(ks).gt.0.) then ! correct for decay (see wetdepo) |
---|
| 565 | drydeposit(ks)=drydeposit(ks)* & |
---|
| 566 | exp(real(abs(ldeltat))*decay(ks)) |
---|
| 567 | endif |
---|
| 568 | else ! no dry deposition |
---|
| 569 | xmass1(j,ks)=xmass1(j,ks)*decfact |
---|
| 570 | endif |
---|
| 571 | |
---|
| 572 | |
---|
| 573 | if (mdomainfill.eq.0) then |
---|
| 574 | if (xmass(npoint(j),ks).gt.0.) & |
---|
| 575 | xmassfract=max(xmassfract,real(npart(npoint(j)))* & |
---|
| 576 | xmass1(j,ks)/xmass(npoint(j),ks)) |
---|
| 577 | else |
---|
| 578 | xmassfract=1. |
---|
| 579 | endif |
---|
| 580 | end do |
---|
| 581 | |
---|
| 582 | if (xmassfract.lt.0.0001) then ! terminate all particles carrying less mass |
---|
| 583 | itra1(j)=-999999999 |
---|
[8a65cb0] | 584 | if (verbosity.gt.0) then |
---|
| 585 | print*,'terminated particle ',j,' for small mass' |
---|
| 586 | endif |
---|
[e200b7a] | 587 | endif |
---|
| 588 | |
---|
| 589 | ! Sabine Eckhardt, June 2008 |
---|
| 590 | ! don't create depofield for backward runs |
---|
| 591 | if (DRYDEP.AND.(ldirect.eq.1)) then |
---|
| 592 | call drydepokernel(nclass(j),drydeposit,real(xtra1(j)), & |
---|
| 593 | real(ytra1(j)),nage,kp) |
---|
| 594 | if (nested_output.eq.1) call drydepokernel_nest( & |
---|
| 595 | nclass(j),drydeposit,real(xtra1(j)),real(ytra1(j)), & |
---|
| 596 | nage,kp) |
---|
| 597 | endif |
---|
| 598 | |
---|
| 599 | ! Terminate trajectories that are older than maximum allowed age |
---|
| 600 | !*************************************************************** |
---|
| 601 | |
---|
| 602 | if (abs(itra1(j)-itramem(j)).ge.lage(nageclass)) then |
---|
[8a65cb0] | 603 | if (linit_cond.ge.1) call initial_cond_calc(itime+lsynctime,j) |
---|
[e200b7a] | 604 | itra1(j)=-999999999 |
---|
[8a65cb0] | 605 | if (verbosity.gt.0) then |
---|
| 606 | print*,'terminated particle ',j,' for age' |
---|
| 607 | endif |
---|
[e200b7a] | 608 | endif |
---|
| 609 | endif |
---|
| 610 | |
---|
| 611 | endif |
---|
| 612 | |
---|
[8a65cb0] | 613 | end do !loop over particles |
---|
| 614 | |
---|
| 615 | ! Counter of "unstable" particle velocity during a time scale of |
---|
| 616 | ! maximumtl=20 minutes (defined in com_mod) |
---|
| 617 | !*************************************************************** |
---|
| 618 | |
---|
| 619 | total_nan_intl=0 |
---|
| 620 | i_nan=i_nan+1 ! added by mc to count nan during a time of maxtl (i.e. maximum tl fixed here to 20 minutes, see com_mod) |
---|
| 621 | sum_nan_count(i_nan)=nan_count |
---|
| 622 | if (i_nan > maxtl/lsynctime) i_nan=1 !lsynctime must be <= maxtl |
---|
| 623 | do ii_nan=1, (maxtl/lsynctime) |
---|
| 624 | total_nan_intl=total_nan_intl+sum_nan_count(ii_nan) |
---|
| 625 | end do |
---|
| 626 | ! Output to keep track of the numerical instabilities in CBL simulation and if |
---|
| 627 | ! they are compromising the final result (or not) |
---|
| 628 | if (cblflag.eq.1) print *,j,itime,'nan_synctime',nan_count,'nan_tl',total_nan_intl |
---|
| 629 | |
---|
| 630 | !!------------------------------------------------------------------------------ |
---|
| 631 | ! these lines below to test the well-mixed condition, modified by mc, not to be |
---|
| 632 | ! included in final release: |
---|
| 633 | ! if (itime.eq.0) then |
---|
| 634 | ! open(551,file='/home/mc/test_cbl/out/WELLMIXEDTEST_CBL_lonlat_9_33_100_3hours_3htp_cd.DAT') |
---|
| 635 | ! open(552,file='/home/mc/test_cbl/out/avg_ol_h_wst_lonlat_9_33_100_3hours_3htp_cd.DAT') |
---|
| 636 | ! end if |
---|
| 637 | ! write(552,'(5F16.7)')itime*1./3600.,avg_wst/well_mixed_norm,avg_ol/well_mixed_norm,avg_h/well_mixed_norm |
---|
| 638 | ! do j=1,25 |
---|
| 639 | ! !write(551,*))itime*1.,h_well/50.*j,well_mixed_vector(j)/well_mixed_norm*50. |
---|
| 640 | ! avg_air_dens(j)=avg_air_dens(j)/well_mixed_vector(j) |
---|
| 641 | ! write(551,'(5F16.7)')itime*1.,h_well/25.*j,well_mixed_vector(j)/well_mixed_norm*25.,avg_air_dens(j),0.04*j |
---|
| 642 | ! end do |
---|
| 643 | !!------------------------------------------------------------------------------ |
---|
[e200b7a] | 644 | |
---|
| 645 | end do |
---|
| 646 | |
---|
| 647 | |
---|
| 648 | ! Complete the calculation of initial conditions for particles not yet terminated |
---|
| 649 | !***************************************************************************** |
---|
| 650 | |
---|
| 651 | do j=1,numpart |
---|
| 652 | if (linit_cond.ge.1) call initial_cond_calc(itime,j) |
---|
| 653 | end do |
---|
| 654 | |
---|
| 655 | if (ipout.eq.2) call partoutput(itime) ! dump particle positions |
---|
| 656 | |
---|
| 657 | if (linit_cond.ge.1) call initial_cond_output(itime) ! dump initial cond. field |
---|
| 658 | |
---|
[78e62dc] | 659 | !close(104) |
---|
[e200b7a] | 660 | |
---|
| 661 | ! De-allocate memory and end |
---|
| 662 | !*************************** |
---|
| 663 | |
---|
| 664 | if (iflux.eq.1) then |
---|
| 665 | deallocate(flux) |
---|
| 666 | endif |
---|
[8a65cb0] | 667 | if (OHREA) then |
---|
| 668 | deallocate(OH_field,OH_hourly,lonOH,latOH,altOH) |
---|
[e200b7a] | 669 | endif |
---|
| 670 | if (ldirect.gt.0) then |
---|
| 671 | deallocate(drygridunc,wetgridunc) |
---|
| 672 | endif |
---|
| 673 | deallocate(gridunc) |
---|
| 674 | deallocate(xpoint1,xpoint2,ypoint1,ypoint2,zpoint1,zpoint2,xmass) |
---|
| 675 | deallocate(ireleasestart,ireleaseend,npart,kindz) |
---|
| 676 | deallocate(xmasssave) |
---|
| 677 | if (nested_output.eq.1) then |
---|
| 678 | deallocate(orooutn, arean, volumen) |
---|
| 679 | if (ldirect.gt.0) then |
---|
| 680 | deallocate(griduncn,drygriduncn,wetgriduncn) |
---|
| 681 | endif |
---|
| 682 | endif |
---|
| 683 | deallocate(outheight,outheighthalf) |
---|
| 684 | deallocate(oroout, area, volume) |
---|
| 685 | |
---|
| 686 | end subroutine timemanager |
---|
| 687 | |
---|