PartMC  2.5.0
aero_state.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2018 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_aero_state module.
7 
8 !> The aero_state_t structure and assocated subroutines.
10 
12  use pmc_aero_sorted
14  use pmc_bin_grid
15  use pmc_aero_data
17  use pmc_aero_dist
18  use pmc_util
19  use pmc_rand
20  use pmc_aero_binned
21  use pmc_mpi
22  use pmc_spec_file
23  use pmc_aero_info
25  use pmc_aero_weight
27 #ifdef PMC_USE_MPI
28  use mpi
29 #endif
30 
31  !> MPI tag for mixing particles between processes.
32  integer, parameter :: aero_state_tag_mix = 4987
33  !> MPI tag for gathering between processes.
34  integer, parameter :: aero_state_tag_gather = 4988
35  !> MPI tag for scattering between processes.
36  integer, parameter :: aero_state_tag_scatter = 4989
37 
38  !> Single flat weighting scheme.
39  integer, parameter :: aero_state_weight_none = 1
40  !> Single flat weighting scheme.
41  integer, parameter :: aero_state_weight_flat = 2
42  !> Power-law weighting scheme.
43  integer, parameter :: aero_state_weight_power = 3
44  !> Coupled number/mass weighting scheme.
45  integer, parameter :: aero_state_weight_nummass = 4
46  !> Flat weighting by source.
47  integer, parameter :: aero_state_weight_flat_source = 5
48  !> Power-law weighting by source.
49  integer, parameter :: aero_state_weight_power_source = 6
50  !> Coupled number/mass weighting by source.
51  integer, parameter :: aero_state_weight_nummass_source = 7
52 
53  !> The current collection of aerosol particles.
54  !!
55  !! The particles in \c aero_state_t are stored in a single flat
56  !! array (the \c apa data member), with a sorting into size bins and
57  !! weight groups/classes possibly stored in the \c aero_sorted data
58  !! member (if \c valid_sort is true).
59  !!
60  !! Every time we remove particles we keep track of the particle ID
61  !! and the action performed in the aero_info_array_t structure. This
62  !! is typically cleared each time we output data to disk.
64  !> Linear array of particles.
65  type(aero_particle_array_t) :: apa
66  !> Sorting of particles into size bins and weight groups/classes.
67  type(aero_sorted_t) :: aero_sorted
68  !> Whether the \c aero_sorted is a correct sorting.
69  logical :: valid_sort
70  !> Weighting functions.
71  type(aero_weight_array_t) :: awa
72  !> Ideal number of computational particles, per weight group and class.
73  real(kind=dp), allocatable :: n_part_ideal(:, :)
74  !> Information on removed particles.
75  type(aero_info_array_t) :: aero_info_array
76  end type aero_state_t
77 
78 contains
79 
80 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
81 
82  !> Return the current number of particles.
83  elemental integer function aero_state_n_part(aero_state)
84 
85  !> Aerosol state.
86  type(aero_state_t), intent(in) :: aero_state
87 
89 
90  end function aero_state_n_part
91 
92 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
93 
94  !> Read the specification for a weighting type from a spec file.
95  subroutine spec_file_read_aero_state_weighting_type(file, weighting_type, &
96  exponent)
97 
98  !> Spec file.
99  type(spec_file_t), intent(inout) :: file
100  !> Aerosol weighting scheme.
101  integer, intent(out) :: weighting_type
102  !> Exponent for power-law weighting (only used if \c weight_type
103  !> is \c AERO_STATE_WEIGHT_POWER).
104  real(kind=dp), intent(out) :: exponent
105 
106  character(len=SPEC_LINE_MAX_VAR_LEN) :: weighting_name
107 
108  !> \page input_format_weight_type Input File Format: Type of aerosol size distribution weighting functions.
109  !!
110  !! The weighting function is specified by the parameters:
111  !! - \b weight_type (string): the type of weighting function ---
112  !! must be one of: \c flat for flat weighting, \c flat_source for
113  !! flat weighting by source, \c power for power weighting,
114  !! \c power_source for power source weighting, \c nummass for number
115  !! and mass weighting, and \c nummass_source for number and mass
116  !! weighting by source. If \c weight_type is \c power or \c
117  !! power_source then the next parameter must also be provided:
118  !! - \b weighting_exponent (real): the exponent for \c power or
119  !! \c power_source. Setting the \c exponent to 0 is equivalent to no
120  !! weighting, while setting the exponent negative uses more
121  !! computational particles at larger diameters and setting the
122  !! exponent positive uses more computaitonal partilces at smaller
123  !! diameters; in practice exponents between 0 and -3 are most useful.
124  !!
125  !! See also:
126  !! - \ref spec_file_format --- the input file text format
127 
128  call spec_file_read_string(file, 'weight_type', weighting_name)
129  if (trim(weighting_name) == 'flat') then
130  weighting_type = aero_state_weight_flat
131  elseif (trim(weighting_name) == 'power') then
132  weighting_type = aero_state_weight_power
133  call spec_file_read_real(file, 'weighting_exponent', exponent)
134  elseif (trim(weighting_name) == 'nummass') then
135  weighting_type = aero_state_weight_nummass
136  elseif (trim(weighting_name) == 'flat_source') then
137  weighting_type = aero_state_weight_flat_source
138  elseif (trim(weighting_name) == 'power_source') then
139  weighting_type = aero_state_weight_power_source
140  call spec_file_read_real(file, 'weighting_exponent', exponent)
141  elseif (trim(weighting_name) == 'nummass_source') then
142  weighting_type = aero_state_weight_nummass_source
143  else
144  call spec_file_die_msg(920321729, file, &
145  "Unknown weighting type: " // trim(weighting_name))
146  end if
147 
148  end subroutine spec_file_read_aero_state_weighting_type
149 
150 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
151 
152  !> Copies weighting information for an \c aero_state.
153  subroutine aero_state_copy_weight(aero_state_from, aero_state_to)
155  !> Reference aerosol.
156  type(aero_state_t), intent(in) :: aero_state_from
157  !> Already allocated.
158  type(aero_state_t), intent(inout) :: aero_state_to
159 
160  aero_state_to%awa = aero_state_from%awa
161 
162  end subroutine aero_state_copy_weight
163 
164 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
165 
166  !> Sets the weighting functions for an \c aero_state.
167  subroutine aero_state_set_weight(aero_state, aero_data, weight_type, &
168  exponent)
170  !> Aerosol to set the weights on.
171  type(aero_state_t), intent(inout) :: aero_state
172  !> Aerosol data.
173  type(aero_data_t), intent(in) :: aero_data
174  !> Type of weighting scheme to use.
175  integer, intent(in) :: weight_type
176  !> Exponent for power-law weighting (only used if \c weight_type
177  !> is \c AERO_STATE_WEIGHT_POWER).
178  real(kind=dp), intent(in), optional :: exponent
179 
180  aero_state%valid_sort = .false.
181  select case(weight_type)
184  call aero_weight_array_set_flat(aero_state%awa, 1)
186  call assert_msg(656670336, present(exponent), &
187  "exponent parameter required for AERO_STATE_WEIGHT_POWER")
188  call aero_weight_array_set_power(aero_state%awa, 1, exponent)
190  call aero_weight_array_set_nummass(aero_state%awa, 1)
192  call aero_weight_array_set_flat(aero_state%awa, &
193  aero_data_n_source(aero_data))
195  call assert_msg(102143848, present(exponent), &
196  "exponent parameter required for AERO_STATE_WEIGHT_POWER")
197  call aero_weight_array_set_power(aero_state%awa, &
198  aero_data_n_source(aero_data), exponent)
200  call aero_weight_array_set_nummass(aero_state%awa, &
201  aero_data_n_source(aero_data))
202  case default
203  call die_msg(969076992, "unknown weight_type: " &
204  // trim(integer_to_string(weight_type)))
205  end select
206 
207  end subroutine aero_state_set_weight
208 
209 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
210 
211  !> Set the ideal number of particles to the given value. The \c
212  !> aero_state%%awa must be already set correctly.
213  subroutine aero_state_set_n_part_ideal(aero_state, n_part)
215  !> Aerosol state (with \c aero_state%%awa set).
216  type(aero_state_t), intent(inout) :: aero_state
217  !> Ideal total number of particles.
218  real(kind=dp), intent(in) :: n_part
219 
220  integer :: n_group, n_class
221 
222  n_group = aero_weight_array_n_group(aero_state%awa)
223  n_class = aero_weight_array_n_class(aero_state%awa)
224  if (allocated(aero_state%n_part_ideal)) then
225  deallocate(aero_state%n_part_ideal)
226  end if
227  allocate(aero_state%n_part_ideal(n_group, n_class))
228  aero_state%n_part_ideal = n_part / real(n_group * n_class, kind=dp)
229 
230  end subroutine aero_state_set_n_part_ideal
231 
232 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
233 
234  !> Determine the appropriate weight class for a source.
235  integer function aero_state_weight_class_for_source(aero_state, source)
237  !> Aerosol state.
238  type(aero_state_t), intent(in) :: aero_state
239  !> Source to find the class for.
240  integer, intent(in) :: source
241 
242  integer :: n_class
243 
244  call assert(932390238, source >= 1)
245  n_class = aero_weight_array_n_class(aero_state%awa)
246  ! we are either using i_class = i_source or always i_class = n_class = 1
247  if (n_class > 1) then
248  call assert(765048788, source <= n_class)
250  else
252  end if
253 
255 
256 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
257 
258  !> Returns the total number of particles in an aerosol distribution.
259  integer function aero_state_total_particles(aero_state, i_group, i_class)
261  !> Aerosol state.
262  type(aero_state_t), intent(in) :: aero_state
263  !> Weight group.
264  integer, optional, intent(in) :: i_group
265  !> Weight class.
266  integer, optional, intent(in) :: i_class
267 
268  integer :: i_part
269 
270  if (present(i_group)) then
271  call assert(908743823, present(i_class))
272  if (aero_state%valid_sort) then
275  aero_state%aero_sorted%group_class%inverse(i_group, i_class))
276  else
277  ! FIXME: should we just sort?
279  do i_part = 1,aero_state_n_part(aero_state)
280  if ((aero_state%apa%particle(i_part)%weight_group == i_group) &
281  .and. &
282  (aero_state%apa%particle(i_part)%weight_class == i_class)) &
283  then
285  end if
286  end do
287  end if
288  else
290  end if
291 
292  end function aero_state_total_particles
293 
294 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
295 
296  !> Returns the total number of particles across all processes.
297  integer function aero_state_total_particles_all_procs(aero_state, i_group, &
298  i_class)
300  !> Aerosol state.
301  type(aero_state_t), intent(in) :: aero_state
302  !> Weight group.
303  integer, optional, intent(in) :: i_group
304  !> Weight class.
305  integer, optional, intent(in) :: i_class
306 
308  aero_state_total_particles(aero_state, i_group, i_class), &
310 
312 
313 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
314 
315  !> Resets an aero_state to have zero particles per bin.
316  subroutine aero_state_zero(aero_state)
318  !> State to zero.
319  type(aero_state_t), intent(inout) :: aero_state
320 
321  integer :: i, n_bin
322 
323  call aero_particle_array_zero(aero_state%apa)
324  aero_state%valid_sort = .false.
325  call aero_info_array_zero(aero_state%aero_info_array)
326 
327  end subroutine aero_state_zero
328 
329 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
330 
331  !> Add the given particle.
332  subroutine aero_state_add_particle(aero_state, aero_particle, aero_data, &
333  allow_resort)
335  !> Aerosol state.
336  type(aero_state_t), intent(inout) :: aero_state
337  !> Particle to add.
338  type(aero_particle_t), intent(in) :: aero_particle
339  !> Aerosol data.
340  type(aero_data_t), intent(in) :: aero_data
341  !> Whether to allow a resort due to the add.
342  logical, optional, intent(in) :: allow_resort
343 
344  if (aero_state%valid_sort) then
345  call aero_sorted_add_particle(aero_state%aero_sorted, aero_state%apa, &
346  aero_particle, aero_data, allow_resort)
347  else
348  call aero_particle_array_add_particle(aero_state%apa, aero_particle)
349  end if
350 
351  end subroutine aero_state_add_particle
352 
353 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
354 
355  !> Remove the given particle without recording it.
356  subroutine aero_state_remove_particle_no_info(aero_state, i_part)
358  !> Aerosol state.
359  type(aero_state_t), intent(inout) :: aero_state
360  !> Index of particle to remove.
361  integer, intent(in) :: i_part
362 
363  if (aero_state%valid_sort) then
364  call aero_sorted_remove_particle(aero_state%aero_sorted, &
365  aero_state%apa, i_part)
366  else
367  call aero_particle_array_remove_particle(aero_state%apa, i_part)
368  end if
369 
371 
372 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
373 
374  !> Remove the given particle and record the removal.
375  subroutine aero_state_remove_particle_with_info(aero_state, i_part, &
376  aero_info)
378  !> Aerosol state.
379  type(aero_state_t), intent(inout) :: aero_state
380  !> Index of particle to remove.
381  integer, intent(in) :: i_part
382  !> Removal info.
383  type(aero_info_t), intent(in) :: aero_info
384 
385  call aero_state_remove_particle_no_info(aero_state, i_part)
386  call aero_info_array_add_aero_info(aero_state%aero_info_array, &
387  aero_info)
388 
390 
391 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
392 
393  !> Remove the given particle and possibly record the removal.
394  subroutine aero_state_remove_particle(aero_state, i_part, &
395  record_removal, aero_info)
397  !> Aerosol state.
398  type(aero_state_t), intent(inout) :: aero_state
399  !> Index of particle to remove.
400  integer, intent(in) :: i_part
401  !> Whether to record the removal in the aero_info_array.
402  logical, intent(in) :: record_removal
403  !> Removal info.
404  type(aero_info_t), intent(in) :: aero_info
405 
406  if (record_removal) then
407  call aero_state_remove_particle_with_info(aero_state, i_part, &
408  aero_info)
409  else
410  call aero_state_remove_particle_no_info(aero_state, i_part)
411  end if
412 
413  end subroutine aero_state_remove_particle
414 
415 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
416 
417  !> Remove a randomly chosen particle from the given bin and return
418  !> it.
419  subroutine aero_state_remove_rand_particle_from_bin(aero_state, &
420  i_bin, i_class, aero_particle)
422  !> Aerosol state.
423  type(aero_state_t), intent(inout) :: aero_state
424  !> Bin number to remove particle from.
425  integer, intent(in) :: i_bin
426  !> Weight class to remove particle from.
427  integer, intent(in) :: i_class
428  !> Removed particle.
429  type(aero_particle_t), intent(inout) :: aero_particle
430 
431  integer :: i_entry, i_part
432 
433  call assert(742996300, aero_state%valid_sort)
434  call assert(392182617, integer_varray_n_entry( &
435  aero_state%aero_sorted%size_class%inverse(i_bin, i_class)) > 0)
437  aero_state%aero_sorted%size_class%inverse(i_bin, i_class)))
438  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
439  i_class)%entry(i_entry)
440  aero_particle = aero_state%apa%particle(i_part)
441  call aero_state_remove_particle_no_info(aero_state, i_part)
442 
444 
445 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
446 
447  !> Add copies or remove a particle, with a given mean number of
448  !> resulting particles.
449  !!
450  !! The particle number \c i_part is either removed, or zero or more
451  !! copies are added, with a random number of copies with the given
452  !! mean \c n_part_mean. The final number of particles is either
453  !! <tt>floor(n_part_mean)</tt> or <tt>ceiling(n_part_mean)</tt>,
454  !! chosen randomly so the mean is \c n_part_mean.
455  subroutine aero_state_dup_particle(aero_state, aero_data, i_part, &
456  n_part_mean, random_weight_group)
458  !> Aerosol state.
459  type(aero_state_t), intent(inout) :: aero_state
460  !> Aerosol data.
461  type(aero_data_t), intent(in) :: aero_data
462  !> Particle number.
463  integer, intent(in) :: i_part
464  !> Mean number of resulting particles.
465  real(kind=dp), intent(in) :: n_part_mean
466  !> Whether particle copies should be placed in a randomly chosen
467  !> weight group.
468  logical, optional, intent(in) :: random_weight_group
469 
470  integer :: n_copies, i_dup, new_group
471  type(aero_particle_t) :: new_aero_particle
472  type(aero_info_t) :: aero_info
473 
474  n_copies = prob_round(n_part_mean)
475  if (n_copies == 0) then
476  aero_info%id = aero_state%apa%particle(i_part)%id
477  aero_info%action = aero_info_weight
478  aero_info%other_id = 0
479  call aero_state_remove_particle_with_info(aero_state, &
480  i_part, aero_info)
481  elseif (n_copies > 1) then
482  do i_dup = 1,(n_copies - 1)
483  new_aero_particle = aero_state%apa%particle(i_part)
484  call aero_particle_new_id(new_aero_particle)
485  if (present(random_weight_group)) then
486  if (random_weight_group) then
487  new_group &
488  = aero_weight_array_rand_group(aero_state%awa, &
489  aero_state%apa%particle(i_part)%weight_class, &
490  aero_particle_radius(aero_state%apa%particle(i_part), &
491  aero_data))
492  call aero_particle_set_weight(new_aero_particle, new_group)
493  end if
494  end if
495  call aero_state_add_particle(aero_state, new_aero_particle, &
496  aero_data)
497  end do
498  end if
499 
500  end subroutine aero_state_dup_particle
501 
502 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
503 
504  !> The number concentration of a single particle (m^{-3}).
505  real(kind=dp) function aero_state_particle_num_conc(aero_state, &
506  aero_particle, aero_data)
508  !> Aerosol state containing the particle.
509  type(aero_state_t), intent(in) :: aero_state
510  !> Aerosol particle.
511  type(aero_particle_t), intent(in) :: aero_particle
512  !> Aerosol data.
513  type(aero_data_t), intent(in) :: aero_data
514 
516  = aero_weight_array_num_conc(aero_state%awa, aero_particle, &
517  aero_data)
518 
519  end function aero_state_particle_num_conc
520 
521 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
522 
523  !> Save the correct number concentrations for later use by
524  !> aero_state_reweight().
525  subroutine aero_state_num_conc_for_reweight(aero_state, aero_data, &
526  reweight_num_conc)
528  !> Aerosol state.
529  type(aero_state_t), intent(in) :: aero_state
530  !> Aerosol data.
531  type(aero_data_t), intent(in) :: aero_data
532  !> Number concentrations for later use by aero_state_reweight().
533  real(kind=dp), intent(out) &
534  :: reweight_num_conc(aero_state_n_part(aero_state))
535 
536  integer :: i_part
537 
538  do i_part = 1,aero_state_n_part(aero_state)
539  reweight_num_conc(i_part) &
540  = aero_weight_array_single_num_conc(aero_state%awa, &
541  aero_state%apa%particle(i_part), aero_data)
542  end do
543 
544  end subroutine aero_state_num_conc_for_reweight
545 
546 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
547 
548  !> Reweight all particles after their constituent volumes have been
549  !> altered.
550  !!
551  !! The pattern for use should be like:
552  !! <pre>
553  !! call aero_state_num_conc_for_reweight(aero_state, aero_data,
554  !! reweight_num_conc)
555  !! ... alter particle species volumes in aero_state ...
556  !! call aero_state_reweight(aero_state, aero_data, reweight_num_conc)
557  !! </pre>
558  subroutine aero_state_reweight(aero_state, aero_data, reweight_num_conc)
560  !> Aerosol state.
561  type(aero_state_t), intent(inout) :: aero_state
562  !> Aerosol data.
563  type(aero_data_t), intent(in) :: aero_data
564  !> Number concentrations previously computed by
565  !> aero_state_num_conc_for_reweight().
566  real(kind=dp), intent(in) &
567  :: reweight_num_conc(aero_state_n_part(aero_state))
568 
569  integer :: i_part, i_group, i_class
570  real(kind=dp) :: n_part_old(size(aero_state%awa%weight, 1), &
571  size(aero_state%awa%weight, 2))
572  real(kind=dp) :: n_part_new(size(aero_state%awa%weight, 1), &
573  size(aero_state%awa%weight, 2))
574  real(kind=dp) :: old_num_conc, new_num_conc, n_part_mean
575 
576  ! find average number of new particles in each weight group, if
577  ! weight is not changed
578  n_part_old = 0d0
579  n_part_new = 0d0
580  do i_part = 1,aero_state_n_part(aero_state)
581  old_num_conc = reweight_num_conc(i_part)
582  new_num_conc = aero_weight_array_single_num_conc(aero_state%awa, &
583  aero_state%apa%particle(i_part), aero_data)
584  n_part_mean = old_num_conc / new_num_conc
585  i_group = aero_state%apa%particle(i_part)%weight_group
586  i_class = aero_state%apa%particle(i_part)%weight_class
587  n_part_new(i_group, i_class) = n_part_new(i_group, i_class) &
588  + n_part_mean
589  n_part_old(i_group, i_class) = n_part_old(i_group, i_class) + 1d0
590  end do
591 
592  ! alter weight to leave the number of computational particles
593  ! per weight bin unchanged
594  do i_group = 1,size(aero_state%awa%weight, 1)
595  do i_class = 1,size(aero_state%awa%weight, 2)
596  if (n_part_old(i_group, i_class) == 0d0) cycle
597  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), &
598  n_part_new(i_group, i_class) / n_part_old(i_group, i_class))
599  end do
600  end do
601 
602  ! work backwards so any additions and removals will only affect
603  ! particles that we've already dealt with
604  do i_part = aero_state_n_part(aero_state),1,-1
605  old_num_conc = reweight_num_conc(i_part)
606  new_num_conc &
607  = aero_weight_array_single_num_conc(aero_state%awa, &
608  aero_state%apa%particle(i_part), aero_data)
609  n_part_mean = old_num_conc / new_num_conc
610  call aero_state_dup_particle(aero_state, aero_data, i_part, &
611  n_part_mean)
612  end do
613 
614  end subroutine aero_state_reweight
615 
616 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
617 
618  !> <tt>aero_state += aero_state_delta</tt>, including combining the
619  !> weights, so the new concentration is the weighted average of the
620  !> two concentrations.
621  subroutine aero_state_add(aero_state, aero_state_delta, aero_data)
623  !> Aerosol state.
624  type(aero_state_t), intent(inout) :: aero_state
625  !> Increment.
626  type(aero_state_t), intent(in) :: aero_state_delta
627  !> Aerosol data.
628  type(aero_data_t), intent(in) :: aero_data
629 
630  call aero_state_add_particles(aero_state, aero_state_delta, aero_data)
631  call aero_weight_array_combine(aero_state%awa, aero_state_delta%awa)
632 
633  end subroutine aero_state_add
634 
635 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
636 
637  !> <tt>aero_state += aero_state_delta</tt>, with the weight
638  !> of \c aero_state left unchanged, so the new concentration is the
639  !> sum of the two concentrations, computed with \c aero_state%%awa.
640  subroutine aero_state_add_particles(aero_state, aero_state_delta, aero_data)
642  !> Aerosol state.
643  type(aero_state_t), intent(inout) :: aero_state
644  !> Increment.
645  type(aero_state_t), intent(in) :: aero_state_delta
646  !> Aerosol data.
647  type(aero_data_t), intent(in) :: aero_data
648 
649  integer :: i_part, i_bin
650 
651  do i_part = 1,aero_state_delta%apa%n_part
652  call aero_state_add_particle(aero_state, &
653  aero_state_delta%apa%particle(i_part), aero_data)
654  end do
655  call aero_info_array_add(aero_state%aero_info_array, &
656  aero_state_delta%aero_info_array)
657 
658  end subroutine aero_state_add_particles
659 
660 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
661 
662  !> Change the weight if necessary to ensure that the addition of
663  !> about \c n_add computational particles will give the correct
664  !> final particle number.
665  subroutine aero_state_prepare_weight_for_add(aero_state, aero_data, &
666  i_group, i_class, n_add, allow_doubling, allow_halving)
668  !> Aero state to add to.
669  type(aero_state_t), intent(inout) :: aero_state
670  !> Aerosol data.
671  type(aero_data_t), intent(in) :: aero_data
672  !> Weight group number to add to.
673  integer, intent(in) :: i_group
674  !> Weight class number to add to.
675  integer, intent(in) :: i_class
676  !> Approximate number of particles to be added at current weighting.
677  real(kind=dp), intent(in) :: n_add
678  !> Whether to allow doubling of the population.
679  logical, intent(in) :: allow_doubling
680  !> Whether to allow halving of the population.
681  logical, intent(in) :: allow_halving
682 
683  integer :: global_n_part, n_group, n_class
684  real(kind=dp) :: mean_n_part, n_part_new, weight_ratio
685  real(kind=dp) :: n_part_ideal_local_group
686 
687  n_group = aero_weight_array_n_group(aero_state%awa)
688  n_class = aero_weight_array_n_class(aero_state%awa)
689  global_n_part = aero_state_total_particles_all_procs(aero_state, &
690  i_group, i_class)
691  mean_n_part = real(global_n_part, kind=dp) / real(pmc_mpi_size(), kind=dp)
692  n_part_new = mean_n_part + n_add
693  if (n_part_new == 0d0) return
694  n_part_ideal_local_group = aero_state%n_part_ideal(i_group, i_class) &
695  / real(pmc_mpi_size(), kind=dp)
696  if ((n_part_new < n_part_ideal_local_group / 2d0) &
697  .or. (n_part_new > n_part_ideal_local_group * 2d0)) &
698  then
699  weight_ratio = n_part_new / n_part_ideal_local_group
700  call aero_state_scale_weight(aero_state, aero_data, i_group, &
701  i_class, weight_ratio, allow_doubling, allow_halving)
702  end if
703 
704  end subroutine aero_state_prepare_weight_for_add
705 
706 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
707 
708  !> Generates a Poisson sample of an \c aero_dist, adding to \c
709  !> aero_state, with the given sample proportion.
710  subroutine aero_state_add_aero_dist_sample(aero_state, aero_data, &
711  aero_dist, sample_prop, create_time, allow_doubling, allow_halving, &
712  n_part_add)
714  !> Aero state to add to.
715  type(aero_state_t), intent(inout) :: aero_state
716  !> Aero data values.
717  type(aero_data_t), intent(in) :: aero_data
718  !> Distribution to sample.
719  type(aero_dist_t), intent(in) :: aero_dist
720  !> Fraction to sample (1).
721  real(kind=dp), intent(in) :: sample_prop
722  !> Creation time for new particles (s).
723  real(kind=dp), intent(in) :: create_time
724  !> Whether to allow doubling of the population.
725  logical, intent(in) :: allow_doubling
726  !> Whether to allow halving of the population.
727  logical, intent(in) :: allow_halving
728  !> Number of particles added.
729  integer, intent(out), optional :: n_part_add
730 
731  real(kind=dp) :: n_samp_avg, radius, total_vol
732  real(kind=dp) :: vols(aero_data_n_spec(aero_data))
733  integer :: n_samp, i_mode, i_samp, i_group, i_class, n_group, n_class
734  type(aero_particle_t) :: aero_particle
735 
736  n_group = size(aero_state%awa%weight, 1)
737  n_class = size(aero_state%awa%weight, 2)
738  if (present(n_part_add)) then
739  n_part_add = 0
740  end if
741  do i_group = 1,n_group
742  do i_mode = 1,aero_dist_n_mode(aero_dist)
743  i_class = aero_state_weight_class_for_source(aero_state, &
744  aero_dist%mode(i_mode)%source)
745 
746  ! adjust weight if necessary
747  n_samp_avg = sample_prop * aero_mode_number(aero_dist%mode(i_mode), &
748  aero_state%awa%weight(i_group, i_class))
749  call aero_state_prepare_weight_for_add(aero_state, aero_data, &
750  i_group, i_class, n_samp_avg, allow_doubling, allow_halving)
751  if (n_samp_avg == 0d0) cycle
752 
753  ! sample particles
754  n_samp_avg = sample_prop * aero_mode_number(aero_dist%mode(i_mode), &
755  aero_state%awa%weight(i_group, i_class))
756  n_samp = rand_poisson(n_samp_avg)
757  if (present(n_part_add)) then
758  n_part_add = n_part_add + n_samp
759  end if
760  do i_samp = 1,n_samp
761  call aero_particle_zero(aero_particle, aero_data)
762  call aero_mode_sample_radius(aero_dist%mode(i_mode), aero_data, &
763  aero_state%awa%weight(i_group, i_class), radius)
764  total_vol = aero_data_rad2vol(aero_data, radius)
765  call aero_mode_sample_vols(aero_dist%mode(i_mode), total_vol, &
766  vols)
767  call aero_particle_set_vols(aero_particle, vols)
768  call aero_particle_new_id(aero_particle)
769  call aero_particle_set_weight(aero_particle, i_group, i_class)
770  call aero_particle_set_create_time(aero_particle, create_time)
771  call aero_particle_set_source(aero_particle, &
772  aero_dist%mode(i_mode)%source)
773  call aero_state_add_particle(aero_state, aero_particle, aero_data)
774  end do
775  end do
776  end do
777 
778  end subroutine aero_state_add_aero_dist_sample
779 
780 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
781 
782  !> Choose a random particle from the aero_state.
783  subroutine aero_state_rand_particle(aero_state, i_part)
785  !> Original state.
786  type(aero_state_t), intent(in) :: aero_state
787  !> Chosen random particle number.
788  integer, intent(out) :: i_part
789 
790  call assert(950725003, aero_state_n_part(aero_state) > 0)
791  i_part = pmc_rand_int(aero_state_n_part(aero_state))
792 
793  end subroutine aero_state_rand_particle
794 
795 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
796 
797  !> Generates a random sample by removing particles from
798  !> aero_state_from and adding them to aero_state_to, which must be
799  !> already allocated (and should have its weight set).
800  !!
801  !! None of the weights are altered by this sampling, making this the
802  !! equivalent of aero_state_add_particles().
803  subroutine aero_state_sample_particles(aero_state_from, aero_state_to, &
804  aero_data, sample_prob, removal_action)
806  !> Original state.
807  type(aero_state_t), intent(inout) :: aero_state_from
808  !> Destination state.
809  type(aero_state_t), intent(inout) :: aero_state_to
810  !> Aerosol data.
811  type(aero_data_t), intent(in) :: aero_data
812  !> Probability of sampling each particle.
813  real(kind=dp), intent(in) :: sample_prob
814  !> Action for removal (see pmc_aero_info module for action
815  !> parameters). Set to AERO_INFO_NONE to not log removal.
816  integer, intent(in) :: removal_action
817 
818  integer :: n_transfer, i_transfer, i_part
819  logical :: do_add, do_remove
820  real(kind=dp) :: num_conc_from, num_conc_to
821  type(aero_info_t) :: aero_info
822 
823  call assert(721006962, (sample_prob >= 0d0) .and. (sample_prob <= 1d0))
824  call aero_state_zero(aero_state_to)
825  call aero_state_copy_weight(aero_state_from, aero_state_to)
826  n_transfer = rand_binomial(aero_state_total_particles(aero_state_from), &
827  sample_prob)
828  i_transfer = 0
829  do while (i_transfer < n_transfer)
830  if (aero_state_total_particles(aero_state_from) <= 0) exit
831  call aero_state_rand_particle(aero_state_from, i_part)
832  num_conc_from = aero_weight_array_num_conc(aero_state_from%awa, &
833  aero_state_from%apa%particle(i_part), aero_data)
834  num_conc_to = aero_weight_array_num_conc(aero_state_to%awa, &
835  aero_state_from%apa%particle(i_part), aero_data)
836 
837  if (num_conc_to == num_conc_from) then ! add and remove
838  do_add = .true.
839  do_remove = .true.
840  elseif (num_conc_to < num_conc_from) then ! always add, maybe remove
841  do_add = .true.
842  do_remove = .false.
843  if (pmc_random() < num_conc_to / num_conc_from) then
844  do_remove = .true.
845  end if
846  else ! always remove, maybe add
847  do_add = .false.
848  if (pmc_random() < num_conc_from / num_conc_to) then
849  do_add = .true.
850  end if
851  do_remove = .true.
852  end if
853  if (do_add) then
854  call aero_state_add_particle(aero_state_to, &
855  aero_state_from%apa%particle(i_part), aero_data)
856  end if
857  if (do_remove) then
858  if (removal_action /= aero_info_none) then
859  aero_info%id = aero_state_from%apa%particle(i_part)%id
860  aero_info%action = removal_action
861  call aero_state_remove_particle_with_info(aero_state_from, &
862  i_part, aero_info)
863  else
864  call aero_state_remove_particle_no_info(aero_state_from, &
865  i_part)
866  end if
867  i_transfer = i_transfer + 1
868  end if
869  end do
870 
871  end subroutine aero_state_sample_particles
872 
873 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
874 
875  !> Generates a random sample by removing particles from
876  !> aero_state_from and adding them to aero_state_to, transfering
877  !> weight as well. This is the equivalent of aero_state_add().
878  subroutine aero_state_sample(aero_state_from, aero_state_to, &
879  aero_data, sample_prob, removal_action)
881  !> Original state.
882  type(aero_state_t), intent(inout) :: aero_state_from
883  !> Destination state (previous contents will be lost).
884  type(aero_state_t), intent(inout) :: aero_state_to
885  !> Aerosol data.
886  type(aero_data_t), intent(in) :: aero_data
887  !> Probability of sampling each particle.
888  real(kind=dp), intent(in) :: sample_prob
889  !> Action for removal (see pmc_aero_info module for action
890  !> parameters). Set to AERO_INFO_NONE to not log removal.
891  integer, intent(in) :: removal_action
892 
893  integer :: n_transfer, i_transfer, i_part
894  logical :: do_add, do_remove, overwrite_to
895  real(kind=dp) :: num_conc_from, num_conc_to
896  type(aero_info_t) :: aero_info
897 
898  call assert(393205561, (sample_prob >= 0d0) .and. (sample_prob <= 1d0))
899  call aero_state_zero(aero_state_to)
900  call aero_state_copy_weight(aero_state_from, aero_state_to)
901  call aero_weight_array_normalize(aero_state_to%awa)
902  n_transfer = rand_binomial(aero_state_total_particles(aero_state_from), &
903  sample_prob)
904  do i_transfer = 1,n_transfer
905  if (aero_state_total_particles(aero_state_from) <= 0) exit
906  call aero_state_rand_particle(aero_state_from, i_part)
907 
908  call aero_state_add_particle(aero_state_to, &
909  aero_state_from%apa%particle(i_part), aero_data)
910  if (removal_action /= aero_info_none) then
911  aero_info%id = aero_state_from%apa%particle(i_part)%id
912  aero_info%action = removal_action
913  call aero_state_remove_particle_with_info(aero_state_from, &
914  i_part, aero_info)
915  else
916  call aero_state_remove_particle_no_info(aero_state_from, &
917  i_part)
918  end if
919  end do
920  overwrite_to = .true.
921  call aero_weight_array_shift(aero_state_from%awa, aero_state_to%awa, &
922  sample_prob, overwrite_to)
923 
924  end subroutine aero_state_sample
925 
926 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
927 
928  !> Create binned number and mass arrays.
929  subroutine aero_state_to_binned(bin_grid, aero_data, aero_state, &
930  aero_binned)
932  !> Bin grid.
933  type(bin_grid_t), intent(in) :: bin_grid
934  !> Aerosol data.
935  type(aero_data_t), intent(in) :: aero_data
936  !> Aerosol state.
937  type(aero_state_t), intent(in) :: aero_state
938  !> Binned distributions.
939  type(aero_binned_t), intent(inout) :: aero_binned
940 
941  integer :: i_part, i_bin
942  real(kind=dp) :: factor
943 
944  aero_binned%num_conc = 0d0
945  aero_binned%vol_conc = 0d0
946  do i_part = 1,aero_state_n_part(aero_state)
947  i_bin = bin_grid_find(bin_grid, &
948  aero_particle_radius(aero_state%apa%particle(i_part), aero_data))
949  if ((i_bin < 1) .or. (i_bin > bin_grid_size(bin_grid))) then
950  call warn_msg(980232449, "particle ID " &
951  // trim(integer_to_string(aero_state%apa%particle(i_part)%id)) &
952  // " outside of bin_grid, discarding")
953  else
954  factor = aero_weight_array_num_conc(aero_state%awa, &
955  aero_state%apa%particle(i_part), aero_data) &
956  / bin_grid%widths(i_bin)
957  aero_binned%vol_conc(i_bin,:) = aero_binned%vol_conc(i_bin,:) &
958  + aero_state%apa%particle(i_part)%vol * factor
959  aero_binned%num_conc(i_bin) = aero_binned%num_conc(i_bin) &
960  + factor
961  end if
962  end do
963 
964  end subroutine aero_state_to_binned
965 
966 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
967 
968  !> Returns the IDs of all particles.
969  function aero_state_ids(aero_state)
971  !> Aerosol state.
972  type(aero_state_t), intent(in) :: aero_state
973 
974  !> Return value.
975  integer :: aero_state_ids(aero_state_n_part(aero_state))
976 
977  integer :: i_part
978 
979  do i_part = 1,aero_state_n_part(aero_state)
980  aero_state_ids(i_part) = aero_state%apa%particle(i_part)%id
981  end do
982 
983  end function aero_state_ids
984 
985 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
986 
987  !> Returns the diameters of all particles.
988  function aero_state_diameters(aero_state, aero_data, include, exclude)
990  !> Aerosol state.
991  type(aero_state_t), intent(in) :: aero_state
992  !> Aerosol data.
993  type(aero_data_t), intent(in) :: aero_data
994  !> Species names to include in the diameter.
995  character(len=*), optional, intent(in) :: include(:)
996  !> Species names to exclude in the diameter.
997  character(len=*), optional, intent(in) :: exclude(:)
998 
999  !> Return diameters array (m).
1000  real(kind=dp) :: aero_state_diameters(aero_state_n_part(aero_state))
1001 
1002  !> Per-particle volume of included components
1003  real(kind=dp) :: volumes(aero_state_n_part(aero_state))
1004 
1005  volumes = aero_state_volumes(aero_state, aero_data, include, exclude)
1006  aero_state_diameters = rad2diam(sphere_vol2rad(volumes))
1007 
1008  end function aero_state_diameters
1009 
1010 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1011 
1012  !> Returns the dry diameters of all particles.
1013  function aero_state_dry_diameters(aero_state, aero_data)
1015  !> Aerosol state.
1016  type(aero_state_t), intent(in) :: aero_state
1017  !> Aerosol data.
1018  type(aero_data_t), intent(in) :: aero_data
1019 
1020  !> Return value (m).
1021  real(kind=dp) :: aero_state_dry_diameters(aero_state_n_part(aero_state))
1022 
1023  aero_state_dry_diameters = aero_particle_dry_diameter( &
1024  aero_state%apa%particle(1:aero_state_n_part(aero_state)), aero_data)
1025 
1026  end function aero_state_dry_diameters
1027 
1028 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1029 
1030  !> Returns the mobility diameters of all particles.
1031  function aero_state_mobility_diameters(aero_state, aero_data, env_state)
1033  !> Aerosol state.
1034  type(aero_state_t), intent(in) :: aero_state
1035  !> Aerosol data.
1036  type(aero_data_t), intent(in) :: aero_data
1037  !> Environment state.
1038  type(env_state_t), intent(in) :: env_state
1039 
1040  !> Return value (m).
1041  real(kind=dp) :: aero_state_mobility_diameters( &
1042  aero_state_n_part(aero_state))
1043 
1044  integer :: i_part
1045 
1046  do i_part = 1,aero_state_n_part(aero_state)
1047  aero_state_mobility_diameters(i_part) &
1049  aero_state%apa%particle(i_part), aero_data, env_state)
1050  end do
1051 
1052  end function aero_state_mobility_diameters
1053 
1054 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1055 
1056  !> Returns the volumes of all particles.
1057  !!
1058  !! If \c include is specified then only those species are included
1059  !! in computing the volumes. If \c exclude is specified then all
1060  !! species except those species are included. If both \c include and
1061  !! \c exclude arguments are specified then only those species in \c
1062  !! include but not in \c exclude are included.
1063  function aero_state_volumes(aero_state, aero_data, include, exclude)
1065  !> Aerosol state.
1066  type(aero_state_t), intent(in) :: aero_state
1067  !> Aerosol data.
1068  type(aero_data_t), optional, intent(in) :: aero_data
1069  !> Species names to include in the mass.
1070  character(len=*), optional, intent(in) :: include(:)
1071  !> Species names to exclude in the mass.
1072  character(len=*), optional, intent(in) :: exclude(:)
1073 
1074  !> Return volumes array (m^3).
1075  real(kind=dp) :: aero_state_volumes(aero_state_n_part(aero_state))
1076 
1077  logical, allocatable :: use_species(:)
1078  integer :: i_name, i_spec
1079 
1080  if ((.not. present(include)) .and. (.not. present(exclude))) then
1081  aero_state_volumes = aero_particle_volume( &
1082  aero_state%apa%particle(1:aero_state_n_part(aero_state)))
1083  else
1084  call assert_msg(599558703, present(aero_data), &
1085  "must provide 'aero_data' if using 'include' or 'exclude'")
1086  allocate(use_species(aero_data_n_spec(aero_data)))
1087  if (present(include)) then
1088  use_species = .false.
1089  do i_name = 1, size(include)
1090  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1091  call assert_msg(111852070, i_spec > 0, &
1092  "unknown species: " // trim(include(i_name)))
1093  use_species(i_spec) = .true.
1094  end do
1095  else
1096  use_species = .true.
1097  end if
1098  if (present(exclude)) then
1099  do i_name = 1, size(exclude)
1100  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1101  call assert_msg(182075590, i_spec > 0, &
1102  "unknown species: " // trim(exclude(i_name)))
1103  use_species(i_spec) = .false.
1104  end do
1105  end if
1106  aero_state_volumes = 0d0
1107  do i_spec = 1,aero_data_n_spec(aero_data)
1108  if (use_species(i_spec)) then
1109  aero_state_volumes = aero_state_volumes &
1111  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1112  i_spec)
1113  end if
1114  end do
1115  end if
1116 
1117  end function aero_state_volumes
1118 
1119 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1120 
1121  !> Returns the masses of all particles.
1122  !!
1123  !! If \c include is specified then only those species are included
1124  !! in computing the masses. If \c exclude is specified then all
1125  !! species except those species are included. If both \c include and
1126  !! \c exclude arguments are specified then only those species in \c
1127  !! include but not in \c exclude are included.
1128  function aero_state_masses(aero_state, aero_data, include, exclude)
1130  !> Aerosol state.
1131  type(aero_state_t), intent(in) :: aero_state
1132  !> Aerosol data.
1133  type(aero_data_t), intent(in) :: aero_data
1134  !> Species names to include in the mass.
1135  character(len=*), optional, intent(in) :: include(:)
1136  !> Species names to exclude in the mass.
1137  character(len=*), optional, intent(in) :: exclude(:)
1138 
1139  !> Return masses array (kg).
1140  real(kind=dp) :: aero_state_masses(aero_state_n_part(aero_state))
1141 
1142  logical :: use_species(aero_data_n_spec(aero_data))
1143  integer :: i_name, i_spec
1144 
1145  if ((.not. present(include)) .and. (.not. present(exclude))) then
1146  aero_state_masses = aero_particle_mass( &
1147  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1148  aero_data)
1149  else
1150  if (present(include)) then
1151  use_species = .false.
1152  do i_name = 1, size(include)
1153  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1154  call assert_msg(963163690, i_spec > 0, &
1155  "unknown species: " // trim(include(i_name)))
1156  use_species(i_spec) = .true.
1157  end do
1158  else
1159  use_species = .true.
1160  end if
1161  if (present(exclude)) then
1162  do i_name = 1, size(exclude)
1163  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1164  call assert_msg(950847713, i_spec > 0, &
1165  "unknown species: " // trim(exclude(i_name)))
1166  use_species(i_spec) = .false.
1167  end do
1168  end if
1169  aero_state_masses = 0d0
1170  do i_spec = 1,aero_data_n_spec(aero_data)
1171  if (use_species(i_spec)) then
1172  aero_state_masses = aero_state_masses &
1174  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1175  i_spec, aero_data)
1176  end if
1177  end do
1178  end if
1179 
1180  end function aero_state_masses
1181 
1182 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1183 
1184  !> Returns the number concentrations of all particles.
1185  function aero_state_num_concs(aero_state, aero_data)
1187  !> Aerosol state.
1188  type(aero_state_t), intent(in) :: aero_state
1189  !> Aerosol data.
1190  type(aero_data_t), intent(in) :: aero_data
1191 
1192  !> Return number concentrations array (m^{-3}).
1193  real(kind=dp) :: aero_state_num_concs(aero_state_n_part(aero_state))
1194 
1195  integer :: i_part
1196 
1197  do i_part = 1,aero_state_n_part(aero_state)
1198  aero_state_num_concs(i_part) &
1199  = aero_state_particle_num_conc(aero_state, &
1200  aero_state%apa%particle(i_part), aero_data)
1201  end do
1202 
1203  end function aero_state_num_concs
1204 
1205 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1206 
1207  !> Returns the total number concentration.
1208  real(kind=dp) function aero_state_total_num_conc(aero_state, aero_data)
1210  !> Aerosol state.
1211  type(aero_state_t), intent(in) :: aero_state
1212  !> Aerosol data.
1213  type(aero_data_t), intent(in) :: aero_data
1214 
1215  integer :: i_part
1216 
1218  do i_part = 1,aero_state_n_part(aero_state)
1220  + aero_state_particle_num_conc(aero_state, &
1221  aero_state%apa%particle(i_part), aero_data)
1222  end do
1223 
1224  end function aero_state_total_num_conc
1225 
1226 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1227 
1228  !> Returns the mass-entropies of all particles.
1229  !!
1230  !! If \c include is specified then only those species are included
1231  !! in computing the entropy. If \c exclude is specified then all
1232  !! species except those species are included. If both \c include and
1233  !! \c exclude arguments are specified then only those species in \c
1234  !! include but not in \c exclude are included. If \c group is
1235  !! specified then the species are divided into two sets, given by
1236  !! those in the group and those not in the group. The entropy is
1237  !! then computed using the total mass of each set.
1238  function aero_state_mass_entropies(aero_state, aero_data, include, exclude, &
1239  group)
1241  !> Aerosol state.
1242  type(aero_state_t), intent(in) :: aero_state
1243  !> Aerosol data.
1244  type(aero_data_t), intent(in) :: aero_data
1245  !> Species names to include in the mass.
1246  character(len=*), optional :: include(:)
1247  !> Species names to exclude in the mass.
1248  character(len=*), optional :: exclude(:)
1249  !> Species names to group together.
1250  character(len=*), optional :: group(:)
1251 
1252  !> Return value.
1253  real(kind=dp) :: aero_state_mass_entropies(aero_state_n_part(aero_state))
1254 
1255  logical :: use_species(aero_data_n_spec(aero_data))
1256  logical :: group_species(aero_data_n_spec(aero_data))
1257  integer :: i_name, i_spec, i_part
1258  real(kind=dp) :: group_mass, non_group_mass, mass
1259 
1260  if (present(include)) then
1261  use_species = .false.
1262  do i_name = 1, size(include)
1263  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1264  call assert_msg(890212002, i_spec > 0, &
1265  "unknown species: " // trim(include(i_name)))
1266  use_species(i_spec) = .true.
1267  end do
1268  else
1269  use_species = .true.
1270  end if
1271  if (present(exclude)) then
1272  do i_name = 1, size(exclude)
1273  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1274  call assert_msg(859945006, i_spec > 0, &
1275  "unknown species: " // trim(exclude(i_name)))
1276  use_species(i_spec) = .false.
1277  end do
1278  end if
1279  if (present(group)) then
1280  group_species = .false.
1281  do i_name = 1, size(group)
1282  i_spec = aero_data_spec_by_name(aero_data, group(i_name))
1283  call assert_msg(376359046, i_spec > 0, &
1284  "unknown species: " // trim(group(i_name)))
1285  group_species(i_spec) = .true.
1286  end do
1287  do i_part = 1,aero_state_n_part(aero_state)
1288  group_mass = 0d0
1289  non_group_mass = 0d0
1290  do i_spec = 1,aero_data_n_spec(aero_data)
1291  if (use_species(i_spec)) then
1292  mass = aero_particle_species_mass( &
1293  aero_state%apa%particle(i_part), i_spec, aero_data)
1294  if (group_species(i_spec)) then
1295  group_mass = group_mass + mass
1296  else
1297  non_group_mass = non_group_mass + mass
1298  end if
1299  end if
1300  end do
1301  aero_state_mass_entropies(i_part) &
1302  = entropy([group_mass, non_group_mass])
1303  end do
1304  else
1305  do i_part = 1,aero_state_n_part(aero_state)
1306  aero_state_mass_entropies(i_part) = entropy(pack( &
1307  aero_particle_species_masses(aero_state%apa%particle(i_part), &
1308  aero_data), use_species))
1309  end do
1310  end if
1311 
1312  end function aero_state_mass_entropies
1313 
1314 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1315 
1316  !> Returns the mixing state metrics of the population.
1317  !!
1318  !! If \c include is specified then only those species are included
1319  !! in computing the entropies. If \c exclude is specified then all
1320  !! species except those species are included. If both \c include and
1321  !! \c exclude arguments are specified then only those species in \c
1322  !! include but not in \c exclude are included. If \c group is
1323  !! specified then the species are divided into two sets, given by
1324  !! those in the group and those not in the group. The entropies are
1325  !! then computed using the total mass of each set.
1326  subroutine aero_state_mixing_state_metrics(aero_state, aero_data, d_alpha, &
1327  d_gamma, chi, include, exclude, group)
1329  !> Aerosol state.
1330  type(aero_state_t), intent(in) :: aero_state
1331  !> Aerosol data.
1332  type(aero_data_t), intent(in) :: aero_data
1333  !> Average particle diversity.
1334  real(kind=dp), intent(out) :: d_alpha
1335  !> Bulk diversity.
1336  real(kind=dp), intent(out) :: d_gamma
1337  !> Mixing state index.
1338  real(kind=dp), intent(out) :: chi
1339  !> Species names to include in the mass.
1340  character(len=*), optional :: include(:)
1341  !> Species names to exclude in the mass.
1342  character(len=*), optional :: exclude(:)
1343  !> Species names to group together.
1344  character(len=*), optional :: group(:)
1345 
1346  real(kind=dp), allocatable :: entropies(:), entropies_of_avg_part(:)
1347  real(kind=dp), allocatable :: masses(:), num_concs(:), &
1348  num_concs_of_avg_part(:), masses_of_avg_part(:)
1349  type(aero_state_t) :: aero_state_averaged
1350  type(bin_grid_t) :: avg_bin_grid
1351 
1352  ! per-particle properties
1353  num_concs = aero_state_num_concs(aero_state, aero_data)
1354  masses = aero_state_masses(aero_state, aero_data, include, exclude)
1355  entropies = aero_state_mass_entropies(aero_state, aero_data, &
1356  include, exclude, group)
1357 
1358  d_alpha = exp(sum(entropies * masses * num_concs) &
1359  / sum(masses * num_concs))
1360 
1361  ! per-particle properties of averaged particles
1362  call bin_grid_make(avg_bin_grid, bin_grid_type_log, 1, 1d-30, 1d10)
1363  aero_state_averaged = aero_state
1364  call aero_state_bin_average_comp(aero_state_averaged, avg_bin_grid, &
1365  aero_data)
1366  num_concs_of_avg_part = aero_state_num_concs(aero_state_averaged, &
1367  aero_data)
1368  masses_of_avg_part = aero_state_masses(aero_state_averaged, &
1369  aero_data, include, exclude)
1370  entropies_of_avg_part = aero_state_mass_entropies(aero_state_averaged, &
1371  aero_data, include, exclude, group)
1372 
1373  d_gamma = exp(sum(entropies_of_avg_part * masses_of_avg_part &
1374  * num_concs_of_avg_part) &
1375  / sum(masses_of_avg_part * num_concs_of_avg_part))
1376 
1377  chi = (d_alpha - 1) / (d_gamma - 1)
1378 
1379  end subroutine aero_state_mixing_state_metrics
1380 
1381 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1382 
1383  !> Returns the approximate critical relative humidity for all particles (1).
1384  function aero_state_approx_crit_rel_humids(aero_state, aero_data, env_state)
1386  !> Aerosol state.
1387  type(aero_state_t), intent(in) :: aero_state
1388  !> Aerosol data.
1389  type(aero_data_t), intent(in) :: aero_data
1390  !> Environment state.
1391  type(env_state_t), intent(in) :: env_state
1392 
1393  !> Return value.
1394  real(kind=dp) &
1395  :: aero_state_approx_crit_rel_humids(aero_state_n_part(aero_state))
1396 
1397  integer :: i_part
1398 
1399  do i_part = 1,aero_state_n_part(aero_state)
1400  aero_state_approx_crit_rel_humids(i_part) = &
1402  aero_state%apa%particle(i_part), aero_data, env_state)
1403  end do
1404 
1406 
1407 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1408 
1409  !> Returns the critical relative humidity for all particles (1).
1410  function aero_state_crit_rel_humids(aero_state, aero_data, env_state)
1412  !> Aerosol state.
1413  type(aero_state_t), intent(in) :: aero_state
1414  !> Aerosol data.
1415  type(aero_data_t), intent(in) :: aero_data
1416  !> Environment state.
1417  type(env_state_t), intent(in) :: env_state
1418 
1419  !> Return value.
1420  real(kind=dp) :: aero_state_crit_rel_humids(aero_state_n_part(aero_state))
1421 
1422  integer :: i_part
1423 
1424  do i_part = 1,aero_state_n_part(aero_state)
1425  aero_state_crit_rel_humids(i_part) = aero_particle_crit_rel_humid( &
1426  aero_state%apa%particle(i_part), aero_data, env_state)
1427  end do
1428 
1429  end function aero_state_crit_rel_humids
1430 
1431 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1432 
1433  !> Does the same thing as aero_state_to_bin() but based on dry radius.
1434  subroutine aero_state_to_binned_dry(bin_grid, aero_data, aero_state, &
1435  aero_binned)
1437  !> Bin grid.
1438  type(bin_grid_t), intent(in) :: bin_grid
1439  !> Aerosol data.
1440  type(aero_data_t), intent(in) :: aero_data
1441  !> Aerosol state.
1442  type(aero_state_t), intent(in) :: aero_state
1443  !> Binned distributions.
1444  type(aero_binned_t), intent(inout) :: aero_binned
1445 
1446  integer :: i_part, i_bin
1447  real(kind=dp) :: factor
1448 
1449  aero_binned%num_conc = 0d0
1450  aero_binned%vol_conc = 0d0
1451  do i_part = 1,aero_state_n_part(aero_state)
1452  i_bin = bin_grid_find(bin_grid, &
1453  aero_particle_solute_radius(aero_state%apa%particle(i_part), &
1454  aero_data))
1455  if ((i_bin < 1) .or. (i_bin > bin_grid_size(bin_grid))) then
1456  call warn_msg(503871022, "particle ID " &
1457  // trim(integer_to_string(aero_state%apa%particle(i_part)%id)) &
1458  // " outside of bin_grid, discarding")
1459  else
1460  factor = aero_weight_array_num_conc(aero_state%awa, &
1461  aero_state%apa%particle(i_part), aero_data) &
1462  / bin_grid%widths(i_bin)
1463  aero_binned%vol_conc(i_bin,:) = aero_binned%vol_conc(i_bin,:) &
1464  + aero_state%apa%particle(i_part)%vol * factor
1465  aero_binned%num_conc(i_bin) = aero_binned%num_conc(i_bin) &
1466  + factor
1467  end if
1468  end do
1469 
1470  end subroutine aero_state_to_binned_dry
1471 
1472 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1473 
1474  !> Doubles number of particles in the given weight group.
1475  subroutine aero_state_double(aero_state, aero_data, i_group, i_class)
1477  !> Aerosol state.
1478  type(aero_state_t), intent(inout) :: aero_state
1479  !> Aerosol data.
1480  type(aero_data_t), intent(in) :: aero_data
1481  !> Weight group to double.
1482  integer, intent(in) :: i_group
1483  !> Weight class to double.
1484  integer, intent(in) :: i_class
1485 
1486  integer :: i_part
1487  type(aero_particle_t) :: aero_particle
1488 
1489  do i_part = 1,aero_state_n_part(aero_state)
1490  if ((aero_state%apa%particle(i_part)%weight_group == i_group) &
1491  .and. (aero_state%apa%particle(i_part)%weight_class == i_class)) &
1492  then
1493  aero_particle = aero_state%apa%particle(i_part)
1494  call aero_particle_new_id(aero_particle)
1495  call aero_state_add_particle(aero_state, aero_particle, aero_data)
1496  end if
1497  end do
1498  aero_state%valid_sort = .false.
1499  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), 0.5d0)
1500 
1501  end subroutine aero_state_double
1502 
1503 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1504 
1505  !> Remove approximately half of the particles in the given weight group.
1506  subroutine aero_state_halve(aero_state, i_group, i_class)
1508  !> Aerosol state.
1509  type(aero_state_t), intent(inout) :: aero_state
1510  !> Weight group to halve.
1511  integer, intent(in) :: i_group
1512  !> Weight class to halve.
1513  integer, intent(in) :: i_class
1514 
1515  integer :: i_part
1516  type(aero_info_t) :: aero_info
1517 
1518  do i_part = aero_state_n_part(aero_state),1,-1
1519  if ((aero_state%apa%particle(i_part)%weight_group == i_group) &
1520  .and. (aero_state%apa%particle(i_part)%weight_class == i_class)) &
1521  then
1522  if (pmc_random() < 0.5d0) then
1523  aero_info%id = aero_state%apa%particle(i_part)%id
1524  aero_info%action = aero_info_halved
1525  call aero_state_remove_particle_with_info(aero_state, i_part, &
1526  aero_info)
1527  end if
1528  end if
1529  end do
1530  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), 2d0)
1531 
1532  end subroutine aero_state_halve
1533 
1534 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1535 
1536  !> Double or halve the particle population in each weight group to
1537  !> maintain close to \c n_part_ideal particles per process,
1538  !> allocated equally amongst the weight groups.
1539  subroutine aero_state_rebalance(aero_state, aero_data, allow_doubling, &
1540  allow_halving, initial_state_warning)
1542  !> Aerosol state.
1543  type(aero_state_t), intent(inout) :: aero_state
1544  !> Aerosol data.
1545  type(aero_data_t), intent(in) :: aero_data
1546  !> Whether to allow doubling of the population.
1547  logical, intent(in) :: allow_doubling
1548  !> Whether to allow halving of the population.
1549  logical, intent(in) :: allow_halving
1550  !> Whether to warn due to initial state doubling/halving.
1551  logical, intent(in) :: initial_state_warning
1552 
1553  integer :: i_group, i_class, n_group, n_class, global_n_part
1554 
1555  n_group = size(aero_state%awa%weight, 1)
1556  n_class = size(aero_state%awa%weight, 2)
1557 
1558  ! if we have less than half the maximum number of particles then
1559  ! double until we fill up the array
1560  if (allow_doubling) then
1561  do i_group = 1,n_group
1562  do i_class = 1,n_class
1563  global_n_part &
1564  = aero_state_total_particles_all_procs(aero_state, i_group, &
1565  i_class)
1566  do while ((real(global_n_part, kind=dp) &
1567  < aero_state%n_part_ideal(i_group, i_class) / 2d0) &
1568  .and. (global_n_part > 0))
1569  if (initial_state_warning) then
1570  call warn_msg(716882783, "doubling particles in initial " &
1571  // "condition")
1572  end if
1573  call aero_state_double(aero_state, aero_data, i_group, i_class)
1574  global_n_part &
1575  = aero_state_total_particles_all_procs(aero_state, &
1576  i_group, i_class)
1577  end do
1578  end do
1579  end do
1580  end if
1581 
1582  ! same for halving if we have too many particles
1583  if (allow_halving) then
1584  do i_group = 1,n_group
1585  do i_class = 1,n_class
1586  do while (real(aero_state_total_particles_all_procs(aero_state, & i_group, i_class), kind=dp) &
1587  > aero_state%n_part_ideal(i_group, i_class) * 2d0)
1588  if (initial_state_warning) then
1589  call warn_msg(661936373, &
1590  "halving particles in initial condition")
1591  end if
1592  call aero_state_halve(aero_state, i_group, i_class)
1593  end do
1594  end do
1595  end do
1596  end if
1597 
1598  end subroutine aero_state_rebalance
1599 
1600 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1601 
1602  !> Scale the weighting of the given group/class by the given ratio,
1603  !> altering particle number as necessary to preserve the number
1604  !> concentration.
1605  subroutine aero_state_scale_weight(aero_state, aero_data, i_group, &
1606  i_class, weight_ratio, allow_doubling, allow_halving)
1607 
1608  !> Aerosol state.
1609  type(aero_state_t), intent(inout) :: aero_state
1610  !> Aerosol data.
1611  type(aero_data_t), intent(in) :: aero_data
1612  !> Weight group number.
1613  integer, intent(in) :: i_group
1614  !> Weight class number.
1615  integer, intent(in) :: i_class
1616  !> Ratio of <tt>new_weight / old_weight</tt>.
1617  real(kind=dp), intent(in) :: weight_ratio
1618  !> Whether to allow doubling of the population.
1619  logical, intent(in) :: allow_doubling
1620  !> Whether to allow halving of the population.
1621  logical, intent(in) :: allow_halving
1622 
1623  real(kind=dp) :: ratio
1624  integer :: i_part, i_remove, n_remove, i_entry, n_part
1625  type(aero_info_t) :: aero_info
1626 
1627  ! We could use the ratio < 1 case unconditionally, but that would
1628  ! have higher variance for the ratio > 1 case than the current
1629  ! scheme.
1630 
1631  call aero_state_sort(aero_state, aero_data)
1632  n_part = integer_varray_n_entry( &
1633  aero_state%aero_sorted%group_class%inverse(i_group, i_class))
1634 
1635  if ((weight_ratio > 1d0) .and. (allow_halving .or. (n_part == 0))) then
1636  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), &
1637  weight_ratio)
1638  n_remove = prob_round(real(n_part, kind=dp) &
1639  * (1d0 - 1d0 / weight_ratio))
1640  do i_remove = 1,n_remove
1642  aero_state%aero_sorted%group_class%inverse(i_group, i_class)))
1643  i_part = aero_state%aero_sorted%group_class%inverse(i_group, &
1644  i_class)%entry(i_entry)
1645  aero_info%id = aero_state%apa%particle(i_part)%id
1646  aero_info%action = aero_info_halved
1647  call aero_state_remove_particle(aero_state, i_part, .true., &
1648  aero_info)
1649  end do
1650  elseif ((weight_ratio < 1d0) &
1651  .and. (allow_doubling .or. (n_part == 0))) then
1652  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), &
1653  weight_ratio)
1654  do i_entry = n_part,1,-1
1655  i_part = aero_state%aero_sorted%group_class%inverse(i_group, &
1656  i_class)%entry(i_entry)
1657  call aero_state_dup_particle(aero_state, aero_data, i_part, &
1658  1d0 / weight_ratio)
1659  end do
1660  end if
1661 
1662  end subroutine aero_state_scale_weight
1663 
1664 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1665 
1666  !> Mix the aero_states between all processes. Currently uses a
1667  !> simple all-to-all diffusion.
1668  subroutine aero_state_mix(aero_state, del_t, mix_timescale, &
1669  aero_data, specify_prob_transfer)
1670 
1671  !> Aerosol state.
1672  type(aero_state_t), intent(inout) :: aero_state
1673  !> Timestep (s).
1674  real(kind=dp), intent(in) :: del_t
1675  !> Mixing timescale (s).
1676  real(kind=dp), intent(in) :: mix_timescale
1677  !> Aero data values.
1678  type(aero_data_t), intent(in) :: aero_data
1679  !> Transfer probability of each particle (0 means no mixing, 1
1680  !> means total mixing).
1681  real(kind=dp), optional, intent(in) :: specify_prob_transfer
1682 
1683 #ifdef PMC_USE_MPI
1684  integer :: rank, n_proc, i_proc, ierr
1685  integer :: buffer_size, buffer_size_check
1686  character, allocatable :: buffer(:)
1687  type(aero_state_t), allocatable :: aero_state_sends(:)
1688  type(aero_state_t), allocatable :: aero_state_recvs(:)
1689  real(kind=dp) :: prob_transfer, prob_not_transferred
1690  real(kind=dp) :: prob_transfer_given_not_transferred
1691 
1692  ! process information
1693  rank = pmc_mpi_rank()
1694  n_proc = pmc_mpi_size()
1695  if (n_proc == 1) then
1696  ! buffer allocation below fails if n_proc == 1
1697  ! so bail out early (nothing to mix anyway)
1698  return
1699  end if
1700 
1701  ! allocate aero_state arrays
1702  allocate(aero_state_sends(n_proc))
1703  allocate(aero_state_recvs(n_proc))
1704 
1705  ! compute the transfer probability
1706  if (present(specify_prob_transfer)) then
1707  prob_transfer = specify_prob_transfer / real(n_proc, kind=dp)
1708  else
1709  prob_transfer = (1d0 - exp(- del_t / mix_timescale)) &
1710  / real(n_proc, kind=dp)
1711  end if
1712 
1713  ! extract particles to send
1714  prob_not_transferred = 1d0
1715  do i_proc = 0,(n_proc - 1)
1716  if (i_proc /= rank) then
1717  ! because we are doing sequential sampling from the aero_state
1718  ! we need to scale up the later transfer probabilities, because
1719  ! the later particles are being transferred conditioned on the
1720  ! fact that they did not transfer already
1721  prob_transfer_given_not_transferred = prob_transfer &
1722  / prob_not_transferred
1723  call aero_state_sample(aero_state, &
1724  aero_state_sends(i_proc + 1), aero_data, &
1725  prob_transfer_given_not_transferred, aero_info_none)
1726  prob_not_transferred = prob_not_transferred - prob_transfer
1727  end if
1728  end do
1729 
1730  ! exchange the particles
1731  call aero_state_mpi_alltoall(aero_state_sends, aero_state_recvs)
1732 
1733  ! process the received particles
1734  do i_proc = 0,(n_proc - 1)
1735  if (i_proc /= rank) then
1736  call aero_state_add(aero_state, aero_state_recvs(i_proc + 1), &
1737  aero_data)
1738  end if
1739  end do
1740 
1741  ! cleanup
1742  deallocate(aero_state_sends)
1743  deallocate(aero_state_recvs)
1744 #endif
1745 
1746  end subroutine aero_state_mix
1747 
1748 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1749 
1750  !> Do an MPI all-to-all transfer of aerosol states.
1751  !!
1752  !! States without particles are not sent.
1753  subroutine aero_state_mpi_alltoall(send, recv)
1754 
1755  !> Array of aero_states to send (one per process).
1756  type(aero_state_t), intent(in) :: send(:)
1757  !> Array of aero_states to receives (one per process).
1758  type(aero_state_t), intent(inout) :: recv(size(send))
1759 
1760 #ifdef PMC_USE_MPI
1761  character, allocatable :: sendbuf(:), recvbuf(:)
1762  integer :: sendcounts(size(send)), sdispls(size(send))
1763  integer :: recvcounts(size(send)), rdispls(size(send))
1764  integer :: i_proc, position, old_position, max_sendbuf_size, ierr
1765 
1766  call assert(978709842, size(send) == pmc_mpi_size())
1767 
1768  max_sendbuf_size = 0
1769  do i_proc = 1,pmc_mpi_size()
1770  if (aero_state_total_particles(send(i_proc)) > 0) then
1771  max_sendbuf_size = max_sendbuf_size &
1772  + pmc_mpi_pack_size_aero_state(send(i_proc))
1773  end if
1774  end do
1775 
1776  allocate(sendbuf(max_sendbuf_size))
1777 
1778  position = 0
1779  do i_proc = 1,pmc_mpi_size()
1780  old_position = position
1781  if (aero_state_total_particles(send(i_proc)) > 0) then
1782  call pmc_mpi_pack_aero_state(sendbuf, position, send(i_proc))
1783  end if
1784  sendcounts(i_proc) = position - old_position
1785  end do
1786  call assert(393267406, position <= max_sendbuf_size)
1787 
1788  call pmc_mpi_alltoall_integer(sendcounts, recvcounts)
1789  allocate(recvbuf(sum(recvcounts)))
1790 
1791  sdispls(1) = 0
1792  rdispls(1) = 0
1793  do i_proc = 2,pmc_mpi_size()
1794  sdispls(i_proc) = sdispls(i_proc - 1) + sendcounts(i_proc - 1)
1795  rdispls(i_proc) = rdispls(i_proc - 1) + recvcounts(i_proc - 1)
1796  end do
1797 
1798  call mpi_alltoallv(sendbuf, sendcounts, sdispls, mpi_character, recvbuf, &
1799  recvcounts, rdispls, mpi_character, mpi_comm_world, ierr)
1800  call pmc_mpi_check_ierr(ierr)
1801 
1802  position = 0
1803  do i_proc = 1,pmc_mpi_size()
1804  call assert(189739257, position == rdispls(i_proc))
1805  if (recvcounts(i_proc) > 0) then
1806  call pmc_mpi_unpack_aero_state(recvbuf, position, recv(i_proc))
1807  end if
1808  end do
1809 
1810  deallocate(sendbuf)
1811  deallocate(recvbuf)
1812 #endif
1813 
1814  end subroutine aero_state_mpi_alltoall
1815 
1816 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1817 
1818  !> Set each aerosol particle to have its original total volume, but
1819  !> species volume ratios given by the total species volume ratio
1820  !> within each bin. This preserves the (weighted) total species
1821  !> volume per bin as well as per-particle total volumes.
1822  subroutine aero_state_bin_average_comp(aero_state, bin_grid, aero_data)
1823 
1824  !> Aerosol state to average.
1825  type(aero_state_t), intent(inout) :: aero_state
1826  !> Bin grid to average within.
1827  type(bin_grid_t), intent(in) :: bin_grid
1828  !> Aerosol data.
1829  type(aero_data_t), intent(in) :: aero_data
1830 
1831  real(kind=dp) :: species_volume_conc(aero_data_n_spec(aero_data))
1832  real(kind=dp) :: total_volume_conc, particle_volume, num_conc
1833  integer :: i_bin, i_class, i_entry, i_part, i_spec
1834 
1835  call aero_state_sort(aero_state, aero_data, bin_grid)
1836 
1837  do i_bin = 1,bin_grid_size(bin_grid)
1838  species_volume_conc = 0d0
1839  total_volume_conc = 0d0
1840  do i_class = 1,size(aero_state%awa%weight, 2)
1841  do i_entry = 1,integer_varray_n_entry( &
1842  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1843  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1844  i_class)%entry(i_entry)
1845  num_conc = aero_weight_array_num_conc(aero_state%awa, &
1846  aero_state%apa%particle(i_part), aero_data)
1847  particle_volume = aero_particle_volume( &
1848  aero_state%apa%particle(i_part))
1849  species_volume_conc = species_volume_conc &
1850  + num_conc * aero_state%apa%particle(i_part)%vol
1851  total_volume_conc = total_volume_conc + num_conc * particle_volume
1852  end do
1853  end do
1854  do i_class = 1,size(aero_state%awa%weight, 2)
1855  do i_entry = 1,integer_varray_n_entry( &
1856  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1857  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1858  i_class)%entry(i_entry)
1859  particle_volume = aero_particle_volume( &
1860  aero_state%apa%particle(i_part))
1861  aero_state%apa%particle(i_part)%vol &
1862  = particle_volume * species_volume_conc / total_volume_conc
1863  end do
1864  end do
1865  end do
1866 
1867  end subroutine aero_state_bin_average_comp
1868 
1869 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1870 
1871  !> Set each aerosol particle to have its original species ratios,
1872  !> but total volume given by the average volume of all particles
1873  !> within each bin.
1874  !!
1875  !! This does not preserve the total species volume
1876  !! per bin. If the \c bin_center parameter is \c .true. then the
1877  !! particles in each bin are set to have the bin center volume,
1878  !! rather than the average volume of the particles in that bin.
1879  !!
1880  !! If the weighting function is not constant (AERO_WEIGHT_TYPE_NONE)
1881  !! then the averaging can be performed in either a number-preserving
1882  !! way or in a volume-preserving way. The volume-preserving way does
1883  !! not preserve species volume ratios in gernal, but will do so if
1884  !! the particle population has already been composition-averaged.
1885  subroutine aero_state_bin_average_size(aero_state, bin_grid, aero_data, &
1886  bin_center, preserve_number)
1887 
1888  !> Aerosol state to average.
1889  type(aero_state_t), intent(inout) :: aero_state
1890  !> Bin grid to average within.
1891  type(bin_grid_t), intent(in) :: bin_grid
1892  !> Aerosol data.
1893  type(aero_data_t), intent(in) :: aero_data
1894  !> Whether to assign the bin center volume (rather than the average
1895  !> volume).
1896  logical, intent(in) :: bin_center
1897  !> Whether to use the number-preserving scheme (otherwise will use
1898  !> the volume-preserving scheme). This parameter has no effect if
1899  !> \c bin_center is \c .true.
1900  logical, intent(in) :: preserve_number
1901 
1902  real(kind=dp) :: total_volume_conc, particle_volume
1903  real(kind=dp) :: new_particle_volume, num_conc, total_num_conc
1904  real(kind=dp) :: lower_volume, upper_volume, center_volume
1905  real(kind=dp) :: lower_function, upper_function, center_function
1906  integer :: i_bin, i_class, i_entry, i_part, i_bisect, n_part
1907  logical :: monotone_increasing, monotone_decreasing
1908 
1909  call aero_state_sort(aero_state, aero_data, bin_grid)
1910 
1911  do i_bin = 1,bin_grid_size(bin_grid)
1912  do i_class = 1,size(aero_state%awa%weight, 2)
1913  if (integer_varray_n_entry( &
1914  aero_state%aero_sorted%size_class%inverse(i_bin, i_class)) &
1915  == 0) then
1916  cycle
1917  end if
1918 
1919  n_part = integer_varray_n_entry( &
1920  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1921  total_num_conc = 0d0
1922  total_volume_conc = 0d0
1923  do i_entry = 1,integer_varray_n_entry( &
1924  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1925  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1926  i_class)%entry(i_entry)
1927  num_conc = aero_weight_array_num_conc(aero_state%awa, &
1928  aero_state%apa%particle(i_part), aero_data)
1929  total_num_conc = total_num_conc + num_conc
1930  particle_volume = aero_particle_volume( &
1931  aero_state%apa%particle(i_part))
1932  total_volume_conc = total_volume_conc &
1933  + num_conc * particle_volume
1934  end do
1935 
1936  ! determine the new_particle_volume for all particles in this bin
1937  if (bin_center) then
1938  new_particle_volume = aero_data_rad2vol(aero_data, &
1939  bin_grid%centers(i_bin))
1940  elseif (aero_weight_array_check_flat(aero_state%awa)) then
1941  num_conc & ! any radius will have the same num_conc
1942  = aero_weight_array_num_conc_at_radius(aero_state%awa, &
1943  i_class, 1d0)
1944  new_particle_volume = total_volume_conc / num_conc &
1945  / real(integer_varray_n_entry( & aero_state%aero_sorted%size_class%inverse(i_bin, i_class)), &
1946  kind=dp)
1947  elseif (preserve_number) then
1948  ! number-preserving scheme: Solve the implicit equation:
1949  ! n_part * W(new_vol) = total_num_conc
1950  !
1951  ! We assume that the weighting function is strictly monotone
1952  ! so this equation has a unique solution and the solution
1953  ! lies between the min and max particle volumes. We use
1954  ! bisection as this doesn't really need to be fast, just
1955  ! robust.
1956 
1957  call aero_weight_array_check_monotonicity(aero_state%awa, &
1958  monotone_increasing, monotone_decreasing)
1959  call assert_msg(214077200, &
1960  monotone_increasing .or. monotone_decreasing, &
1961  "monotone weight function required for averaging")
1962 
1963  ! initialize to min and max particle volumes
1964  do i_entry = 1,n_part
1965  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1966  i_class)%entry(i_entry)
1967  particle_volume = aero_particle_volume( &
1968  aero_state%apa%particle(i_part))
1969  if (i_part == 1) then
1970  lower_volume = particle_volume
1971  upper_volume = particle_volume
1972  end if
1973  lower_volume = min(lower_volume, particle_volume)
1974  upper_volume = max(upper_volume, particle_volume)
1975  end do
1976 
1977  lower_function = real(n_part, kind=dp) &
1978  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1979  i_class, aero_data_vol2rad(aero_data, lower_volume)) &
1980  - total_num_conc
1981  upper_function = real(n_part, kind=dp) &
1982  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1983  i_class, aero_data_vol2rad(aero_data, upper_volume)) &
1984  - total_num_conc
1985 
1986  ! do 50 rounds of bisection (2^50 = 10^15)
1987  do i_bisect = 1,50
1988  center_volume = (lower_volume + upper_volume) / 2d0
1989  center_function = real(n_part, kind=dp) &
1990  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1991  i_class, aero_data_vol2rad(aero_data, center_volume)) &
1992  - total_num_conc
1993  if ((lower_function > 0d0 .and. center_function > 0d0) &
1994  .or. (lower_function < 0d0 .and. center_function < 0d0)) &
1995  then
1996  lower_volume = center_volume
1997  lower_function = center_function
1998  else
1999  upper_volume = center_volume
2000  upper_function = center_function
2001  end if
2002  end do
2003 
2004  new_particle_volume = center_volume
2005  else
2006  ! volume-preserving scheme: Solve the implicit equation:
2007  ! n_part * W(new_vol) * new_vol = total_volume_conc
2008  !
2009  ! We assume that the weighting function is strictly monotone
2010  ! so this equation has a unique solution and the solution
2011  ! lies between the min and max particle volumes. We use
2012  ! bisection as this doesn't really need to be fast, just
2013  ! robust.
2014 
2015  call aero_weight_array_check_monotonicity(aero_state%awa, &
2016  monotone_increasing, monotone_decreasing)
2017  call assert_msg(483078128, &
2018  monotone_increasing .or. monotone_decreasing, &
2019  "monotone weight function required for averaging")
2020 
2021  ! initialize to min and max particle volumes
2022  do i_entry = 1,n_part
2023  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
2024  i_class)%entry(i_entry)
2025  particle_volume = aero_particle_volume( &
2026  aero_state%apa%particle(i_part))
2027  if (i_part == 1) then
2028  lower_volume = particle_volume
2029  upper_volume = particle_volume
2030  end if
2031  lower_volume = min(lower_volume, particle_volume)
2032  upper_volume = max(upper_volume, particle_volume)
2033  end do
2034 
2035  lower_function = real(n_part, kind=dp) &
2036  * aero_weight_array_num_conc_at_radius( &
2037  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2038  lower_volume)) * lower_volume - total_volume_conc
2039  upper_function = real(n_part, kind=dp) &
2040  * aero_weight_array_num_conc_at_radius( &
2041  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2042  upper_volume)) * upper_volume - total_volume_conc
2043 
2044  ! do 50 rounds of bisection (2^50 = 10^15)
2045  do i_bisect = 1,50
2046  center_volume = (lower_volume + upper_volume) / 2d0
2047  center_function = real(n_part, kind=dp) &
2048  * aero_weight_array_num_conc_at_radius( &
2049  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2050  center_volume)) * center_volume - total_volume_conc
2051  if ((lower_function > 0d0 .and. center_function > 0d0) &
2052  .or. (lower_function < 0d0 .and. center_function < 0d0)) &
2053  then
2054  lower_volume = center_volume
2055  lower_function = center_function
2056  else
2057  upper_volume = center_volume
2058  upper_function = center_function
2059  end if
2060  end do
2061 
2062  new_particle_volume = center_volume
2063  end if
2064 
2065  do i_entry = 1,n_part
2066  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
2067  i_class)%entry(i_entry)
2068  particle_volume = aero_particle_volume( &
2069  aero_state%apa%particle(i_part))
2070  aero_state%apa%particle(i_part)%vol &
2071  = aero_state%apa%particle(i_part)%vol &
2072  / particle_volume * new_particle_volume
2073  end do
2074  end do
2075  end do
2076 
2077  end subroutine aero_state_bin_average_size
2078 
2079 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2080 
2081  !> Make all particles dry (water set to zero).
2082  subroutine aero_state_make_dry(aero_state, aero_data)
2083 
2084  !> Aerosol state to make dry.
2085  type(aero_state_t), intent(inout) :: aero_state
2086  !> Aerosol data.
2087  type(aero_data_t), intent(in) :: aero_data
2088 
2089  integer :: i_part
2090  real(kind=dp) :: reweight_num_conc(aero_state%apa%n_part)
2091 
2092  ! We're modifying particle diameters, so bin sorting is now invalid
2093  aero_state%valid_sort = .false.
2094 
2095  call aero_state_num_conc_for_reweight(aero_state, aero_data, &
2096  reweight_num_conc)
2097  if (aero_data%i_water > 0) then
2098  do i_part = 1,aero_state_n_part(aero_state)
2099  aero_state%apa%particle(i_part)%vol(aero_data%i_water) = 0d0
2100  end do
2101  aero_state%valid_sort = .false.
2102  end if
2103  ! adjust particles to account for weight changes
2104  call aero_state_reweight(aero_state, aero_data, reweight_num_conc)
2105 
2106  end subroutine aero_state_make_dry
2107 
2108 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2109 
2110  !> Determines the number of bytes required to pack the given value.
2111  integer function pmc_mpi_pack_size_aero_state(val)
2112 
2113  !> Value to pack.
2114  type(aero_state_t), intent(in) :: val
2115 
2116  integer :: total_size, i_group
2117 
2118  total_size = 0
2119  total_size = total_size + pmc_mpi_pack_size_apa(val%apa)
2120  total_size = total_size + pmc_mpi_pack_size_aero_weight_array(val%awa)
2121  total_size = total_size + pmc_mpi_pack_size_real_array_2d(val%n_part_ideal)
2122  total_size = total_size + pmc_mpi_pack_size_aia(val%aero_info_array)
2123  pmc_mpi_pack_size_aero_state = total_size
2124 
2125  end function pmc_mpi_pack_size_aero_state
2126 
2127 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2128 
2129  !> Packs the given value into the buffer, advancing position.
2130  subroutine pmc_mpi_pack_aero_state(buffer, position, val)
2131 
2132  !> Memory buffer.
2133  character, intent(inout) :: buffer(:)
2134  !> Current buffer position.
2135  integer, intent(inout) :: position
2136  !> Value to pack.
2137  type(aero_state_t), intent(in) :: val
2138 
2139 #ifdef PMC_USE_MPI
2140  integer :: prev_position, i_group
2141 
2142  prev_position = position
2143  call pmc_mpi_pack_aero_particle_array(buffer, position, val%apa)
2144  call pmc_mpi_pack_aero_weight_array(buffer,position,val%awa)
2145  call pmc_mpi_pack_real_array_2d(buffer, position, val%n_part_ideal)
2146  call pmc_mpi_pack_aero_info_array(buffer, position, val%aero_info_array)
2147  call assert(850997402, &
2148  position - prev_position <= pmc_mpi_pack_size_aero_state(val))
2149 #endif
2150 
2151  end subroutine pmc_mpi_pack_aero_state
2152 
2153 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2154 
2155  !> Unpacks the given value from the buffer, advancing position.
2156  subroutine pmc_mpi_unpack_aero_state(buffer, position, val)
2157 
2158  !> Memory buffer.
2159  character, intent(inout) :: buffer(:)
2160  !> Current buffer position.
2161  integer, intent(inout) :: position
2162  !> Value to pack.
2163  type(aero_state_t), intent(inout) :: val
2164 
2165 #ifdef PMC_USE_MPI
2166  integer :: prev_position, i_group, n_group
2167 
2168  val%valid_sort = .false.
2169  prev_position = position
2170  call pmc_mpi_unpack_aero_particle_array(buffer, position, val%apa)
2171  call pmc_mpi_unpack_aero_weight_array(buffer,position,val%awa)
2172  call pmc_mpi_unpack_real_array_2d(buffer, position, val%n_part_ideal)
2173  call pmc_mpi_unpack_aero_info_array(buffer, position, val%aero_info_array)
2174  call assert(132104747, &
2175  position - prev_position <= pmc_mpi_pack_size_aero_state(val))
2176 #endif
2177 
2178  end subroutine pmc_mpi_unpack_aero_state
2179 
2180 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2181 
2182  !> Gathers data from all processes into one aero_state on process 0.
2183  subroutine aero_state_mpi_gather(aero_state, aero_state_total, aero_data)
2184 
2185  !> Local aero_state.
2186  type(aero_state_t), intent(in) :: aero_state
2187  !> Centralized aero_state (only on process 0).
2188  type(aero_state_t), intent(inout) :: aero_state_total
2189  !> Aero data values.
2190  type(aero_data_t), intent(in) :: aero_data
2191 
2192 #ifdef PMC_USE_MPI
2193  type(aero_state_t) :: aero_state_transfer
2194  integer :: n_proc, ierr, status(mpi_status_size)
2195  integer :: buffer_size, max_buffer_size, i_proc, position
2196  character, allocatable :: buffer(:)
2197 #endif
2198 
2199  if (pmc_mpi_rank() == 0) then
2200  aero_state_total = aero_state
2201  end if
2202 
2203 #ifdef PMC_USE_MPI
2204 
2205  if (pmc_mpi_rank() /= 0) then
2206  ! send data from remote processes
2207  max_buffer_size = 0
2208  max_buffer_size = max_buffer_size &
2209  + pmc_mpi_pack_size_aero_state(aero_state)
2210  allocate(buffer(max_buffer_size))
2211  position = 0
2212  call pmc_mpi_pack_aero_state(buffer, position, aero_state)
2213  call assert(542772170, position <= max_buffer_size)
2214  buffer_size = position
2215  call mpi_send(buffer, buffer_size, mpi_character, 0, &
2216  aero_state_tag_gather, mpi_comm_world, ierr)
2217  call pmc_mpi_check_ierr(ierr)
2218  deallocate(buffer)
2219  else
2220  ! root process receives data from each remote proc
2221  n_proc = pmc_mpi_size()
2222  do i_proc = 1,(n_proc - 1)
2223  ! determine buffer size at root process
2224  call mpi_probe(i_proc, aero_state_tag_gather, mpi_comm_world, &
2225  status, ierr)
2226  call pmc_mpi_check_ierr(ierr)
2227  call mpi_get_count(status, mpi_character, buffer_size, ierr)
2228  call pmc_mpi_check_ierr(ierr)
2229 
2230  ! get buffer at root process
2231  allocate(buffer(buffer_size))
2232  call mpi_recv(buffer, buffer_size, mpi_character, i_proc, &
2233  aero_state_tag_gather, mpi_comm_world, status, ierr)
2234 
2235  ! unpack it
2236  position = 0
2237  call pmc_mpi_unpack_aero_state(buffer, position, &
2238  aero_state_transfer)
2239  call assert(518174881, position == buffer_size)
2240  deallocate(buffer)
2241 
2242  call aero_state_add(aero_state_total, aero_state_transfer, &
2243  aero_data)
2244  end do
2245  end if
2246 
2247 #endif
2248 
2249  end subroutine aero_state_mpi_gather
2250 
2251 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2252 
2253  !> Write the aero particle dimension to the given NetCDF file if it
2254  !> is not already present and in any case return the associated
2255  !> dimid.
2256  subroutine aero_state_netcdf_dim_aero_particle(aero_state, ncid, &
2257  dimid_aero_particle)
2258 
2259  !> aero_state structure.
2260  type(aero_state_t), intent(in) :: aero_state
2261  !> NetCDF file ID, in data mode.
2262  integer, intent(in) :: ncid
2263  !> Dimid of the aero particle dimension.
2264  integer, intent(out) :: dimid_aero_particle
2265 
2266  integer :: status, i_part
2267  integer :: varid_aero_particle
2268  integer :: aero_particle_centers(aero_state_n_part(aero_state))
2269 
2270  ! try to get the dimension ID
2271  status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
2272  if (status == nf90_noerr) return
2273  if (status /= nf90_ebaddim) call pmc_nc_check(status)
2274 
2275  ! dimension not defined, so define now define it
2276  call pmc_nc_check(nf90_redef(ncid))
2277 
2278  call pmc_nc_check(nf90_def_dim(ncid, "aero_particle", &
2279  aero_state_n_part(aero_state), dimid_aero_particle))
2280 
2281  call pmc_nc_check(nf90_enddef(ncid))
2282 
2283  do i_part = 1,aero_state_n_part(aero_state)
2284  aero_particle_centers(i_part) = i_part
2285  end do
2286  call pmc_nc_write_integer_1d(ncid, aero_particle_centers, &
2287  "aero_particle", (/ dimid_aero_particle /), &
2288  description="dummy dimension variable (no useful value)")
2289 
2291 
2292 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2293 
2294  !> Write the aero removed dimension to the given NetCDF file if it
2295  !> is not already present and in any case return the associated
2296  !> dimid.
2297  subroutine aero_state_netcdf_dim_aero_removed(aero_state, ncid, &
2298  dimid_aero_removed)
2299 
2300  !> aero_state structure.
2301  type(aero_state_t), intent(in) :: aero_state
2302  !> NetCDF file ID, in data mode.
2303  integer, intent(in) :: ncid
2304  !> Dimid of the aero removed dimension.
2305  integer, intent(out) :: dimid_aero_removed
2306 
2307  integer :: status, i_remove, dim_size
2308  integer :: varid_aero_removed
2309  integer :: aero_removed_centers( &
2310  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2311 
2312  ! try to get the dimension ID
2313  status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed)
2314  if (status == nf90_noerr) return
2315  if (status /= nf90_ebaddim) call pmc_nc_check(status)
2316 
2317  ! dimension not defined, so define now define it
2318  call pmc_nc_check(nf90_redef(ncid))
2319 
2320  dim_size = max(1, aero_info_array_n_item(aero_state%aero_info_array))
2321  call pmc_nc_check(nf90_def_dim(ncid, "aero_removed", &
2322  dim_size, dimid_aero_removed))
2323 
2324  call pmc_nc_check(nf90_enddef(ncid))
2325 
2326  do i_remove = 1,dim_size
2327  aero_removed_centers(i_remove) = i_remove
2328  end do
2329  call pmc_nc_write_integer_1d(ncid, aero_removed_centers, &
2330  "aero_removed", (/ dimid_aero_removed /), &
2331  description="dummy dimension variable (no useful value)")
2332 
2333  end subroutine aero_state_netcdf_dim_aero_removed
2334 
2335 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2336 
2337  !> Write full state.
2338  subroutine aero_state_output_netcdf(aero_state, ncid, aero_data, &
2339  record_removals, record_optical)
2340 
2341  !> aero_state to write.
2342  type(aero_state_t), intent(in) :: aero_state
2343  !> NetCDF file ID, in data mode.
2344  integer, intent(in) :: ncid
2345  !> aero_data structure.
2346  type(aero_data_t), intent(in) :: aero_data
2347  !> Whether to output particle removal info.
2348  logical, intent(in) :: record_removals
2349  !> Whether to output aerosol optical properties.
2350  logical, intent(in) :: record_optical
2351 
2352  integer :: dimid_aero_particle, dimid_aero_species, dimid_aero_source
2353  integer :: dimid_aero_removed
2354  integer :: i_part, i_remove
2355  real(kind=dp) :: aero_particle_mass(aero_state_n_part(aero_state), &
2356  aero_data_n_spec(aero_data))
2357  integer :: aero_n_orig_part(aero_state_n_part(aero_state), &
2358  aero_data_n_source(aero_data))
2359  integer :: aero_particle_weight_group(aero_state_n_part(aero_state))
2360  integer :: aero_particle_weight_class(aero_state_n_part(aero_state))
2361  real(kind=dp) :: aero_absorb_cross_sect(aero_state_n_part(aero_state))
2362  real(kind=dp) :: aero_scatter_cross_sect(aero_state_n_part(aero_state))
2363  real(kind=dp) :: aero_asymmetry(aero_state_n_part(aero_state))
2364  real(kind=dp) :: aero_refract_shell_real(aero_state_n_part(aero_state))
2365  real(kind=dp) :: aero_refract_shell_imag(aero_state_n_part(aero_state))
2366  real(kind=dp) :: aero_refract_core_real(aero_state_n_part(aero_state))
2367  real(kind=dp) :: aero_refract_core_imag(aero_state_n_part(aero_state))
2368  real(kind=dp) :: aero_core_vol(aero_state_n_part(aero_state))
2369  integer :: aero_water_hyst_leg(aero_state_n_part(aero_state))
2370  real(kind=dp) :: aero_num_conc(aero_state_n_part(aero_state))
2371  integer :: aero_id(aero_state_n_part(aero_state))
2372  real(kind=dp) :: aero_least_create_time(aero_state_n_part(aero_state))
2373  real(kind=dp) :: aero_greatest_create_time(aero_state_n_part(aero_state))
2374  integer :: aero_removed_id( &
2375  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2376  integer :: aero_removed_action( &
2377  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2378  integer :: aero_removed_other_id( &
2379  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2380 
2381  !> \page output_format_aero_state Output File Format: Aerosol Particle State
2382  !!
2383  !! The aerosol state consists of a set of individual aerosol
2384  !! particles, each with its own individual properties. The
2385  !! properties of all particles are stored in arrays, one per
2386  !! property. For example, <tt>aero_absorb_cross_sect(i)</tt> gives
2387  !! the absorption cross section of particle number \c i, while
2388  !! <tt>aero_particle_mass(i,s)</tt> gives the mass of species \c s
2389  !! in particle \c i. The aerosol species are described in \ref
2390  !! output_format_aero_data.
2391  !!
2392  !! Each aerosol particle \c i represents a number concentration
2393  !! given by <tt>aero_num_conc(i)</tt>. Multiplying a per-particle
2394  !! quantity by the respective number concentration gives the
2395  !! concentration of that quantity contributed by the particle. For
2396  !! example, summing <tt>aero_particle_mass(i,s) *
2397  !! aero_num_conc(i)</tt> over all \c i gives the total mass
2398  !! concentration of species \c s in kg/m^3. Similarly, summing
2399  !! <tt>aero_absorb_cross_sect(i) * aero_num_conc(i)</tt> over all
2400  !! \c i will give the concentration of scattering cross section in
2401  !! m^2/m^3.
2402  !!
2403  !! FIXME: the aero_weight is also output
2404  !!
2405  !! The aerosol state uses the \c aero_species NetCDF dimension as
2406  !! specified in the \ref output_format_aero_data section, as well
2407  !! as the NetCDF dimension:
2408  !! - \b aero_particle: number of aerosol particles
2409  !!
2410  !! The aerosol state NetCDF variables are:
2411  !! - \b aero_particle (dim \c aero_particle): dummy dimension variable
2412  !! (no useful value)
2413  !! - \b aero_particle_mass (unit kg,
2414  !! dim <tt>aero_particle x aero_species</tt>): constituent masses of
2415  !! each aerosol particle - <tt>aero_particle_mass(i,s)</tt> gives the
2416  !! mass of species \c s in particle \c i
2417  !! - \b aero_n_orig_part (dim <tt>aero_particle x
2418  !! aero_source</tt>): number of original particles from each
2419  !! source that formed each aerosol particle -
2420  !! <tt>aero_n_orig_part(i,s)</tt> is the number of particles
2421  !! from source \c s that contributed to particle \c i - when
2422  !! particle \c i first enters the simulation (by emissions,
2423  !! dilution, etc.) it has <tt>aero_n_orig_part(i,s) = 1</tt>
2424  !! for the source number \c s it came from (otherwise zero)
2425  !! and when two particles coagulate, their values of \c
2426  !! aero_n_orig_part are added (the number of coagulation
2427  !! events that formed each particle is thus
2428  !! <tt>sum(aero_n_orig_part(i,:)) - 1</tt>)
2429  !! - \b aero_particle_weight_group (dim <tt>aero_particle</tt>):
2430  !! weight group number (1 to <tt>aero_weight_group</tt>) of
2431  !! each aerosol particle
2432  !! - \b aero_particle_weight_class (dim <tt>aero_particle</tt>):
2433  !! weight class number (1 to <tt>aero_weight_class</tt>) of each
2434  !! aerosol particle
2435  !! - \b aero_absorb_cross_sect (unit m^2, dim \c aero_particle):
2436  !! optical absorption cross sections of each aerosol particle
2437  !! - \b aero_scatter_cross_sect (unit m^2, dim \c aero_particle):
2438  !! optical scattering cross sections of each aerosol particle
2439  !! - \b aero_asymmetry (dimensionless, dim \c aero_particle): optical
2440  !! asymmetry parameters of each aerosol particle
2441  !! - \b aero_refract_shell_real (dimensionless, dim \c aero_particle):
2442  !! real part of the refractive indices of the shell of each
2443  !! aerosol particle
2444  !! - \b aero_refract_shell_imag (dimensionless, dim \c aero_particle):
2445  !! imaginary part of the refractive indices of the shell of each
2446  !! aerosol particle
2447  !! - \b aero_refract_core_real (dimensionless, dim \c aero_particle):
2448  !! real part of the refractive indices of the core of each
2449  !! aerosol particle
2450  !! - \b aero_refract_core_imag (dimensionless, dim \c aero_particle):
2451  !! imaginary part of the refractive indices of the core of each
2452  !! aerosol particle
2453  !! - \b aero_core_vol (unit m^3, dim \c aero_particle): volume of the
2454  !! optical cores of each aerosol particle
2455  !! - \b aero_water_hyst_leg (dim \c aero_particle): integers
2456  !! specifying which leg of the water hysteresis curve each
2457  !! particle is on, using the MOSAIC numbering convention
2458  !! - \b aero_num_conc (unit m^{-3}, dim \c aero_particle): number
2459  !! concentration associated with each particle
2460  !! - \b aero_id (dim \c aero_particle): unique ID number of each
2461  !! aerosol particle
2462  !! - \b aero_least_create_time (unit s, dim \c aero_particle): least
2463  !! (earliest) creation time of any original constituent particles
2464  !! that coagulated to form each particle, measured from the start
2465  !! of the simulation - a particle is said to be created when it
2466  !! first enters the simulation (by emissions, dilution, etc.)
2467  !! - \b aero_greatest_create_time (unit s, dim \c
2468  !! aero_particle): greatest (latest) creation time of any
2469  !! original constituent particles that coagulated to form each
2470  !! particle, measured from the start of the simulation - a
2471  !! particle is said to be created when it first enters the
2472  !! simulation (by emissions, dilution, etc.)
2473 
2474  call aero_weight_array_output_netcdf(aero_state%awa, ncid)
2475 
2476  call aero_data_netcdf_dim_aero_species(aero_data, ncid, &
2477  dimid_aero_species)
2478  call aero_data_netcdf_dim_aero_source(aero_data, ncid, &
2479  dimid_aero_source)
2480 
2481  if (aero_state_n_part(aero_state) > 0) then
2482  call aero_state_netcdf_dim_aero_particle(aero_state, ncid, &
2483  dimid_aero_particle)
2484 
2485  do i_part = 1,aero_state_n_part(aero_state)
2486  aero_particle_mass(i_part, :) &
2487  = aero_state%apa%particle(i_part)%vol * aero_data%density
2488  aero_n_orig_part(i_part, :) &
2489  = aero_state%apa%particle(i_part)%n_orig_part
2490  aero_particle_weight_group(i_part) &
2491  = aero_state%apa%particle(i_part)%weight_group
2492  aero_particle_weight_class(i_part) &
2493  = aero_state%apa%particle(i_part)%weight_class
2494  aero_water_hyst_leg(i_part) &
2495  = aero_state%apa%particle(i_part)%water_hyst_leg
2496  aero_num_conc(i_part) &
2497  = aero_state_particle_num_conc(aero_state, &
2498  aero_state%apa%particle(i_part), aero_data)
2499  aero_id(i_part) = aero_state%apa%particle(i_part)%id
2500  aero_least_create_time(i_part) &
2501  = aero_state%apa%particle(i_part)%least_create_time
2502  aero_greatest_create_time(i_part) &
2503  = aero_state%apa%particle(i_part)%greatest_create_time
2504  if (record_optical) then
2505  aero_absorb_cross_sect(i_part) &
2506  = aero_state%apa%particle(i_part)%absorb_cross_sect
2507  aero_scatter_cross_sect(i_part) &
2508  = aero_state%apa%particle(i_part)%scatter_cross_sect
2509  aero_asymmetry(i_part) = aero_state%apa%particle(i_part)%asymmetry
2510  aero_refract_shell_real(i_part) &
2511  = real(aero_state%apa%particle(i_part)%refract_shell)
2512  aero_refract_shell_imag(i_part) &
2513  = aimag(aero_state%apa%particle(i_part)%refract_shell)
2514  aero_refract_core_real(i_part) &
2515  = real(aero_state%apa%particle(i_part)%refract_core)
2516  aero_refract_core_imag(i_part) &
2517  = aimag(aero_state%apa%particle(i_part)%refract_core)
2518  aero_core_vol(i_part) = aero_state%apa%particle(i_part)%core_vol
2519  end if
2520  end do
2521  call pmc_nc_write_real_2d(ncid, aero_particle_mass, &
2522  "aero_particle_mass", (/ dimid_aero_particle, &
2523  dimid_aero_species /), unit="kg", &
2524  long_name="constituent masses of each aerosol particle")
2525  call pmc_nc_write_integer_2d(ncid, aero_n_orig_part, &
2526  "aero_n_orig_part", (/ dimid_aero_particle, &
2527  dimid_aero_source /), &
2528  long_name="number of original constituent particles from " &
2529  // "each source that coagulated to form each aerosol particle")
2530  call pmc_nc_write_integer_1d(ncid, aero_particle_weight_group, &
2531  "aero_particle_weight_group", (/ dimid_aero_particle /), &
2532  long_name="weight group number of each aerosol particle")
2533  call pmc_nc_write_integer_1d(ncid, aero_particle_weight_class, &
2534  "aero_particle_weight_class", (/ dimid_aero_particle /), &
2535  long_name="weight class number of each aerosol particle")
2536  call pmc_nc_write_integer_1d(ncid, aero_water_hyst_leg, &
2537  "aero_water_hyst_leg", (/ dimid_aero_particle /), &
2538  long_name="leg of the water hysteresis curve leg of each "&
2539  // "aerosol particle")
2540  call pmc_nc_write_real_1d(ncid, aero_num_conc, &
2541  "aero_num_conc", (/ dimid_aero_particle /), unit="m^{-3}", &
2542  long_name="number concentration for each particle")
2543  call pmc_nc_write_integer_1d(ncid, aero_id, &
2544  "aero_id", (/ dimid_aero_particle /), &
2545  long_name="unique ID number of each aerosol particle")
2546  call pmc_nc_write_real_1d(ncid, aero_least_create_time, &
2547  "aero_least_create_time", (/ dimid_aero_particle /), unit="s", &
2548  long_name="least creation time of each aerosol particle", &
2549  description="least (earliest) creation time of any original " &
2550  // "constituent particles that coagulated to form each " &
2551  // "particle, measured from the start of the simulation")
2552  call pmc_nc_write_real_1d(ncid, aero_greatest_create_time, &
2553  "aero_greatest_create_time", (/ dimid_aero_particle /), &
2554  unit="s", &
2555  long_name="greatest creation time of each aerosol particle", &
2556  description="greatest (latest) creation time of any original " &
2557  // "constituent particles that coagulated to form each " &
2558  // "particle, measured from the start of the simulation")
2559  if (record_optical) then
2560  call pmc_nc_write_real_1d(ncid, aero_absorb_cross_sect, &
2561  "aero_absorb_cross_sect", (/ dimid_aero_particle /), &
2562  unit="m^2", &
2563  long_name="optical absorption cross sections of each " &
2564  // "aerosol particle")
2565  call pmc_nc_write_real_1d(ncid, aero_scatter_cross_sect, &
2566  "aero_scatter_cross_sect", (/ dimid_aero_particle /), &
2567  unit="m^2", &
2568  long_name="optical scattering cross sections of each " &
2569  // "aerosol particle")
2570  call pmc_nc_write_real_1d(ncid, aero_asymmetry, &
2571  "aero_asymmetry", (/ dimid_aero_particle /), unit="1", &
2572  long_name="optical asymmetry parameters of each " &
2573  // "aerosol particle")
2574  call pmc_nc_write_real_1d(ncid, aero_refract_shell_real, &
2575  "aero_refract_shell_real", (/ dimid_aero_particle /), &
2576  unit="1", &
2577  long_name="real part of the refractive indices of the " &
2578  // "shell of each aerosol particle")
2579  call pmc_nc_write_real_1d(ncid, aero_refract_shell_imag, &
2580  "aero_refract_shell_imag", (/ dimid_aero_particle /), &
2581  unit="1", &
2582  long_name="imaginary part of the refractive indices of " &
2583  // "the shell of each aerosol particle")
2584  call pmc_nc_write_real_1d(ncid, aero_refract_core_real, &
2585  "aero_refract_core_real", (/ dimid_aero_particle /), &
2586  unit="1", &
2587  long_name="real part of the refractive indices of the core " &
2588  // "of each aerosol particle")
2589  call pmc_nc_write_real_1d(ncid, aero_refract_core_imag, &
2590  "aero_refract_core_imag", (/ dimid_aero_particle /), &
2591  unit="1", &
2592  long_name="imaginary part of the refractive indices of " &
2593  // "the core of each aerosol particle")
2594  call pmc_nc_write_real_1d(ncid, aero_core_vol, &
2595  "aero_core_vol", (/ dimid_aero_particle /), unit="m^3", &
2596  long_name="volume of the optical cores of each " &
2597  // "aerosol particle")
2598  end if
2599  end if
2600 
2601  ! FIXME: move this to aero_info_array.F90, together with
2602  ! aero_state_netcdf_dim_aero_removed() ?
2603  if (record_removals) then
2604  call aero_state_netcdf_dim_aero_removed(aero_state, ncid, &
2605  dimid_aero_removed)
2606  if (aero_info_array_n_item(aero_state%aero_info_array) >= 1) then
2607  do i_remove = 1,aero_info_array_n_item(aero_state%aero_info_array)
2608  aero_removed_id(i_remove) = &
2609  aero_state%aero_info_array%aero_info(i_remove)%id
2610  aero_removed_action(i_remove) = &
2611  aero_state%aero_info_array%aero_info(i_remove)%action
2612  aero_removed_other_id(i_remove) = &
2613  aero_state%aero_info_array%aero_info(i_remove)%other_id
2614  end do
2615  else
2616  aero_removed_id(1) = 0
2617  aero_removed_action(1) = aero_info_none
2618  aero_removed_other_id(1) = 0
2619  end if
2620  call pmc_nc_write_integer_1d(ncid, aero_removed_id, &
2621  "aero_removed_id", (/ dimid_aero_removed /), &
2622  long_name="ID of removed particles")
2623  call pmc_nc_write_integer_1d(ncid, aero_removed_action, &
2624  "aero_removed_action", (/ dimid_aero_removed /), &
2625  long_name="reason for particle removal", &
2626  description="valid is 0 (invalid entry), 1 (removed due to " &
2627  // "dilution), 2 (removed due to coagulation -- combined " &
2628  // "particle ID is in \c aero_removed_other_id), 3 (removed " &
2629  // "due to populating halving), or 4 (removed due to " &
2630  // "weighting changes")
2631  call pmc_nc_write_integer_1d(ncid, aero_removed_other_id, &
2632  "aero_removed_other_id", (/ dimid_aero_removed /), &
2633  long_name="ID of other particle involved in removal", &
2634  description="if <tt>aero_removed_action(i)</tt> is 2 " &
2635  // "(due to coagulation), then " &
2636  // "<tt>aero_removed_other_id(i)</tt> is the ID of the " &
2637  // "resulting combined particle, or 0 if the new particle " &
2638  // "was not created")
2639  end if
2640 
2641  end subroutine aero_state_output_netcdf
2642 
2643  ! this belongs in the subroutine above, but is outside because
2644  ! Doxygen 1.8.7 doesn't resolve references when multiple \page
2645  ! blocks are in one subroutine
2646 
2647  !> \page output_format_aero_removed Output File Format: Aerosol Particle Removal Information
2648  !!
2649  !! When an aerosol particle is introduced into the simulation it
2650  !! is assigned a unique ID number. This ID number will persist
2651  !! over time, allowing tracking of a paticular particle's
2652  !! evolution. If the \c record_removals variable in the input spec
2653  !! file is \c yes, then the every time a particle is removed from
2654  !! the simulation its removal will be recorded in the removal
2655  !! information.
2656  !!
2657  !! The removal information written at timestep \c n contains
2658  !! information about every particle ID that is present at time
2659  !! <tt>(n - 1)</tt> but not present at time \c n.
2660  !!
2661  !! The removal information is always written in the output files,
2662  !! even if no particles were removed in the previous
2663  !! timestep. Unfortunately, NetCDF files cannot contain arrays of
2664  !! length 0. In the case of no particles being removed, the \c
2665  !! aero_removed dimension will be set to 1 and
2666  !! <tt>aero_removed_action(1)</tt> will be 0 (\c AERO_INFO_NONE).
2667  !!
2668  !! When two particles coagulate, the ID number of the combined
2669  !! particle will be the ID particle of the largest constituent, if
2670  !! possible (weighting functions can make this impossible to
2671  !! achieve). A given particle ID may thus be lost due to
2672  !! coagulation (if the resulting combined particle has a different
2673  !! ID), or the ID may be preserved (as the ID of the combined
2674  !! particle). Only if the ID is lost will the particle be recorded
2675  !! in the removal information, and in this case
2676  !! <tt>aero_removed_action(i)</tt> will be 2 (\c AERO_INFO_COAG)
2677  !! and <tt>aero_removed_other_id(i)</tt> will be the ID number of
2678  !! the combined particle.
2679  !!
2680  !! The aerosol removal information NetCDF dimensions are:
2681  !! - \b aero_removed: number of aerosol particles removed from the
2682  !! simulation during the previous timestep (or 1, as described
2683  !! above)
2684  !!
2685  !! The aerosol removal information NetCDF variables are:
2686  !! - \b aero_removed (dim \c aero_removed): dummy dimension variable
2687  !! (no useful value)
2688  !! - \b aero_removed_id (dim \c aero_removed): the ID number of each
2689  !! removed particle
2690  !! - \b aero_removed_action (dim \c aero_removed): the reasons for
2691  !! removal for each particle, with values:
2692  !! - 0 (\c AERO_INFO_NONE): no information (invalid entry)
2693  !! - 1 (\c AERO_INFO_DILUTION): particle was removed due to dilution
2694  !! with outside air
2695  !! - 2 (\c AERO_INFO_COAG): particle was removed due to coagulation
2696  !! - 3 (\c AERO_INFO_HALVED): particle was removed due to halving of
2697  !! the aerosol population
2698  !! - 4 (\c AERO_INFO_WEIGHT): particle was removed due to adjustments
2699  !! in the particle's weighting function
2700  !! - \b aero_removed_other_id (dim \c aero_removed): the ID number of
2701  !! the combined particle formed by coagulation, if the removal reason
2702  !! was coagulation (2, \c AERO_INFO_COAG). May be 0, if the new
2703  !! coagulated particle was not created due to weighting.
2704 
2705 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2706 
2707  !> Read full state.
2708  subroutine aero_state_input_netcdf(aero_state, ncid, aero_data)
2709 
2710  !> aero_state to read.
2711  type(aero_state_t), intent(inout) :: aero_state
2712  !> NetCDF file ID, in data mode.
2713  integer, intent(in) :: ncid
2714  !> aero_data structure.
2715  type(aero_data_t), intent(in) :: aero_data
2716 
2717  integer :: dimid_aero_particle, dimid_aero_removed, n_info_item, n_part
2718  integer :: i_bin, i_part_in_bin, i_part, i_remove, status
2719  type(aero_particle_t) :: aero_particle
2720  character(len=1000) :: name
2721 
2722  real(kind=dp), allocatable :: aero_particle_mass(:,:)
2723  integer, allocatable :: aero_n_orig_part(:,:)
2724  integer, allocatable :: aero_particle_weight_group(:)
2725  integer, allocatable :: aero_particle_weight_class(:)
2726  real(kind=dp), allocatable :: aero_absorb_cross_sect(:)
2727  real(kind=dp), allocatable :: aero_scatter_cross_sect(:)
2728  real(kind=dp), allocatable :: aero_asymmetry(:)
2729  real(kind=dp), allocatable :: aero_refract_shell_real(:)
2730  real(kind=dp), allocatable :: aero_refract_shell_imag(:)
2731  real(kind=dp), allocatable :: aero_refract_core_real(:)
2732  real(kind=dp), allocatable :: aero_refract_core_imag(:)
2733  real(kind=dp), allocatable :: aero_core_vol(:)
2734  integer, allocatable :: aero_water_hyst_leg(:)
2735  real(kind=dp), allocatable :: aero_num_conc(:)
2736  integer, allocatable :: aero_id(:)
2737  real(kind=dp), allocatable :: aero_least_create_time(:)
2738  real(kind=dp), allocatable :: aero_greatest_create_time(:)
2739  integer, allocatable :: aero_removed_id(:)
2740  integer, allocatable :: aero_removed_action(:)
2741  integer, allocatable :: aero_removed_other_id(:)
2742 
2743  call aero_state_zero(aero_state)
2744 
2745  status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
2746  if (status == nf90_ebaddim) then
2747  ! no aero_particle dimension means no particles present
2748  call aero_weight_array_input_netcdf(aero_state%awa, ncid)
2749  return
2750  end if
2751  call pmc_nc_check(status)
2752  call pmc_nc_check(nf90_inquire_dimension(ncid, dimid_aero_particle, &
2753  name, n_part))
2754 
2755  call pmc_nc_read_real_2d(ncid, aero_particle_mass, &
2756  "aero_particle_mass")
2757  call pmc_nc_read_integer_2d(ncid, aero_n_orig_part, &
2758  "aero_n_orig_part")
2759  call pmc_nc_read_integer_1d(ncid, aero_particle_weight_group, &
2760  "aero_particle_weight_group")
2761  call pmc_nc_read_integer_1d(ncid, aero_particle_weight_class, &
2762  "aero_particle_weight_class")
2763  call pmc_nc_read_real_1d(ncid, aero_absorb_cross_sect, &
2764  "aero_absorb_cross_sect", must_be_present=.false.)
2765  call pmc_nc_read_real_1d(ncid, aero_scatter_cross_sect, &
2766  "aero_scatter_cross_sect", must_be_present=.false.)
2767  call pmc_nc_read_real_1d(ncid, aero_asymmetry, &
2768  "aero_asymmetry", must_be_present=.false.)
2769  call pmc_nc_read_real_1d(ncid, aero_refract_shell_real, &
2770  "aero_refract_shell_real", must_be_present=.false.)
2771  call pmc_nc_read_real_1d(ncid, aero_refract_shell_imag, &
2772  "aero_refract_shell_imag", must_be_present=.false.)
2773  call pmc_nc_read_real_1d(ncid, aero_refract_core_real, &
2774  "aero_refract_core_real", must_be_present=.false.)
2775  call pmc_nc_read_real_1d(ncid, aero_refract_core_imag, &
2776  "aero_refract_core_imag", must_be_present=.false.)
2777  call pmc_nc_read_real_1d(ncid, aero_core_vol, &
2778  "aero_core_vol", must_be_present=.false.)
2779  call pmc_nc_read_integer_1d(ncid, aero_water_hyst_leg, &
2780  "aero_water_hyst_leg")
2781  call pmc_nc_read_real_1d(ncid, aero_num_conc, &
2782  "aero_num_conc")
2783  call pmc_nc_read_integer_1d(ncid, aero_id, &
2784  "aero_id")
2785  call pmc_nc_read_real_1d(ncid, aero_least_create_time, &
2786  "aero_least_create_time")
2787  call pmc_nc_read_real_1d(ncid, aero_greatest_create_time, &
2788  "aero_greatest_create_time")
2789 
2790  call aero_weight_array_input_netcdf(aero_state%awa, ncid)
2791  call aero_state_set_n_part_ideal(aero_state, 0d0)
2792 
2793  do i_part = 1,n_part
2794  aero_particle%vol = aero_particle_mass(i_part, :) / aero_data%density
2795  aero_particle%n_orig_part = aero_n_orig_part(i_part, :)
2796  aero_particle%weight_group = aero_particle_weight_group(i_part)
2797  aero_particle%weight_class = aero_particle_weight_class(i_part)
2798  if (size(aero_absorb_cross_sect) == n_part) then
2799  aero_particle%absorb_cross_sect = aero_absorb_cross_sect(i_part)
2800  end if
2801  if (size(aero_scatter_cross_sect) == n_part) then
2802  aero_particle%scatter_cross_sect = aero_scatter_cross_sect(i_part)
2803  end if
2804  if (size(aero_asymmetry) == n_part) then
2805  aero_particle%asymmetry = aero_asymmetry(i_part)
2806  end if
2807  if ((size(aero_refract_shell_real) == n_part) &
2808  .and. (size(aero_refract_shell_imag) == n_part)) then
2809  aero_particle%refract_shell = &
2810  cmplx(aero_refract_shell_real(i_part), &
2811  aero_refract_shell_imag(i_part), kind=dc)
2812  end if
2813  if ((size(aero_refract_core_real) == n_part) &
2814  .and. (size(aero_refract_core_imag) == n_part)) then
2815  aero_particle%refract_core = cmplx(aero_refract_core_real(i_part), &
2816  aero_refract_core_imag(i_part), kind=dc)
2817  end if
2818  if (size(aero_core_vol) == n_part) then
2819  aero_particle%core_vol = aero_core_vol(i_part)
2820  end if
2821  aero_particle%water_hyst_leg = aero_water_hyst_leg(i_part)
2822  aero_particle%id = aero_id(i_part)
2823  aero_particle%least_create_time = aero_least_create_time(i_part)
2824  aero_particle%greatest_create_time = aero_greatest_create_time(i_part)
2825 
2826  call assert(314368871, almost_equal(aero_num_conc(i_part), &
2827  aero_weight_array_num_conc(aero_state%awa, aero_particle, &
2828  aero_data)))
2829 
2830  call aero_state_add_particle(aero_state, aero_particle, aero_data)
2831  end do
2832 
2833  call pmc_nc_read_integer_1d(ncid, aero_removed_id, &
2834  "aero_removed_id", must_be_present=.false.)
2835  call pmc_nc_read_integer_1d(ncid, aero_removed_action, &
2836  "aero_removed_action", must_be_present=.false.)
2837  call pmc_nc_read_integer_1d(ncid, aero_removed_other_id, &
2838  "aero_removed_other_id", must_be_present=.false.)
2839 
2840  n_info_item = size(aero_removed_id)
2841  if (n_info_item >= 1) then
2842  if ((n_info_item > 1) &
2843  .or. ((n_info_item == 1) .and. (aero_removed_id(1) /= 0))) then
2844  call aero_info_array_enlarge_to(aero_state%aero_info_array, &
2845  n_info_item)
2846  do i_remove = 1,n_info_item
2847  aero_state%aero_info_array%aero_info(i_remove)%id &
2848  = aero_removed_id(i_remove)
2849  aero_state%aero_info_array%aero_info(i_remove)%action &
2850  = aero_removed_action(i_remove)
2851  aero_state%aero_info_array%aero_info(i_remove)%other_id &
2852  = aero_removed_other_id(i_remove)
2853  end do
2854  end if
2855  end if
2856 
2857  end subroutine aero_state_input_netcdf
2858 
2859 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2860 
2861  !> Sorts the particles if necessary.
2862  subroutine aero_state_sort(aero_state, aero_data, bin_grid, all_procs_same)
2863 
2864  !> Aerosol state to sort.
2865  type(aero_state_t), intent(inout) :: aero_state
2866  !> Aerosol data.
2867  type(aero_data_t), intent(in) :: aero_data
2868  !> Bin grid.
2869  type(bin_grid_t), optional, intent(in) :: bin_grid
2870  !> Whether all processors should use the same bin grid.
2871  logical, optional, intent(in) :: all_procs_same
2872 
2873  call aero_sorted_remake_if_needed(aero_state%aero_sorted, aero_state%apa, &
2874  aero_data, aero_state%valid_sort, size(aero_state%awa%weight, 1), &
2875  size(aero_state%awa%weight, 2), bin_grid, all_procs_same)
2876  aero_state%valid_sort = .true.
2877 
2878  end subroutine aero_state_sort
2879 
2880 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2881 
2882  !> Check that aerosol state data is consistent.
2883  subroutine aero_state_check(aero_state, aero_data)
2884 
2885  !> Aerosol state to check.
2886  type(aero_state_t), intent(in) :: aero_state
2887  !> Aerosol data.
2888  type(aero_data_t), intent(in) :: aero_data
2889 
2890  logical, parameter :: continue_on_error = .false.
2891 
2892  call aero_particle_array_check(aero_state%apa, aero_data, &
2893  continue_on_error)
2894  if (aero_state%valid_sort) then
2895  call aero_sorted_check(aero_state%aero_sorted, aero_state%apa, &
2896  aero_data, size(aero_state%awa%weight, 1), &
2897  size(aero_state%awa%weight, 2), continue_on_error)
2898  end if
2899 
2900  end subroutine aero_state_check
2901 
2902 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2903 
2904 end module pmc_aero_state
2905 
subroutine aero_state_add_particles(aero_state, aero_state_delta, aero_data)
aero_state += aero_state_delta, with the weight of aero_state left unchanged, so the new concentratio...
Definition: aero_state.F90:641
subroutine aero_state_sample(aero_state_from, aero_state_to, aero_data, sample_prob, removal_action)
Generates a random sample by removing particles from aero_state_from and adding them to aero_state_to...
Definition: aero_state.F90:880
integer, parameter aero_info_halved
Particle was removed due to halving of the aerosol population.
Definition: aero_info.F90:25
subroutine aero_state_input_netcdf(aero_state, ncid, aero_data)
Read full state.
elemental real(kind=dp) function aero_particle_radius(aero_particle, aero_data)
Total radius of the particle (m).
subroutine aero_state_halve(aero_state, i_group, i_class)
Remove approximately half of the particles in the given weight group.
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:77
subroutine aero_state_prepare_weight_for_add(aero_state, aero_data, i_group, i_class, n_add, allow_doubling, allow_halving)
Change the weight if necessary to ensure that the addition of about n_add computational particles wil...
Definition: aero_state.F90:667
integer function aero_state_total_particles_all_procs(aero_state, i_group, i_class)
Returns the total number of particles across all processes.
Definition: aero_state.F90:299
real(kind=dp) function aero_weight_array_single_num_conc(aero_weight_array, aero_particle, aero_data)
Compute the number concentration for a particle (m^{-3}).
subroutine aero_state_netcdf_dim_aero_particle(aero_state, ncid, dimid_aero_particle)
Write the aero particle dimension to the given NetCDF file if it is not already present and in any ca...
subroutine aero_state_to_binned_dry(bin_grid, aero_data, aero_state, aero_binned)
Does the same thing as aero_state_to_bin() but based on dry radius.
real(kind=dp) function entropy(p)
Compute the entropy of a probability mass function (non necessarily normalized).
Definition: util.F90:1630
elemental real(kind=dp) function aero_particle_mass(aero_particle, aero_data)
Total mass of the particle (kg).
subroutine aero_weight_array_check_monotonicity(aero_weight_array, monotone_increasing, monotone_decreasing)
Determine whether all weight functions in an array are monotone increasing, monotone decreasing...
subroutine aero_state_mix(aero_state, del_t, mix_timescale, aero_data, specify_prob_transfer)
Mix the aero_states between all processes. Currently uses a simple all-to-all diffusion.
integer, parameter aero_info_weight
Particle was removed due to adjustments in the particle&#39;s weighting function.
Definition: aero_info.F90:28
subroutine aero_state_scale_weight(aero_state, aero_data, i_group, i_class, weight_ratio, allow_doubling, allow_halving)
Scale the weighting of the given group/class by the given ratio, altering particle number as necessar...
Random number generators.
Definition: rand.F90:9
An input file with extra data for printing messages.
Definition: spec_file.F90:59
integer function rand_poisson(mean)
Generate a Poisson-distributed random number with the given mean.
Definition: rand.F90:253
elemental integer function aero_particle_array_n_part(aero_particle_array)
Return the current number of particles.
subroutine aero_state_sort(aero_state, aero_data, bin_grid, all_procs_same)
Sorts the particles if necessary.
real(kind=dp) function aero_particle_solute_radius(aero_particle, aero_data)
Returns the total solute radius (m).
real(kind=dp) function, dimension( aero_state_n_part(aero_state)) aero_state_mobility_diameters(aero_state, aero_data, env_state)
Returns the mobility diameters of all particles.
1-D array of particles, used by aero_state to store the particles.
elemental integer function aero_info_array_n_item(aero_info_array)
Return the current number of items.
real(kind=dp) function aero_particle_mobility_diameter(aero_particle, aero_data, env_state)
Mobility diameter of the particle (m).
subroutine aero_state_rand_particle(aero_state, i_part)
Choose a random particle from the aero_state.
Definition: aero_state.F90:784
subroutine aero_particle_set_source(aero_particle, i_source)
Sets the aerosol particle source.
integer, parameter aero_state_weight_flat_source
Flat weighting by source.
Definition: aero_state.F90:47
subroutine aero_state_dup_particle(aero_state, aero_data, i_part, n_part_mean, random_weight_group)
Add copies or remove a particle, with a given mean number of resulting particles. ...
Definition: aero_state.F90:457
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
The aero_weight_array_t structure and associated subroutines.
subroutine pmc_mpi_alltoall_integer(send, recv)
Does an all-to-all transfer of integers.
Definition: mpi.F90:1469
integer function pmc_mpi_pack_size_apa(val)
Determines the number of bytes required to pack the given value.
subroutine aero_info_array_add_aero_info(aero_info_array, aero_info)
Adds the given aero_info to the end of the array.
subroutine aero_state_rebalance(aero_state, aero_data, allow_doubling, allow_halving, initial_state_warning)
Double or halve the particle population in each weight group to maintain close to n_part_ideal partic...
subroutine aero_sorted_remake_if_needed(aero_sorted, aero_particle_array, aero_data, valid_sort, n_group, n_class, bin_grid, all_procs_same)
Remake a sorting if particles are getting too close to the edges.
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_diameters(aero_state, aero_data, include, exclude)
Returns the diameters of all particles.
Definition: aero_state.F90:989
real(kind=dp) elemental function sphere_vol2rad(v)
Convert mass-equivalent volume (m^3) to geometric radius (m) for spherical particles.
Definition: util.F90:240
subroutine aero_weight_array_set_nummass(aero_weight_array, n_class)
Allocates an aero_weight_array as joint flat/power-3 weightings..
integer, parameter aero_state_tag_scatter
MPI tag for scattering between processes.
Definition: aero_state.F90:36
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
subroutine aero_info_array_zero(aero_info_array)
Sets an aero_info_array to contain zero data.
real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:100
real(kind=dp) function aero_state_particle_num_conc(aero_state, aero_particle, aero_data)
The number concentration of a single particle (m^{-3}).
Definition: aero_state.F90:507
The aero_particle_array_t structure and assoicated subroutines.
subroutine warn_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition: util.F90:37
real(kind=dp) function aero_particle_crit_rel_humid(aero_particle, aero_data, env_state)
Returns the critical relative humidity (1).
real(kind=dp) function, dimension(aero_data_n_spec(aero_data)) aero_particle_species_masses(aero_particle, aero_data)
Mass of all species in the particle (kg).
elemental integer function integer_varray_n_entry(integer_varray)
Return the current number of entries.
The aero_dist_t structure and associated subroutines.
Definition: aero_dist.F90:18
subroutine aero_info_array_add(aero_info_array, aero_info_array_delta)
Adds aero_info_array_delta to the end of aero_info_array.
integer function pmc_mpi_pack_size_aero_weight_array(val)
Determines the number of bytes required to pack the given value.
subroutine aero_data_netcdf_dim_aero_species(aero_data, ncid, dimid_aero_species)
Write the aero species dimension to the given NetCDF file if it is not already present and in any cas...
Definition: aero_data.F90:550
integer, parameter aero_state_tag_gather
MPI tag for gathering between processes.
Definition: aero_state.F90:34
subroutine pmc_mpi_pack_aero_info_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
subroutine aero_particle_array_remove_particle(aero_particle_array, index)
Removes the particle at the given index.
real(kind=dp) function aero_particle_approx_crit_rel_humid(aero_particle, aero_data, env_state)
Returns the approximate critical relative humidity (1).
subroutine aero_particle_set_weight(aero_particle, i_group, i_class)
Sets the aerosol particle weight group.
subroutine aero_state_add_aero_dist_sample(aero_state, aero_data, aero_dist, sample_prop, create_time, allow_doubling, allow_halving, n_part_add)
Generates a Poisson sample of an aero_dist, adding to aero_state, with the given sample proportion...
Definition: aero_state.F90:713
subroutine aero_state_mpi_gather(aero_state, aero_state_total, aero_data)
Gathers data from all processes into one aero_state on process 0.
The aero_particle_t structure and associated subroutines.
elemental real(kind=dp) function aero_particle_species_volume(aero_particle, i_spec)
Volume of a single species in the particle (m^3).
subroutine aero_info_array_enlarge_to(aero_info_array, n)
Possibly enlarges the given array, ensuring that it is at least of size n.
elemental integer function aero_state_n_part(aero_state)
Return the current number of particles.
Definition: aero_state.F90:84
subroutine aero_particle_new_id(aero_particle)
Assigns a globally-unique new ID number to the particle.
The aero_sorted_t structure and assocated subroutines.
Definition: aero_sorted.F90:9
subroutine pmc_mpi_unpack_aero_weight_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
subroutine aero_weight_array_input_netcdf(aero_weight_array, ncid)
Read full aero_weight_array.
real(kind=dp) function aero_weight_array_num_conc_at_radius(aero_weight_array, i_class, radius)
Compute the total number concentration at a given radius (m^3).
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:139
subroutine spec_file_read_real(file, name, var)
Read a real number from a spec file that must have the given name.
Definition: spec_file.F90:562
integer function pmc_mpi_pack_size_aero_state(val)
Determines the number of bytes required to pack the given value.
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
subroutine aero_particle_array_zero(aero_particle_array)
Resets an aero_particle_array to contain zero particles.
subroutine aero_particle_set_create_time(aero_particle, create_time)
Sets the creation times for the particle.
integer function pmc_mpi_pack_size_real_array_2d(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:530
subroutine aero_state_to_binned(bin_grid, aero_data, aero_state, aero_binned)
Create binned number and mass arrays.
Definition: aero_state.F90:931
subroutine aero_state_set_n_part_ideal(aero_state, n_part)
Set the ideal number of particles to the given value. The aero_state%awa must be already set correctl...
Definition: aero_state.F90:214
A complete aerosol distribution, consisting of several modes.
Definition: aero_dist.F90:33
subroutine aero_state_double(aero_state, aero_data, i_group, i_class)
Doubles number of particles in the given weight group.
1-D arrays of aero_info_t structure.
subroutine aero_particle_array_check(aero_particle_array, aero_data, continue_on_error)
Check that the particle array data is consistent.
subroutine aero_state_copy_weight(aero_state_from, aero_state_to)
Copies weighting information for an aero_state.
Definition: aero_state.F90:154
subroutine aero_state_reweight(aero_state, aero_data, reweight_num_conc)
Reweight all particles after their constituent volumes have been altered.
Definition: aero_state.F90:559
subroutine aero_state_remove_particle_no_info(aero_state, i_part)
Remove the given particle without recording it.
Definition: aero_state.F90:357
integer, parameter aero_state_tag_mix
MPI tag for mixing particles between processes.
Definition: aero_state.F90:32
subroutine aero_weight_array_set_flat(aero_weight_array, n_class)
Allocates an aero_weight_array as flat weightings.
subroutine aero_state_remove_particle_with_info(aero_state, i_part, aero_info)
Remove the given particle and record the removal.
Definition: aero_state.F90:377
real(kind=dp) elemental function rad2diam(r)
Convert radius (m) to diameter (m).
Definition: util.F90:252
An array of aerosol size distribution weighting functions.
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_num_concs(aero_state, aero_data)
Returns the number concentrations of all particles.
real(kind=dp) function aero_weight_array_num_conc(aero_weight_array, aero_particle, aero_data)
Compute the number concentration for a particle (m^{-3}).
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:134
integer function aero_state_total_particles(aero_state, i_group, i_class)
Returns the total number of particles in an aerosol distribution.
Definition: aero_state.F90:260
subroutine aero_weight_array_set_power(aero_weight_array, n_class, exponent)
Allocates an aero_weight_array as power weightings.
elemental integer function bin_grid_size(bin_grid)
Return the number of bins in the grid, or -1 if the bin grid is not allocated.
Definition: bin_grid.F90:51
integer, parameter aero_info_none
No information.
Definition: aero_info.F90:19
subroutine pmc_mpi_pack_aero_state(buffer, position, val)
Packs the given value into the buffer, advancing position.
subroutine aero_sorted_check(aero_sorted, aero_particle_array, aero_data, n_group, n_class, continue_on_error)
Check sorting.
subroutine aero_particle_set_vols(aero_particle, vols)
Sets the aerosol particle volumes.
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_masses(aero_state, aero_data, include, exclude)
Returns the masses of all particles.
The current collection of aerosol particles.
Definition: aero_state.F90:63
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:74
integer function pmc_rand_int(n)
Returns a random integer between 1 and n.
Definition: rand.F90:173
subroutine pmc_mpi_unpack_aero_state(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
The integer_varray_t structure and assocated subroutines.
elemental subroutine aero_weight_scale(aero_weight, factor)
Scale the weight by the given fraction, so new_weight = old_weight * factor.
Definition: aero_weight.F90:61
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
subroutine aero_state_netcdf_dim_aero_removed(aero_state, ncid, dimid_aero_removed)
Write the aero removed dimension to the given NetCDF file if it is not already present and in any cas...
integer, parameter aero_state_weight_nummass_source
Coupled number/mass weighting by source.
Definition: aero_state.F90:51
Single aerosol particle data structure.
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
subroutine aero_state_bin_average_size(aero_state, bin_grid, aero_data, bin_center, preserve_number)
Set each aerosol particle to have its original species ratios, but total volume given by the average ...
subroutine aero_state_make_dry(aero_state, aero_data)
Make all particles dry (water set to zero).
subroutine aero_state_remove_rand_particle_from_bin(aero_state, i_bin, i_class, aero_particle)
Remove a randomly chosen particle from the given bin and return it.
Definition: aero_state.F90:421
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_dry_diameters(aero_state, aero_data)
Returns the dry diameters of all particles.
subroutine bin_grid_make(bin_grid, type, n_bin, min, max)
Generates the bin grid given the range and number of bins.
Definition: bin_grid.F90:84
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_approx_crit_rel_humids(aero_state, aero_data, env_state)
Returns the approximate critical relative humidity for all particles (1).
subroutine aero_weight_array_combine(aero_weight_array, aero_weight_array_delta)
Combine aero_weight_array_delta into aero_weight_array with a harmonic mean.
integer function aero_state_weight_class_for_source(aero_state, source)
Determine the appropriate weight class for a source.
Definition: aero_state.F90:236
subroutine aero_state_mixing_state_metrics(aero_state, aero_data, d_alpha, d_gamma, chi, include, exclude, group)
Returns the mixing state metrics of the population.
Reading formatted text input.
Definition: spec_file.F90:43
subroutine aero_state_sample_particles(aero_state_from, aero_state_to, aero_data, sample_prob, removal_action)
Generates a random sample by removing particles from aero_state_from and adding them to aero_state_to...
Definition: aero_state.F90:805
integer, parameter aero_state_weight_nummass
Coupled number/mass weighting scheme.
Definition: aero_state.F90:45
integer, parameter aero_state_weight_power
Power-law weighting scheme.
Definition: aero_state.F90:43
integer function aero_weight_array_rand_group(aero_weight_array, i_class, radius)
Choose a random group at the given radius, with probability inversely proportional to group weight at...
subroutine pmc_mpi_pack_aero_weight_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
logical function aero_weight_array_check_flat(aero_weight_array)
Check whether a given aero_weight array is flat in total.
integer function bin_grid_find(bin_grid, val)
Find the bin number that contains a given value.
Definition: bin_grid.F90:144
integer function prob_round(val)
Round val to floor(val) or ceiling(val) with probability proportional to the relative distance from v...
Definition: rand.F90:215
subroutine aero_sorted_remove_particle(aero_sorted, aero_particle_array, i_part)
Remove a particle from both an aero_sorted and the corresponding aero_particle_array.
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_mass_entropies(aero_state, aero_data, include, exclude, group)
Returns the mass-entropies of all particles.
elemental real(kind=dp) function aero_particle_volume(aero_particle)
Total volume of the particle (m^3).
Sorting of particles into bins.
Definition: aero_sorted.F90:47
subroutine pmc_mpi_check_ierr(ierr)
Dies if ierr is not ok.
Definition: mpi.F90:40
subroutine aero_state_num_conc_for_reweight(aero_state, aero_data, reweight_num_conc)
Save the correct number concentrations for later use by aero_state_reweight().
Definition: aero_state.F90:527
subroutine aero_particle_zero(aero_particle, aero_data)
Resets an aero_particle to be zero.
elemental subroutine aero_weight_array_normalize(aero_weight_array)
Normalizes the aero_weight_array to a non-zero value.
subroutine aero_data_netcdf_dim_aero_source(aero_data, ncid, dimid_aero_source)
Write the aero source dimension to the given NetCDF file if it is not already present and in any case...
Definition: aero_data.F90:607
subroutine aero_state_add(aero_state, aero_state_delta, aero_data)
aero_state += aero_state_delta, including combining the weights, so the new concentration is the weig...
Definition: aero_state.F90:622
subroutine aero_state_zero(aero_state)
Resets an aero_state to have zero particles per bin.
Definition: aero_state.F90:317
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
subroutine pmc_mpi_unpack_real_array_2d(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1042
integer function pmc_mpi_pack_size_aia(val)
Determines the number of bytes required to pack the given value.
subroutine pmc_mpi_pack_aero_particle_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
integer, parameter aero_state_weight_flat
Single flat weighting scheme.
Definition: aero_state.F90:41
integer, parameter aero_state_weight_none
Single flat weighting scheme.
Definition: aero_state.F90:39
The aero_info_t structure and associated subroutines.
Definition: aero_info.F90:9
subroutine pmc_mpi_unpack_aero_particle_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
elemental integer function aero_data_n_source(aero_data)
Return the number of aerosol sources, or -1 if uninitialized.
Definition: aero_data.F90:230
The aero_info_array_t structure and assoicated subroutines.
subroutine spec_file_read_string(file, name, var)
Read a string from a spec file that must have a given name.
Definition: spec_file.F90:604
subroutine aero_state_bin_average_comp(aero_state, bin_grid, aero_data)
Set each aerosol particle to have its original total volume, but species volume ratios given by the t...
elemental integer function aero_dist_n_mode(aero_dist)
Return the number of modes.
Definition: aero_dist.F90:44
integer, parameter bin_grid_type_log
Logarithmically spaced bin grid.
Definition: bin_grid.F90:23
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Definition: aero_data.F90:214
subroutine aero_state_check(aero_state, aero_data)
Check that aerosol state data is consistent.
integer function aero_weight_array_n_group(aero_weight_array)
Return the number of weight groups.
subroutine pmc_mpi_pack_real_array_2d(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:784
subroutine aero_state_set_weight(aero_state, aero_data, weight_type, exponent)
Sets the weighting functions for an aero_state.
Definition: aero_state.F90:169
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:766
integer, parameter aero_state_weight_power_source
Power-law weighting by source.
Definition: aero_state.F90:49
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_crit_rel_humids(aero_state, aero_data, env_state)
Returns the critical relative humidity for all particles (1).
integer function aero_weight_array_n_class(aero_weight_array)
Return the number of weight classes.
subroutine aero_weight_array_shift(aero_weight_array_from, aero_weight_array_to, sample_prop, overwrite_to)
Adjust source and destination weights to reflect moving sample_prop proportion of particles from aero...
integer function pmc_mpi_size()
Returns the total number of processes.
Definition: mpi.F90:134
Aerosol material properties and associated data.
Definition: aero_data.F90:41
subroutine aero_state_add_particle(aero_state, aero_particle, aero_data, allow_resort)
Add the given particle.
Definition: aero_state.F90:334
subroutine aero_state_mpi_alltoall(send, recv)
Do an MPI all-to-all transfer of aerosol states.
Common utility subroutines.
Definition: util.F90:9
The aero_binned_t structure and associated subroutines.
Definition: aero_binned.F90:9
real(kind=dp) function aero_state_total_num_conc(aero_state, aero_data)
Returns the total number concentration.
Wrapper functions for MPI.
Definition: mpi.F90:13
integer function, dimension(aero_state_n_part(aero_state)) aero_state_ids(aero_state)
Returns the IDs of all particles.
Definition: aero_state.F90:970
integer, parameter dc
Kind of a double precision complex number.
Definition: constants.F90:14
Information about removed particles describing the sink.
Definition: aero_info.F90:48
integer function rand_binomial(n, p)
Generate a Binomial-distributed random number with the given parameters.
Definition: rand.F90:320
subroutine aero_state_remove_particle(aero_state, i_part, record_removal, aero_info)
Remove the given particle and possibly record the removal.
Definition: aero_state.F90:396
The aero_weight_t structure and associated subroutines.
Definition: aero_weight.F90:9
elemental real(kind=dp) function aero_particle_dry_diameter(aero_particle, aero_data)
Total dry diameter of the particle (m).
integer function aero_data_spec_by_name(aero_data, name)
Returns the number of the species in aero_data with the given name, or returns 0 if there is no such ...
Definition: aero_data.F90:247
subroutine aero_sorted_add_particle(aero_sorted, aero_particle_array, aero_particle, aero_data, allow_resort)
Add a new particle to both an aero_sorted and the corresponding aero_particle_array.
subroutine pmc_mpi_allreduce_sum_integer(val, val_sum)
Computes the sum of val across all processes, storing the result in val_sum on all processes...
Definition: mpi.F90:1203
elemental real(kind=dp) function aero_particle_species_mass(aero_particle, i_spec, aero_data)
Mass of a single species in the particle (kg).
Aerosol number and volume distributions stored per bin.
Definition: aero_binned.F90:37
subroutine pmc_mpi_unpack_aero_info_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
subroutine aero_particle_array_add_particle(aero_particle_array, aero_particle)
Adds the given particle to the end of the array.
logical function almost_equal(d1, d2)
Tests whether two real numbers are almost equal using only a relative tolerance.
Definition: util.F90:318
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_volumes(aero_state, aero_data, include, exclude)
Returns the volumes of all particles.