PartMC  2.3.0
run_sect.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 Nicole Riemer and Matthew West
2 ! Copyright (C) Andreas Bott
3 ! Licensed under the GNU General Public License version 2 or (at your
4 ! option) any later version. See the file COPYING for details.
5 
6 !> \file
7 !> The pmc_run_sect module.
8 
9 !> 1D sectional simulation.
10 !!
11 !! Sectional code based on \c coad1d.f by Andreas Bott
12 !! - http://www.meteo.uni-bonn.de/mitarbeiter/ABott/
13 !! - Released under the GPL to Nicole Riemer (personal communication)
14 !! - A. Bott, A flux method for the numerical solution of the
15 !! stochastic collection equation, J. Atmos. Sci. 55, 2284-2293,
16 !! 1998.
18 
19  use pmc_bin_grid
20  use pmc_aero_binned
21  use pmc_util
22  use pmc_aero_dist
23  use pmc_scenario
24  use pmc_env_state
25  use pmc_aero_data
26  use pmc_coag_kernel
27  use pmc_output
28  use pmc_gas_data
29  use pmc_gas_state
30 
31  !> Options controlling the operation of run_sect().
33  !> Final time (s).
34  real(kind=dp) :: t_max
35  !> Timestep for coagulation (s).
36  real(kind=dp) :: del_t
37  !> Output interval (0 disables) (s).
38  real(kind=dp) :: t_output
39  !> Progress interval (0 disables) (s).
40  real(kind=dp) :: t_progress
41  !> Whether to do coagulation.
42  logical :: do_coagulation
43  !> Output prefix.
44  character(len=300) :: prefix
45  !> Type of coagulation kernel.
46  integer :: coag_kernel_type
47  !> UUID of the simulation.
48  character(len=PMC_UUID_LEN) :: uuid
49  end type run_sect_opt_t
50 
51 contains
52 
53 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
54 
55  !> Run a sectional simulation.
56  subroutine run_sect(bin_grid, gas_data, aero_data, aero_dist, &
57  scenario, env_state, run_sect_opt)
58 
59  !> Bin grid.
60  type(bin_grid_t), intent(in) :: bin_grid
61  !> Gas data.
62  type(gas_data_t), intent(in) :: gas_data
63  !> Aerosol data.
64  type(aero_data_t), intent(in) :: aero_data
65  !> Aerosol distribution.
66  type(aero_dist_t), intent(inout) :: aero_dist
67  !> Environment data.
68  type(scenario_t), intent(inout) :: scenario
69  !> Environment state.
70  type(env_state_t), intent(inout) :: env_state
71  !> Options.
72  type(run_sect_opt_t), intent(in) :: run_sect_opt
73 
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
84  type(env_state_t) :: old_env_state
85  type(aero_binned_t) :: aero_binned
86  type(gas_state_t) :: gas_state
87 
88  integer i, j, i_time, num_t, i_summary
89  logical do_output, do_progress
90 
91  call check_time_multiple("t_max", run_sect_opt%t_max, &
92  "del_t", run_sect_opt%del_t)
93  call check_time_multiple("t_output", run_sect_opt%t_output, &
94  "del_t", run_sect_opt%del_t)
95  call check_time_multiple("t_progress", run_sect_opt%t_progress, &
96  "del_t", run_sect_opt%del_t)
97 
98  ! g : spectral mass distribution (mg/cm^3)
99  ! e : droplet mass grid (mg)
100  ! r : droplet radius grid (um)
101  ! log_width : constant grid distance of logarithmic grid
102 
103  if (aero_data%n_spec /= 1) then
104  call die_msg(844211192, &
105  'run_sect() can only use one aerosol species')
106  end if
107 
108  ! output data structure
109  call aero_binned_allocate_size(aero_binned, bin_grid%n_bin, &
110  aero_data%n_spec)
111  aero_binned%vol_conc = 0d0
112  call gas_state_allocate_size(gas_state, gas_data%n_spec)
113 
114  ! mass and radius grid
115  do i = 1,bin_grid%n_bin
116  r(i) = bin_grid%centers(i) * 1d6 ! radius in m to um
117  e(i) = rad2vol(bin_grid%centers(i)) &
118  * aero_data%density(1) * 1d6 ! vol in m^3 to mass in mg
119  end do
120 
121  ! initial mass distribution
122  call aero_binned_add_aero_dist(aero_binned, bin_grid, aero_data, &
123  aero_dist)
124 
125  call courant(bin_grid%n_bin, bin_grid%widths(1), e, ima, c)
126 
127  ! initialize time
128  last_progress_time = 0d0
129  time = 0d0
130  i_summary = 1
131 
132  ! precompute kernel values for all pairs of bins
133  call bin_kernel(bin_grid%n_bin, bin_grid%centers, aero_data, &
134  run_sect_opt%coag_kernel_type, env_state, k_bin)
135  call smooth_bin_kernel(bin_grid%n_bin, k_bin, ck)
136  do i = 1,bin_grid%n_bin
137  do j = 1,bin_grid%n_bin
138  ck(i,j) = ck(i,j) * 1d6 ! m^3/s to cm^3/s
139  end do
140  end do
141 
142  ! multiply kernel with constant timestep and logarithmic grid distance
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)
146  end do
147  end do
148 
149  ! initial output
150  call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
151  last_output_time, do_output)
152  if (do_output) then
153  call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
154  aero_binned, gas_data, gas_state, env_state, i_summary, &
155  time, run_sect_opt%t_output, run_sect_opt%uuid)
156  end if
157 
158  ! main time-stepping loop
159  num_t = nint(run_sect_opt%t_max / run_sect_opt%del_t)
160  call env_state_allocate(old_env_state)
161  do i_time = 1, num_t
162 
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) &
169  / rad2vol(bin_grid%centers)
170  end if
171 
172  time = run_sect_opt%t_max * real(i_time, kind=dp) &
173  / real(num_t, kind=dp)
174 
175  call env_state_copy(env_state, old_env_state)
176  call scenario_update_env_state(scenario, env_state, time)
177  call scenario_update_gas_state(scenario, run_sect_opt%del_t, &
178  env_state, old_env_state, gas_data, gas_state)
179  call scenario_update_aero_binned(scenario, run_sect_opt%del_t, &
180  env_state, old_env_state, bin_grid, aero_data, aero_binned)
181 
182  ! print output
183  call check_event(time, run_sect_opt%del_t, run_sect_opt%t_output, &
184  last_output_time, do_output)
185  if (do_output) then
186  i_summary = i_summary + 1
187  call output_sectional(run_sect_opt%prefix, bin_grid, aero_data, &
188  aero_binned, gas_data, gas_state, env_state, i_summary, &
189  time, run_sect_opt%t_output, run_sect_opt%uuid)
190  end if
191 
192  ! print progress to stdout
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
198  end if
199  end do
200 
201  call env_state_deallocate(old_env_state)
202  call aero_binned_deallocate(aero_binned)
203  call gas_state_deallocate(gas_state)
204 
205  end subroutine run_sect
206 
207 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
208 
209  !> Collision subroutine, exponential approach.
210  subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, &
211  c, ima, g, r, e, ck, ec)
212 
213  integer n_bin
214  real(kind=dp) dt
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)
228 
229  real(kind=dp), parameter :: gmin = 1d-60
230 
231  integer i, i0, i1, j, k, kp
232  real(kind=dp) x0, gsi, gsj, gsk, gk, x1, flux
233 
234  do i = 1,n_bin
235  prod(i) = 0d0
236  ploss(i) = 0d0
237  end do
238 
239  ! lower and upper integration limit i0,i1
240  do i0 = 1,(n_bin - 1)
241  if (g(i0) .gt. gmin) goto 2000
242  end do
243 2000 continue
244  do i1 = (n_bin - 1),1,-1
245  if (g(i1) .gt. gmin) goto 2010
246  end do
247 2010 continue
248 
249  do i = i0,i1
250  do j = i,i1
251  k = ima(i,j) ! k = 0 means that i + j goes nowhere
252  kp = k + 1
253 
254  x0 = ck(i,j) * g(i) * g(j)
255  x0 = min(x0, g(i) * e(j))
256 
257  if (j .ne. k) x0 = min(x0, g(j) * e(i))
258  gsi = x0 / e(j)
259  gsj = x0 / e(i)
260  gsk = gsi + gsj
261 
262  ! loss from positions i, j
263  ploss(i) = ploss(i) + gsi
264  ploss(j) = ploss(j) + gsj
265  g(i) = g(i) - gsi
266  g(j) = g(j) - gsj
267 
268  if (k > 0) then ! do we have a valid bin for the coagulation result?
269  gk = g(k) + gsk
270 
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))))
275  flux = min(flux, gk)
276  g(k) = gk - flux
277  g(kp) = g(kp) + flux
278  ! gain for positions i, j
279  prod(k) = prod(k) + gsk - flux
280  prod(kp) = prod(kp) + flux
281  end if
282  end if
283  end do
284  end do
285 
286  end subroutine coad
287 
288 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
289 
290  !> Determines the Courant number for each bin pair.
291  subroutine courant(n_bin, log_width, e, ima, c)
292 
293  !> Number of bins.
294  integer, intent(in) :: n_bin
295  !> Bin scale factor.
296  real(kind=dp), intent(in) :: log_width
297  !> Droplet mass grid (mg).
298  real(kind=dp), intent(in) :: e(n_bin)
299  !> i + j goes in bin ima(i,j).
300  integer, intent(out) :: ima(n_bin,n_bin)
301  !> Courant number for bin pairs.
302  real(kind=dp), intent(out) :: c(n_bin,n_bin)
303 
304  integer i, j, k, kk
305  real(kind=dp) x0
306 
307  c = 0d0 ! added to avoid uninitialized access errors
308  ima = 0 ! ima(i,j) = 0 means that particles i + j go nowhere
309  do i = 1,n_bin
310  do j = i,n_bin
311  x0 = e(i) + e(j)
312  ! this is basically the same as particle_in_bin(), but that
313  ! is actually slightly different than what was always done
314  ! here
315 
316  ! MW 2011-04-28: I think the above comment no longer
317  ! applies, and we can make this just
318  ! bin_grid_find(). FIXME.
319  k = find_1d(n_bin, e, x0)
320  if (k < n_bin) then
321  k = k + 1
322  if (c(i,j) .lt. 1d0 - 1d-08) then
323  kk = k - 1
324  c(i,j) = log(x0 / e(k-1)) / (3d0 * log_width)
325  else
326  c(i,j) = 0d0
327  kk = k
328  end if
329  ima(i,j) = min(n_bin - 1, kk)
330  end if
331  c(j,i) = c(i,j)
332  ima(j,i) = ima(i,j)
333  end do
334  end do
335 
336  end subroutine courant
337 
338 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
339 
340  !> Smooths kernel values for bin pairs, and halves the self-rate.
341  subroutine smooth_bin_kernel(n_bin, k, k_smooth)
342 
343  !> Number of bins.
344  integer, intent(in) :: n_bin
345  !> Kernel values.
346  real(kind=dp), intent(in) :: k(n_bin,n_bin)
347  !> Smoothed kernel values.
348  real(kind=dp), intent(out) :: k_smooth(n_bin,n_bin)
349 
350  integer i, j, im, ip, jm, jp
351 
352  do i = 1,n_bin
353  do j = 1,n_bin
354  jm = max0(j - 1, 1)
355  im = max0(i - 1, 1)
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)) &
360  + 0.5d0 * k(i,j)
361  if (i .eq. j) then
362  k_smooth(i,j) = 0.5d0 * k_smooth(i,j)
363  end if
364  end do
365  end do
366 
367  end subroutine smooth_bin_kernel
368 
369 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
370 
371 end module pmc_run_sect
subroutine aero_binned_allocate_size(aero_binned, n_bin, n_spec)
Allocate an aero_binned_t of the given size.
Definition: aero_binned.F90:60
The scenario_t structure and associated subroutines.
Definition: scenario.F90:9
1D sectional simulation.
Definition: run_sect.F90:17
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
real(kind=dp) elemental function rad2vol(r)
Convert radius (m) to volume (m^3).
Definition: util.F90:274
The aero_dist_t structure and associated subroutines.
Definition: aero_dist.F90:18
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
Definition: util.F90:577
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
subroutine aero_binned_deallocate(aero_binned)
Free internal memory in an aero_binned_t structure.
Definition: aero_binned.F90:79
subroutine scenario_update_env_state(scenario, env_state, time)
Update time-dependent contents of the environment. scenario_init_env_state() should have been called ...
Definition: scenario.F90:301
The gas_data_t structure and associated subroutines.
Definition: gas_data.F90:9
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...
Definition: util.F90:368
A complete aerosol distribution, consisting of several modes.
Definition: aero_dist.F90:33
Common utility subroutines.
Definition: util.F90:9
Current environment state.
Definition: env_state.F90:26
subroutine env_state_deallocate(env_state)
Free all storage.
Definition: env_state.F90:78
Options controlling the operation of run_sect().
Definition: run_sect.F90:32
subroutine gas_state_deallocate(gas_state)
Free all storage.
Definition: gas_state.F90:63
subroutine run_sect(bin_grid, gas_data, aero_data, aero_dist, scenario, env_state, run_sect_opt)
Run a sectional simulation.
Definition: run_sect.F90:56
The gas_state_t structure and associated subroutines.
Definition: gas_state.F90:9
Scenario data.
Definition: scenario.F90:51
subroutine scenario_update_gas_state(scenario, delta_t, env_state, old_env_state, gas_data, gas_state)
Do gas emissions and background dilution.
Definition: scenario.F90:375
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
subroutine env_state_copy(env_from, env_to)
env_to = env_from
Definition: env_state.F90:138
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
The aero_binned_t structure and associated subroutines.
Definition: aero_binned.F90:9
subroutine check_event(time, timestep, interval, last_time, do_event)
Computes whether an event is scheduled to take place.
Definition: util.F90:399
Current state of the gas mixing ratios in the system.
Definition: gas_state.F90:26
subroutine coad(n_bin, dt, taug, taup, taul, tauu, prod, ploss, c, ima, g, r, e, ck, ec)
Collision subroutine, exponential approach.
Definition: run_sect.F90:210
Constant gas data.
Definition: gas_data.F90:29
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.
Definition: scenario.F90:504
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.
Definition: gas_state.F90:48
subroutine courant(n_bin, log_width, e, ima, c)
Determines the Courant number for each bin pair.
Definition: run_sect.F90:291
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.
Definition: env_state.F90:56
subroutine smooth_bin_kernel(n_bin, k, k_smooth)
Smooths kernel values for bin pairs, and halves the self-rate.
Definition: run_sect.F90:341
Aerosol material properties and associated data.
Definition: aero_data.F90:40
Write data in NetCDF format.
Definition: output.F90:68
Generic coagulation kernel.
Definition: coag_kernel.F90:9
Aerosol number and volume distributions stored per bin.
Definition: aero_binned.F90:33
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.
Definition: output.F90:676