34 real(kind=dp) :: t_max
36 real(kind=dp) :: del_t
38 real(kind=dp) :: t_output
40 real(kind=dp) :: t_progress
42 logical :: do_coagulation
44 character(len=300) :: prefix
46 integer :: coag_kernel_type
48 character(len=PMC_UUID_LEN) :: uuid
56 subroutine run_sect(bin_grid, gas_data, aero_data, aero_dist, &
57 scenario, env_state, run_sect_opt)
74 real(kind=dp) c(bin_grid%n_bin,bin_grid%n_bin)
75 integer ima(bin_grid%n_bin,bin_grid%n_bin)
76 real(kind=dp) g(bin_grid%n_bin), r(bin_grid%n_bin), e(bin_grid%n_bin)
77 real(kind=dp) k_bin(bin_grid%n_bin,bin_grid%n_bin)
78 real(kind=dp) ck(bin_grid%n_bin,bin_grid%n_bin)
79 real(kind=dp) ec(bin_grid%n_bin,bin_grid%n_bin)
80 real(kind=dp) taug(bin_grid%n_bin), taup(bin_grid%n_bin)
81 real(kind=dp) taul(bin_grid%n_bin), tauu(bin_grid%n_bin)
82 real(kind=dp) prod(bin_grid%n_bin), ploss(bin_grid%n_bin)
83 real(kind=dp) time, last_output_time, last_progress_time
88 integer i, j, i_time, num_t, i_summary
89 logical do_output, do_progress
92 "del_t", run_sect_opt%del_t)
94 "del_t", run_sect_opt%del_t)
96 "del_t", run_sect_opt%del_t)
103 if (aero_data%n_spec /= 1)
then
105 'run_sect() can only use one aerosol species')
111 aero_binned%vol_conc = 0d0
115 do i = 1,bin_grid%n_bin
116 r(i) = bin_grid%centers(i) * 1d6
117 e(i) =
rad2vol(bin_grid%centers(i)) &
118 * aero_data%density(1) * 1d6
125 call
courant(bin_grid%n_bin, bin_grid%widths(1), e, ima, c)
128 last_progress_time = 0d0
133 call
bin_kernel(bin_grid%n_bin, bin_grid%centers, aero_data, &
134 run_sect_opt%coag_kernel_type, env_state, k_bin)
136 do i = 1,bin_grid%n_bin
137 do j = 1,bin_grid%n_bin
138 ck(i,j) = ck(i,j) * 1d6
143 do i = 1,bin_grid%n_bin
144 do j = 1,bin_grid%n_bin
145 ck(i,j) = ck(i,j) * run_sect_opt%del_t * bin_grid%widths(i)
150 call
check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
151 last_output_time, do_output)
154 aero_binned, gas_data, gas_state, env_state, i_summary, &
155 time, run_sect_opt%t_output, run_sect_opt%uuid)
159 num_t = nint(run_sect_opt%t_max / run_sect_opt%del_t)
163 if (run_sect_opt%do_coagulation)
then
164 g = aero_binned%vol_conc(:,1) * aero_data%density(1)
165 call
coad(bin_grid%n_bin, run_sect_opt%del_t, taug, taup, taul, &
166 tauu, prod, ploss, c, ima, g, r, e, ck, ec)
167 aero_binned%vol_conc(:,1) = g / aero_data%density(1)
168 aero_binned%num_conc = aero_binned%vol_conc(:,1) &
172 time = run_sect_opt%t_max *
real(i_time, kind=dp) &
173 /
real(num_t, kind=dp)
178 env_state, old_env_state, gas_data, gas_state)
180 env_state, old_env_state, bin_grid, aero_data, aero_binned)
183 call
check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
184 last_output_time, do_output)
186 i_summary = i_summary + 1
188 aero_binned, gas_data, gas_state, env_state, i_summary, &
189 time, run_sect_opt%t_output, run_sect_opt%uuid)
193 call
check_event(time, run_sect_opt%del_t, run_sect_opt%t_progress, &
194 last_progress_time, do_progress)
195 if (do_progress)
then
196 write(*,
'(a6,a8)')
'step',
'time'
197 write(*,
'(i6,f8.1)') i_time, time
210 subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, &
211 c, ima, g, r, e, ck, ec)
215 real(kind=dp) taug(n_bin)
216 real(kind=dp) taup(n_bin)
217 real(kind=dp) taul(n_bin)
218 real(kind=dp) tauu(n_bin)
219 real(kind=dp) prod(n_bin)
220 real(kind=dp) ploss(n_bin)
221 real(kind=dp) c(n_bin,n_bin)
222 integer ima(n_bin,n_bin)
223 real(kind=dp) g(n_bin)
224 real(kind=dp) r(n_bin)
225 real(kind=dp) e(n_bin)
226 real(kind=dp) ck(n_bin,n_bin)
227 real(kind=dp) ec(n_bin,n_bin)
229 real(kind=dp),
parameter :: gmin = 1d-60
231 integer i, i0, i1, j, k, kp
232 real(kind=dp) x0, gsi, gsj, gsk, gk, x1, flux
240 do i0 = 1,(n_bin - 1)
241 if (g(i0) .gt. gmin) goto 2000
244 do i1 = (n_bin - 1),1,-1
245 if (g(i1) .gt. gmin) goto 2010
254 x0 = ck(i,j) * g(i) * g(j)
255 x0 = min(x0, g(i) * e(j))
257 if (j .ne. k) x0 = min(x0, g(j) * e(i))
263 ploss(i) = ploss(i) + gsi
264 ploss(j) = ploss(j) + gsj
271 if (gk .gt. gmin)
then
272 x1 = log(g(kp) / gk + 1d-60)
273 flux = gsk / x1 * (exp(0.5d0 * x1) &
274 - exp(x1 * (0.5d0 - c(i,j))))
279 prod(k) = prod(k) + gsk - flux
280 prod(kp) = prod(kp) + flux
291 subroutine courant(n_bin, log_width, e, ima, c)
294 integer,
intent(in) :: n_bin
296 real(kind=dp),
intent(in) :: log_width
298 real(kind=dp),
intent(in) :: e(n_bin)
300 integer,
intent(out) :: ima(n_bin,n_bin)
302 real(kind=dp),
intent(out) :: c(n_bin,n_bin)
322 if (c(i,j) .lt. 1d0 - 1d-08)
then
324 c(i,j) = log(x0 / e(k-1)) / (3d0 * log_width)
329 ima(i,j) = min(n_bin - 1, kk)
344 integer,
intent(in) :: n_bin
346 real(kind=dp),
intent(in) :: k(n_bin,n_bin)
348 real(kind=dp),
intent(out) :: k_smooth(n_bin,n_bin)
350 integer i, j, im, ip, jm, jp
356 jp = min0(j + 1, n_bin)
357 ip = min0(i + 1, n_bin)
358 k_smooth(i,j) = 0.125d0 * (k(i,jm) + k(im,j) &
359 + k(ip,j) + k(i,jp)) &
362 k_smooth(i,j) = 0.5d0 * k_smooth(i,j)
subroutine aero_binned_allocate_size(aero_binned, n_bin, n_spec)
Allocate an aero_binned_t of the given size.
The scenario_t structure and associated subroutines.
The aero_data_t structure and associated subroutines.
subroutine die_msg(code, error_msg)
Error immediately.
real(kind=dp) elemental function rad2vol(r)
Convert radius (m) to volume (m^3).
The aero_dist_t structure and associated subroutines.
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
The env_state_t structure and associated subroutines.
subroutine aero_binned_deallocate(aero_binned)
Free internal memory in an aero_binned_t structure.
subroutine scenario_update_env_state(scenario, env_state, time)
Update time-dependent contents of the environment. scenario_init_env_state() should have been called ...
The gas_data_t structure and associated subroutines.
subroutine check_time_multiple(first_name, first_time, second_name, second_time)
Check that the first time interval is close to an integer multiple of the second, and warn if it is n...
A complete aerosol distribution, consisting of several modes.
Common utility subroutines.
Current environment state.
subroutine env_state_deallocate(env_state)
Free all storage.
Options controlling the operation of run_sect().
subroutine gas_state_deallocate(gas_state)
Free all storage.
subroutine run_sect(bin_grid, gas_data, aero_data, aero_dist, scenario, env_state, run_sect_opt)
Run a sectional simulation.
The gas_state_t structure and associated subroutines.
subroutine scenario_update_gas_state(scenario, delta_t, env_state, old_env_state, gas_data, gas_state)
Do gas emissions and background dilution.
1D grid, either logarithmic or linear.
subroutine env_state_copy(env_from, env_to)
env_to = env_from
The bin_grid_t structure and associated subroutines.
The aero_binned_t structure and associated subroutines.
subroutine check_event(time, timestep, interval, last_time, do_event)
Computes whether an event is scheduled to take place.
Current state of the gas mixing ratios in the system.
subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, c, ima, g, r, e, ck, ec)
Collision subroutine, exponential approach.
subroutine scenario_update_aero_binned(scenario, delta_t, env_state, old_env_state, bin_grid, aero_data, aero_binned)
Do emissions and background dilution from the environment for a binned aerosol distribution.
subroutine bin_kernel(n_bin, bin_r, aero_data, coag_kernel_type, env_state, k)
Computes an array of kernel values for each bin pair. k(i,j) is the kernel value at the centers of bi...
subroutine gas_state_allocate_size(gas_state, n_spec)
Allocate storage for gas species of the given size.
subroutine courant(n_bin, log_width, e, ima, c)
Determines the Courant number for each bin pair.
subroutine aero_binned_add_aero_dist(aero_binned, bin_grid, aero_data, aero_dist)
Add an aero_dist_t to an aero_binned_t.
subroutine env_state_allocate(env_state)
Allocate an empty environment.
subroutine smooth_bin_kernel(n_bin, k, k_smooth)
Smooths kernel values for bin pairs, and halves the self-rate.
Aerosol material properties and associated data.
Write data in NetCDF format.
Generic coagulation kernel.
Aerosol number and volume distributions stored per bin.
subroutine output_sectional(prefix, bin_grid, aero_data, aero_binned, gas_data, gas_state, env_state, index, time, del_t, uuid)
Write the current sectional data.