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