56 real(kind=
dp),
allocatable :: temp_time(:)
58 real(kind=
dp),
allocatable :: temp(:)
61 real(kind=
dp),
allocatable :: pressure_time(:)
63 real(kind=
dp),
allocatable :: pressure(:)
66 real(kind=
dp),
allocatable :: height_time(:)
68 real(kind=
dp),
allocatable :: height(:)
71 real(kind=
dp),
allocatable :: gas_emission_time(:)
73 real(kind=
dp),
allocatable :: gas_emission_rate_scale(:)
78 real(kind=
dp),
allocatable :: gas_dilution_time(:)
80 real(kind=
dp),
allocatable :: gas_dilution_rate(:)
85 real(kind=
dp),
allocatable :: aero_emission_time(:)
87 real(kind=
dp),
allocatable :: aero_emission_rate_scale(:)
92 real(kind=
dp),
allocatable :: aero_dilution_time(:)
94 real(kind=
dp),
allocatable :: aero_dilution_rate(:)
99 integer :: loss_function_type
117 real(kind=
dp),
intent(in) :: time
119 env_state%temp =
interp_1d(scenario%temp_time, scenario%temp, time)
120 env_state%pressure =
interp_1d(scenario%pressure_time, &
121 scenario%pressure, time)
122 env_state%height =
interp_1d(scenario%height_time, scenario%height, time)
123 env_state%elapsed_time = time
125 env_state%solar_zenith_angle = 0d0
140 real(kind=
dp),
intent(in) :: time
143 real(kind=
dp) :: pmv_old, pmv_new
145 real(kind=
dp) :: pressure_old
147 real(kind=
dp) :: temp_old
153 pressure_old = env_state%pressure
154 temp_old = env_state%temp
156 env_state%temp =
interp_1d(scenario%temp_time, scenario%temp, time)
157 env_state%pressure =
interp_1d(scenario%pressure_time, &
158 scenario%pressure, time)
160 pmv_new = pmv_old * env_state%pressure / pressure_old
163 env_state%height =
interp_1d(scenario%height_time, scenario%height, time)
164 env_state%elapsed_time = time
208 old_env_state, gas_data, gas_state)
213 real(kind=
dp),
intent(in) :: delta_t
223 real(kind=
dp) :: emission_rate_scale, dilution_rate, p
227 if (
size(scenario%gas_emission) > 0)
then
229 scenario%gas_emission_time, scenario%gas_emission_rate_scale, &
230 env_state%elapsed_time, emissions, emission_rate_scale)
232 p = emission_rate_scale * delta_t / env_state%height
238 scenario%gas_dilution_time, scenario%gas_dilution_rate, &
239 env_state%elapsed_time, background, dilution_rate)
240 p = exp(- dilution_rate * delta_t)
241 if (env_state%height > old_env_state%height)
then
242 p = p * old_env_state%height / env_state%height
260 old_env_state, aero_data, aero_state, n_emit, n_dil_in, n_dil_out, &
261 allow_doubling, allow_halving)
266 real(kind=
dp),
intent(in) :: delta_t
276 integer,
intent(out) :: n_emit
278 integer,
intent(out) :: n_dil_in
280 integer,
intent(out) :: n_dil_out
282 logical,
intent(in) :: allow_doubling
284 logical,
intent(in) :: allow_halving
286 real(kind=
dp),
parameter :: sample_timescale = 3600.0d0
287 real(kind=
dp) :: characteristic_factor, sample_timescale_effective
288 real(kind=
dp) :: emission_rate_scale, dilution_rate, p
292 sample_timescale_effective = max(1.0, min(sample_timescale, &
293 env_state%elapsed_time))
294 characteristic_factor = sample_timescale_effective / delta_t
297 if (
size(scenario%aero_emission) > 0)
then
299 scenario%aero_emission_time, scenario%aero_emission_rate_scale, &
300 env_state%elapsed_time, emissions, emission_rate_scale)
301 p = emission_rate_scale * delta_t / env_state%height
303 emissions, p, characteristic_factor, env_state%elapsed_time, &
304 allow_doubling, allow_halving, n_emit)
309 scenario%aero_dilution_time, scenario%aero_dilution_rate, &
310 env_state%elapsed_time, background, dilution_rate)
311 p = exp(- dilution_rate * delta_t)
312 if (env_state%height > old_env_state%height)
then
313 p = p * old_env_state%height / env_state%height
319 aero_data, 1d0 - p, aero_info_dilution)
323 background, 1d0 - p, characteristic_factor, env_state%elapsed_time, &
324 allow_doubling, allow_halving, n_dil_in)
332 call aero_weight_array_scale(aero_state%awa, &
333 old_env_state%inverse_density * (1.0d0 / env_state%inverse_density))
336 call aero_weight_array_scale(aero_state%awa, &
337 old_env_state%temp * env_state%pressure &
338 / (env_state%temp * old_env_state%pressure))
350 old_env_state, bin_grid, aero_data, aero_binned)
355 real(kind=
dp),
intent(in) :: delta_t
365 type(aero_binned_t),
intent(inout) :: aero_binned
367 real(kind=
dp) :: emission_rate_scale, dilution_rate, p
369 type(aero_binned_t) :: emissions_binned, background_binned
373 scenario%aero_emission_time, scenario%aero_emission_rate_scale, &
374 env_state%elapsed_time, emissions, emission_rate_scale)
375 call aero_binned_add_aero_dist(emissions_binned, bin_grid, aero_data, &
377 p = emission_rate_scale * delta_t / env_state%height
378 call aero_binned_add_scaled(aero_binned, emissions_binned, p)
382 scenario%aero_dilution_time, scenario%aero_dilution_rate, &
383 env_state%elapsed_time, background, dilution_rate)
384 call aero_binned_add_aero_dist(background_binned, bin_grid, aero_data, &
386 p = exp(- dilution_rate * delta_t)
387 if (env_state%height > old_env_state%height)
then
388 p = p * old_env_state%height / env_state%height
390 call aero_binned_scale(aero_binned, p)
391 call aero_binned_add_scaled(aero_binned, background_binned, 1d0 - p)
399 aero_data, env_state)
404 real(kind=
dp),
intent(in) :: vol
406 real(kind=
dp),
intent(in) :: density
425 aero_data, env_state)
429 aero_data, env_state) &
431 aero_data, env_state)
433 call die_msg(201594391,
"Unknown loss function id: " &
448 real(kind=
dp),
intent(in) :: vol
450 real(kind=
dp),
intent(in) :: density
459 real(kind=
dp) :: density_air
460 real(kind=
dp) :: visc_d, visc_k
461 real(kind=
dp) :: gas_speed, gas_mean_free_path
462 real(kind=
dp) :: knud, cunning
463 real(kind=
dp) :: grav
464 real(kind=
dp) :: r_s, r_a
465 real(kind=
dp) :: alpha, beta, gamma, a, eps_0
466 real(kind=
dp) :: diff_p
467 real(kind=
dp) :: von_karman
468 real(kind=
dp) :: st, sc, u_star
469 real(kind=
dp) :: e_b, e_im, e_in, r1
470 real(kind=
dp) :: u_mean, z_ref, z_rough
485 density_air = (
const%air_molec_weight * env_state%pressure) &
486 / (
const%univ_gas_const * env_state%temp)
488 visc_d = 1.8325d-5 * (416.16 / (env_state%temp + 120.0d0)) &
489 * (env_state%temp / 296.16)**1.5d0
491 visc_k = visc_d / density_air
494 sqrt((8.0d0 *
const%boltzmann * env_state%temp *
const%avagadro) / &
497 gas_mean_free_path = (2.0d0 * visc_d) / (density_air * gas_speed)
499 knud = (2.0d0 * gas_mean_free_path) / d_p
501 cunning = 1.0d0 + knud * (1.257d0 + 0.4d0 * exp(-1.1d0 / knud))
505 v_s = (density * d_p**2.0d0 * grav * cunning) / (18.0d0 * visc_d)
509 u_star = .4d0 * u_mean / log(z_ref / z_rough)
510 r_a = (1.0d0 / (.4d0 * u_star)) * log(z_ref / z_rough)
512 diff_p = (
const%boltzmann * env_state%temp * cunning) / &
513 (3.d0 *
const%pi * visc_d * d_p)
519 e_in = .5d0 * (d_p / a)**2.0d0
522 st = (v_s * u_star) / (grav * a)
523 e_im = (st / (alpha + st))**beta
530 r_s = 1.0d0 / (eps_0 * u_star * (e_b + e_in + e_im) * r1)
533 v_d = v_s + (1.0d0 / (r_a + r_s + r_a * r_s * v_s))
549 real(kind=
dp),
intent(in) :: vol
556 integer,
parameter :: n_sample = 3
558 real(kind=
dp) :: d, d_min, d_max, loss
561 d_min = minval(aero_data%density)
562 d_max = maxval(aero_data%density)
595 integer,
parameter :: n_sample = 3
597 real(kind=
dp),
parameter :: over_scale = 2d0
599 real(kind=
dp) :: v_low, v_high, vol, r, r_max
609 r_max = max(r_max, r)
611 loss_max(b) = r_max * over_scale
632 real(kind=
dp),
intent(in) :: delta_t
640 integer :: c, b, s, i_part
641 real(kind=
dp) :: over_rate, over_prob, rand_real, rand_geom
648 do i_part = aero_state%apa%n_part, 1, -1
650 aero_data, aero_state, env_state, i_part, 1d0)
657 if (.not. aero_state%aero_sorted%removal_rate_bounds_valid)
then
659 aero_state%aero_sorted%bin_grid, aero_data, env_state, &
660 aero_state%aero_sorted%removal_rate_max)
661 aero_state%aero_sorted%removal_rate_bounds_valid = .true.
664 do c = 1,aero_sorted_n_class(aero_state%aero_sorted)
665 do b = 1,aero_sorted_n_bin(aero_state%aero_sorted)
666 over_rate = aero_state%aero_sorted%removal_rate_max(b)
667 if (delta_t * over_rate <= 0d0) cycle
668 over_prob = 1d0 - exp(-delta_t * over_rate)
671 do s = aero_state%aero_sorted%size_class%inverse(b, c)%n_entry, &
674 aero_state%aero_sorted%size_class%inverse(b, c)%entry(s)
676 aero_data, aero_state, env_state, i_part, 1d0)
680 s = aero_state%aero_sorted%size_class%inverse(b, c)%n_entry + 1
683 if (rand_real <= 0d0)
exit
684 rand_geom = -log(rand_real) / (delta_t * over_rate) + 1d0
685 if (rand_geom >= real(s, kind=
dp))
exit
686 s = s - floor(rand_geom)
692 aero_state%aero_sorted%size_class%inverse(b, c)%entry(s)
694 aero_data, aero_state, env_state, i_part, over_prob)
712 aero_data, aero_state, env_state, i_part, over_prob)
717 real(kind=
dp),
intent(in) :: delta_t
725 integer,
intent(in) :: i_part
727 real(kind=
dp),
intent(in) :: over_prob
729 real(kind=
dp) :: prob, rate, vol, density
730 type(aero_info_t) :: aero_info
732 vol = aero_particle_volume(aero_state%apa%particle(i_part))
733 density = aero_particle_density(aero_state%apa%particle(i_part), aero_data)
735 prob = 1d0 - exp(-delta_t * rate)
737 "particle loss upper bound estimation is too tight: " &
742 aero_info%id = aero_state%apa%particle(i_part)%id
743 aero_info%action = aero_info_dilution
744 aero_info%other_id = 0
758 integer,
intent(in) :: aero_mode_type
771 subroutine spec_file_read_scenario(file, gas_data, aero_data, &
772 read_aero_weight_classes, scenario)
781 logical,
intent(in) :: read_aero_weight_classes
785 character(len=PMC_MAX_FILENAME_LEN) :: sub_filename
787 character(len=SPEC_LINE_MAX_VAR_LEN) :: function_name
833 scenario%temp_time, scenario%temp)
840 scenario%pressure_time, scenario%pressure)
847 scenario%height_time, scenario%height)
853 call spec_file_read_gas_states_times_rates(sub_file, gas_data, &
854 scenario%gas_emission_time, scenario%gas_emission_rate_scale, &
855 scenario%gas_emission)
861 call spec_file_read_gas_states_times_rates(sub_file, gas_data, &
862 scenario%gas_dilution_time, scenario%gas_dilution_rate, &
863 scenario%gas_background)
869 call spec_file_read_aero_dists_times_rates(sub_file, aero_data, &
870 read_aero_weight_classes, scenario%aero_emission_time, &
871 scenario%aero_emission_rate_scale, scenario%aero_emission)
877 call spec_file_read_aero_dists_times_rates(sub_file, aero_data, &
878 read_aero_weight_classes, scenario%aero_dilution_time, &
879 scenario%aero_dilution_rate, scenario%aero_background)
884 if (trim(function_name) ==
'none')
then
886 else if (trim(function_name) ==
'constant')
then
888 else if (trim(function_name) ==
'volume')
then
890 else if (trim(function_name) ==
'drydep')
then
892 else if (trim(function_name) ==
'chamber')
then
894 call spec_file_read_chamber(file, scenario%chamber)
897 "Unknown loss function type: " // trim(function_name))
900 end subroutine spec_file_read_scenario
1011 integer :: total_size, i, n
1030 if (
allocated(val%gas_emission_time))
then
1031 do i = 1,
size(val%gas_emission)
1032 total_size = total_size &
1036 if (
allocated(val%gas_dilution_time))
then
1037 do i = 1,
size(val%gas_background)
1038 total_size = total_size &
1042 if (
allocated(val%aero_emission_time))
then
1043 do i = 1,
size(val%aero_emission)
1044 total_size = total_size &
1048 if (
allocated(val%aero_dilution_time))
then
1049 do i = 1,
size(val%aero_background)
1050 total_size = total_size &
1065 character,
intent(inout) :: buffer(:)
1067 integer,
intent(inout) :: position
1072 integer :: prev_position, i
1074 prev_position = position
1087 val%aero_emission_rate_scale)
1092 if (
allocated(val%gas_emission_time))
then
1093 do i = 1,
size(val%gas_emission)
1097 if (
allocated(val%gas_dilution_time))
then
1098 do i = 1,
size(val%gas_background)
1102 if (
allocated(val%aero_emission_time))
then
1103 do i = 1,
size(val%aero_emission)
1107 if (
allocated(val%aero_dilution_time))
then
1108 do i = 1,
size(val%aero_background)
1124 character,
intent(inout) :: buffer(:)
1126 integer,
intent(inout) :: position
1131 integer :: prev_position, i
1133 prev_position = position
1142 val%gas_emission_rate_scale)
1147 val%aero_emission_rate_scale)
1152 if (
allocated(val%gas_emission))
deallocate(val%gas_emission)
1153 if (
allocated(val%gas_background))
deallocate(val%gas_background)
1154 if (
allocated(val%aero_emission))
deallocate(val%aero_emission)
1155 if (
allocated(val%aero_background))
deallocate(val%aero_background)
1156 if (
allocated(val%gas_emission_time))
then
1157 allocate(val%gas_emission(
size(val%gas_emission_time)))
1158 do i = 1,
size(val%gas_emission)
1162 if (
allocated(val%gas_dilution_time))
then
1163 allocate(val%gas_background(
size(val%gas_dilution_time)))
1164 do i = 1,
size(val%gas_background)
1166 val%gas_background(i))
1169 if (
allocated(val%aero_emission_time))
then
1170 allocate(val%aero_emission(
size(val%aero_emission_time)))
1171 do i = 1,
size(val%aero_emission)
1175 if (
allocated(val%aero_dilution_time))
then
1176 allocate(val%aero_background(
size(val%aero_dilution_time)))
1177 do i = 1,
size(val%aero_background)
1179 val%aero_background(i))