PartMC  2.4.0
aero_state.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2017 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)
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 
995  !> Return diameters array (m).
996  real(kind=dp) :: aero_state_diameters(aero_state_n_part(aero_state))
997 
998  aero_state_diameters = aero_particle_diameter( &
999  aero_state%apa%particle(1:aero_state_n_part(aero_state)), aero_data)
1000 
1001  end function aero_state_diameters
1002 
1003 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1004 
1005  !> Returns the dry diameters of all particles.
1006  function aero_state_dry_diameters(aero_state, aero_data)
1008  !> Aerosol state.
1009  type(aero_state_t), intent(in) :: aero_state
1010  !> Aerosol data.
1011  type(aero_data_t), intent(in) :: aero_data
1012 
1013  !> Return value (m).
1014  real(kind=dp) :: aero_state_dry_diameters(aero_state_n_part(aero_state))
1015 
1016  aero_state_dry_diameters = aero_particle_dry_diameter( &
1017  aero_state%apa%particle(1:aero_state_n_part(aero_state)), aero_data)
1018 
1019  end function aero_state_dry_diameters
1020 
1021 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1022 
1023  !> Returns the mobility diameters of all particles.
1024  function aero_state_mobility_diameters(aero_state, aero_data, env_state)
1026  !> Aerosol state.
1027  type(aero_state_t), intent(in) :: aero_state
1028  !> Aerosol data.
1029  type(aero_data_t), intent(in) :: aero_data
1030  !> Environment state.
1031  type(env_state_t), intent(in) :: env_state
1032 
1033  !> Return value (m).
1034  real(kind=dp) :: aero_state_mobility_diameters( &
1035  aero_state_n_part(aero_state))
1036 
1037  integer :: i_part
1038 
1039  do i_part = 1,aero_state_n_part(aero_state)
1040  aero_state_mobility_diameters(i_part) &
1042  aero_state%apa%particle(i_part), aero_data, env_state)
1043  end do
1044 
1045  end function aero_state_mobility_diameters
1046 
1047 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1048 
1049  !> Returns the volumes of all particles.
1050  !!
1051  !! If \c include is specified then only those species are included
1052  !! in computing the volumes. If \c exclude is specified then all
1053  !! species except those species are included. If both \c include and
1054  !! \c exclude arguments are specified then only those species in \c
1055  !! include but not in \c exclude are included.
1056  function aero_state_volumes(aero_state, aero_data, include, exclude)
1058  !> Aerosol state.
1059  type(aero_state_t), intent(in) :: aero_state
1060  !> Aerosol data.
1061  type(aero_data_t), optional, intent(in) :: aero_data
1062  !> Species names to include in the mass.
1063  character(len=*), optional, intent(in) :: include(:)
1064  !> Species names to exclude in the mass.
1065  character(len=*), optional, intent(in) :: exclude(:)
1066 
1067  !> Return volumes array (m^3).
1068  real(kind=dp) :: aero_state_volumes(aero_state_n_part(aero_state))
1069 
1070  logical, allocatable :: use_species(:)
1071  integer :: i_name, i_spec
1072 
1073  if ((.not. present(include)) .and. (.not. present(exclude))) then
1074  aero_state_volumes = aero_particle_volume( &
1075  aero_state%apa%particle(1:aero_state_n_part(aero_state)))
1076  else
1077  call assert_msg(599558703, present(aero_data), &
1078  "must provide 'aero_data' if using 'include' or 'exclude'")
1079  allocate(use_species(aero_data_n_spec(aero_data)))
1080  if (present(include)) then
1081  use_species = .false.
1082  do i_name = 1, size(include)
1083  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1084  call assert_msg(111852070, i_spec > 0, &
1085  "unknown species: " // trim(include(i_name)))
1086  use_species(i_spec) = .true.
1087  end do
1088  else
1089  use_species = .true.
1090  end if
1091  if (present(exclude)) then
1092  do i_name = 1, size(exclude)
1093  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1094  call assert_msg(182075590, i_spec > 0, &
1095  "unknown species: " // trim(exclude(i_name)))
1096  use_species(i_spec) = .false.
1097  end do
1098  end if
1099  aero_state_volumes = 0d0
1100  do i_spec = 1,aero_data_n_spec(aero_data)
1101  if (use_species(i_spec)) then
1102  aero_state_volumes = aero_state_volumes &
1104  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1105  i_spec)
1106  end if
1107  end do
1108  end if
1109 
1110  end function aero_state_volumes
1111 
1112 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1113 
1114  !> Returns the masses of all particles.
1115  !!
1116  !! If \c include is specified then only those species are included
1117  !! in computing the masses. If \c exclude is specified then all
1118  !! species except those species are included. If both \c include and
1119  !! \c exclude arguments are specified then only those species in \c
1120  !! include but not in \c exclude are included.
1121  function aero_state_masses(aero_state, aero_data, include, exclude)
1123  !> Aerosol state.
1124  type(aero_state_t), intent(in) :: aero_state
1125  !> Aerosol data.
1126  type(aero_data_t), intent(in) :: aero_data
1127  !> Species names to include in the mass.
1128  character(len=*), optional, intent(in) :: include(:)
1129  !> Species names to exclude in the mass.
1130  character(len=*), optional, intent(in) :: exclude(:)
1131 
1132  !> Return masses array (kg).
1133  real(kind=dp) :: aero_state_masses(aero_state_n_part(aero_state))
1134 
1135  logical :: use_species(aero_data_n_spec(aero_data))
1136  integer :: i_name, i_spec
1137 
1138  if ((.not. present(include)) .and. (.not. present(exclude))) then
1139  aero_state_masses = aero_particle_mass( &
1140  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1141  aero_data)
1142  else
1143  if (present(include)) then
1144  use_species = .false.
1145  do i_name = 1, size(include)
1146  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1147  call assert_msg(963163690, i_spec > 0, &
1148  "unknown species: " // trim(include(i_name)))
1149  use_species(i_spec) = .true.
1150  end do
1151  else
1152  use_species = .true.
1153  end if
1154  if (present(exclude)) then
1155  do i_name = 1, size(exclude)
1156  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1157  call assert_msg(950847713, i_spec > 0, &
1158  "unknown species: " // trim(exclude(i_name)))
1159  use_species(i_spec) = .false.
1160  end do
1161  end if
1162  aero_state_masses = 0d0
1163  do i_spec = 1,aero_data_n_spec(aero_data)
1164  if (use_species(i_spec)) then
1165  aero_state_masses = aero_state_masses &
1167  aero_state%apa%particle(1:aero_state_n_part(aero_state)), &
1168  i_spec, aero_data)
1169  end if
1170  end do
1171  end if
1172 
1173  end function aero_state_masses
1174 
1175 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1176 
1177  !> Returns the number concentrations of all particles.
1178  function aero_state_num_concs(aero_state, aero_data)
1180  !> Aerosol state.
1181  type(aero_state_t), intent(in) :: aero_state
1182  !> Aerosol data.
1183  type(aero_data_t), intent(in) :: aero_data
1184 
1185  !> Return number concentrations array (m^{-3}).
1186  real(kind=dp) :: aero_state_num_concs(aero_state_n_part(aero_state))
1187 
1188  integer :: i_part
1189 
1190  do i_part = 1,aero_state_n_part(aero_state)
1191  aero_state_num_concs(i_part) &
1192  = aero_state_particle_num_conc(aero_state, &
1193  aero_state%apa%particle(i_part), aero_data)
1194  end do
1195 
1196  end function aero_state_num_concs
1197 
1198 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1199 
1200  !> Returns the total number concentration.
1201  real(kind=dp) function aero_state_total_num_conc(aero_state, aero_data)
1203  !> Aerosol state.
1204  type(aero_state_t), intent(in) :: aero_state
1205  !> Aerosol data.
1206  type(aero_data_t), intent(in) :: aero_data
1207 
1208  integer :: i_part
1209 
1211  do i_part = 1,aero_state_n_part(aero_state)
1213  + aero_state_particle_num_conc(aero_state, &
1214  aero_state%apa%particle(i_part), aero_data)
1215  end do
1216 
1217  end function aero_state_total_num_conc
1218 
1219 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1220 
1221  !> Returns the mass-entropies of all particles.
1222  !!
1223  !! If \c include is specified then only those species are included
1224  !! in computing the entropy. If \c exclude is specified then all
1225  !! species except those species are included. If both \c include and
1226  !! \c exclude arguments are specified then only those species in \c
1227  !! include but not in \c exclude are included. If \c group is
1228  !! specified then the species are divided into two sets, given by
1229  !! those in the group and those not in the group. The entropy is
1230  !! then computed using the total mass of each set.
1231  function aero_state_mass_entropies(aero_state, aero_data, include, exclude, &
1232  group)
1234  !> Aerosol state.
1235  type(aero_state_t), intent(in) :: aero_state
1236  !> Aerosol data.
1237  type(aero_data_t), intent(in) :: aero_data
1238  !> Species names to include in the mass.
1239  character(len=*), optional :: include(:)
1240  !> Species names to exclude in the mass.
1241  character(len=*), optional :: exclude(:)
1242  !> Species names to group together.
1243  character(len=*), optional :: group(:)
1244 
1245  !> Return value.
1246  real(kind=dp) :: aero_state_mass_entropies(aero_state_n_part(aero_state))
1247 
1248  logical :: use_species(aero_data_n_spec(aero_data))
1249  logical :: group_species(aero_data_n_spec(aero_data))
1250  integer :: i_name, i_spec, i_part
1251  real(kind=dp) :: group_mass, non_group_mass, mass
1252 
1253  if (present(include)) then
1254  use_species = .false.
1255  do i_name = 1, size(include)
1256  i_spec = aero_data_spec_by_name(aero_data, include(i_name))
1257  call assert_msg(890212002, i_spec > 0, &
1258  "unknown species: " // trim(include(i_name)))
1259  use_species(i_spec) = .true.
1260  end do
1261  else
1262  use_species = .true.
1263  end if
1264  if (present(exclude)) then
1265  do i_name = 1, size(exclude)
1266  i_spec = aero_data_spec_by_name(aero_data, exclude(i_name))
1267  call assert_msg(859945006, i_spec > 0, &
1268  "unknown species: " // trim(exclude(i_name)))
1269  use_species(i_spec) = .false.
1270  end do
1271  end if
1272  if (present(group)) then
1273  group_species = .false.
1274  do i_name = 1, size(group)
1275  i_spec = aero_data_spec_by_name(aero_data, group(i_name))
1276  call assert_msg(376359046, i_spec > 0, &
1277  "unknown species: " // trim(group(i_name)))
1278  group_species(i_spec) = .true.
1279  end do
1280  do i_part = 1,aero_state_n_part(aero_state)
1281  group_mass = 0d0
1282  non_group_mass = 0d0
1283  do i_spec = 1,aero_data_n_spec(aero_data)
1284  if (use_species(i_spec)) then
1285  mass = aero_particle_species_mass( &
1286  aero_state%apa%particle(i_part), i_spec, aero_data)
1287  if (group_species(i_spec)) then
1288  group_mass = group_mass + mass
1289  else
1290  non_group_mass = non_group_mass + mass
1291  end if
1292  end if
1293  end do
1294  aero_state_mass_entropies(i_part) &
1295  = entropy([group_mass, non_group_mass])
1296  end do
1297  else
1298  do i_part = 1,aero_state_n_part(aero_state)
1299  aero_state_mass_entropies(i_part) = entropy(pack( &
1300  aero_particle_species_masses(aero_state%apa%particle(i_part), &
1301  aero_data), use_species))
1302  end do
1303  end if
1304 
1305  end function aero_state_mass_entropies
1306 
1307 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1308 
1309  !> Returns the mixing state metrics of the population.
1310  !!
1311  !! If \c include is specified then only those species are included
1312  !! in computing the entropies. If \c exclude is specified then all
1313  !! species except those species are included. If both \c include and
1314  !! \c exclude arguments are specified then only those species in \c
1315  !! include but not in \c exclude are included. If \c group is
1316  !! specified then the species are divided into two sets, given by
1317  !! those in the group and those not in the group. The entropies are
1318  !! then computed using the total mass of each set.
1319  subroutine aero_state_mixing_state_metrics(aero_state, aero_data, d_alpha, &
1320  d_gamma, chi, include, exclude, group)
1322  !> Aerosol state.
1323  type(aero_state_t), intent(in) :: aero_state
1324  !> Aerosol data.
1325  type(aero_data_t), intent(in) :: aero_data
1326  !> Average particle diversity.
1327  real(kind=dp), intent(out) :: d_alpha
1328  !> Bulk diversity.
1329  real(kind=dp), intent(out) :: d_gamma
1330  !> Mixing state index.
1331  real(kind=dp), intent(out) :: chi
1332  !> Species names to include in the mass.
1333  character(len=*), optional :: include(:)
1334  !> Species names to exclude in the mass.
1335  character(len=*), optional :: exclude(:)
1336  !> Species names to group together.
1337  character(len=*), optional :: group(:)
1338 
1339  real(kind=dp), allocatable :: entropies(:), entropies_of_avg_part(:)
1340  real(kind=dp), allocatable :: masses(:), num_concs(:), &
1341  num_concs_of_avg_part(:), masses_of_avg_part(:)
1342  type(aero_state_t) :: aero_state_averaged
1343  type(bin_grid_t) :: avg_bin_grid
1344 
1345  ! per-particle properties
1346  num_concs = aero_state_num_concs(aero_state, aero_data)
1347  masses = aero_state_masses(aero_state, aero_data, include, exclude)
1348  entropies = aero_state_mass_entropies(aero_state, aero_data, &
1349  include, exclude, group)
1350 
1351  d_alpha = exp(sum(entropies * masses * num_concs) &
1352  / sum(masses * num_concs))
1353 
1354  ! per-particle properties of averaged particles
1355  call bin_grid_make(avg_bin_grid, bin_grid_type_log, 1, 1d-30, 1d10)
1356  aero_state_averaged = aero_state
1357  call aero_state_bin_average_comp(aero_state_averaged, avg_bin_grid, &
1358  aero_data)
1359  num_concs_of_avg_part = aero_state_num_concs(aero_state_averaged, &
1360  aero_data)
1361  masses_of_avg_part = aero_state_masses(aero_state_averaged, &
1362  aero_data, include, exclude)
1363  entropies_of_avg_part = aero_state_mass_entropies(aero_state_averaged, &
1364  aero_data, include, exclude, group)
1365 
1366  d_gamma = exp(sum(entropies_of_avg_part * masses_of_avg_part &
1367  * num_concs_of_avg_part) &
1368  / sum(masses_of_avg_part * num_concs_of_avg_part))
1369 
1370  chi = (d_alpha - 1) / (d_gamma - 1)
1371 
1372  end subroutine aero_state_mixing_state_metrics
1373 
1374 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1375 
1376  !> Returns the approximate critical relative humidity for all particles (1).
1377  function aero_state_approx_crit_rel_humids(aero_state, aero_data, env_state)
1379  !> Aerosol state.
1380  type(aero_state_t), intent(in) :: aero_state
1381  !> Aerosol data.
1382  type(aero_data_t), intent(in) :: aero_data
1383  !> Environment state.
1384  type(env_state_t), intent(in) :: env_state
1385 
1386  !> Return value.
1387  real(kind=dp) &
1388  :: aero_state_approx_crit_rel_humids(aero_state_n_part(aero_state))
1389 
1390  integer :: i_part
1391 
1392  do i_part = 1,aero_state_n_part(aero_state)
1393  aero_state_approx_crit_rel_humids(i_part) = &
1395  aero_state%apa%particle(i_part), aero_data, env_state)
1396  end do
1397 
1399 
1400 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1401 
1402  !> Returns the critical relative humidity for all particles (1).
1403  function aero_state_crit_rel_humids(aero_state, aero_data, env_state)
1405  !> Aerosol state.
1406  type(aero_state_t), intent(in) :: aero_state
1407  !> Aerosol data.
1408  type(aero_data_t), intent(in) :: aero_data
1409  !> Environment state.
1410  type(env_state_t), intent(in) :: env_state
1411 
1412  !> Return value.
1413  real(kind=dp) :: aero_state_crit_rel_humids(aero_state_n_part(aero_state))
1414 
1415  integer :: i_part
1416 
1417  do i_part = 1,aero_state_n_part(aero_state)
1418  aero_state_crit_rel_humids(i_part) = aero_particle_crit_rel_humid( &
1419  aero_state%apa%particle(i_part), aero_data, env_state)
1420  end do
1421 
1422  end function aero_state_crit_rel_humids
1423 
1424 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1425 
1426  !> Does the same thing as aero_state_to_bin() but based on dry radius.
1427  subroutine aero_state_to_binned_dry(bin_grid, aero_data, aero_state, &
1428  aero_binned)
1430  !> Bin grid.
1431  type(bin_grid_t), intent(in) :: bin_grid
1432  !> Aerosol data.
1433  type(aero_data_t), intent(in) :: aero_data
1434  !> Aerosol state.
1435  type(aero_state_t), intent(in) :: aero_state
1436  !> Binned distributions.
1437  type(aero_binned_t), intent(inout) :: aero_binned
1438 
1439  integer :: i_part, i_bin
1440  real(kind=dp) :: factor
1441 
1442  aero_binned%num_conc = 0d0
1443  aero_binned%vol_conc = 0d0
1444  do i_part = 1,aero_state_n_part(aero_state)
1445  i_bin = bin_grid_find(bin_grid, &
1446  aero_particle_solute_radius(aero_state%apa%particle(i_part), &
1447  aero_data))
1448  if ((i_bin < 1) .or. (i_bin > bin_grid_size(bin_grid))) then
1449  call warn_msg(503871022, "particle ID " &
1450  // trim(integer_to_string(aero_state%apa%particle(i_part)%id)) &
1451  // " outside of bin_grid, discarding")
1452  else
1453  factor = aero_weight_array_num_conc(aero_state%awa, &
1454  aero_state%apa%particle(i_part), aero_data) &
1455  / bin_grid%widths(i_bin)
1456  aero_binned%vol_conc(i_bin,:) = aero_binned%vol_conc(i_bin,:) &
1457  + aero_state%apa%particle(i_part)%vol * factor
1458  aero_binned%num_conc(i_bin) = aero_binned%num_conc(i_bin) &
1459  + factor
1460  end if
1461  end do
1462 
1463  end subroutine aero_state_to_binned_dry
1464 
1465 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1466 
1467  !> Doubles number of particles in the given weight group.
1468  subroutine aero_state_double(aero_state, aero_data, i_group, i_class)
1470  !> Aerosol state.
1471  type(aero_state_t), intent(inout) :: aero_state
1472  !> Aerosol data.
1473  type(aero_data_t), intent(in) :: aero_data
1474  !> Weight group to double.
1475  integer, intent(in) :: i_group
1476  !> Weight class to double.
1477  integer, intent(in) :: i_class
1478 
1479  integer :: i_part
1480  type(aero_particle_t) :: aero_particle
1481 
1482  do i_part = 1,aero_state_n_part(aero_state)
1483  if ((aero_state%apa%particle(i_part)%weight_group == i_group) &
1484  .and. (aero_state%apa%particle(i_part)%weight_class == i_class)) &
1485  then
1486  aero_particle = aero_state%apa%particle(i_part)
1487  call aero_particle_new_id(aero_particle)
1488  call aero_state_add_particle(aero_state, aero_particle, aero_data)
1489  end if
1490  end do
1491  aero_state%valid_sort = .false.
1492  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), 0.5d0)
1493 
1494  end subroutine aero_state_double
1495 
1496 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1497 
1498  !> Remove approximately half of the particles in the given weight group.
1499  subroutine aero_state_halve(aero_state, i_group, i_class)
1501  !> Aerosol state.
1502  type(aero_state_t), intent(inout) :: aero_state
1503  !> Weight group to halve.
1504  integer, intent(in) :: i_group
1505  !> Weight class to halve.
1506  integer, intent(in) :: i_class
1507 
1508  integer :: i_part
1509  type(aero_info_t) :: aero_info
1510 
1511  do i_part = aero_state_n_part(aero_state),1,-1
1512  if ((aero_state%apa%particle(i_part)%weight_group == i_group) &
1513  .and. (aero_state%apa%particle(i_part)%weight_class == i_class)) &
1514  then
1515  if (pmc_random() < 0.5d0) then
1516  aero_info%id = aero_state%apa%particle(i_part)%id
1517  aero_info%action = aero_info_halved
1518  call aero_state_remove_particle_with_info(aero_state, i_part, &
1519  aero_info)
1520  end if
1521  end if
1522  end do
1523  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), 2d0)
1524 
1525  end subroutine aero_state_halve
1526 
1527 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1528 
1529  !> Double or halve the particle population in each weight group to
1530  !> maintain close to \c n_part_ideal particles per process,
1531  !> allocated equally amongst the weight groups.
1532  subroutine aero_state_rebalance(aero_state, aero_data, allow_doubling, &
1533  allow_halving, initial_state_warning)
1535  !> Aerosol state.
1536  type(aero_state_t), intent(inout) :: aero_state
1537  !> Aerosol data.
1538  type(aero_data_t), intent(in) :: aero_data
1539  !> Whether to allow doubling of the population.
1540  logical, intent(in) :: allow_doubling
1541  !> Whether to allow halving of the population.
1542  logical, intent(in) :: allow_halving
1543  !> Whether to warn due to initial state doubling/halving.
1544  logical, intent(in) :: initial_state_warning
1545 
1546  integer :: i_group, i_class, n_group, n_class, global_n_part
1547 
1548  n_group = size(aero_state%awa%weight, 1)
1549  n_class = size(aero_state%awa%weight, 2)
1550 
1551  ! if we have less than half the maximum number of particles then
1552  ! double until we fill up the array
1553  if (allow_doubling) then
1554  do i_group = 1,n_group
1555  do i_class = 1,n_class
1556  global_n_part &
1557  = aero_state_total_particles_all_procs(aero_state, i_group, &
1558  i_class)
1559  do while ((real(global_n_part, kind=dp) &
1560  < aero_state%n_part_ideal(i_group, i_class) / 2d0) &
1561  .and. (global_n_part > 0))
1562  if (initial_state_warning) then
1563  call warn_msg(716882783, "doubling particles in initial " &
1564  // "condition")
1565  end if
1566  call aero_state_double(aero_state, aero_data, i_group, i_class)
1567  global_n_part &
1568  = aero_state_total_particles_all_procs(aero_state, &
1569  i_group, i_class)
1570  end do
1571  end do
1572  end do
1573  end if
1574 
1575  ! same for halving if we have too many particles
1576  if (allow_halving) then
1577  do i_group = 1,n_group
1578  do i_class = 1,n_class
1579  do while (real(aero_state_total_particles_all_procs(aero_state, & i_group, i_class), kind=dp) &
1580  > aero_state%n_part_ideal(i_group, i_class) * 2d0)
1581  if (initial_state_warning) then
1582  call warn_msg(661936373, &
1583  "halving particles in initial condition")
1584  end if
1585  call aero_state_halve(aero_state, i_group, i_class)
1586  end do
1587  end do
1588  end do
1589  end if
1590 
1591  end subroutine aero_state_rebalance
1592 
1593 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1594 
1595  !> Scale the weighting of the given group/class by the given ratio,
1596  !> altering particle number as necessary to preserve the number
1597  !> concentration.
1598  subroutine aero_state_scale_weight(aero_state, aero_data, i_group, &
1599  i_class, weight_ratio, allow_doubling, allow_halving)
1600 
1601  !> Aerosol state.
1602  type(aero_state_t), intent(inout) :: aero_state
1603  !> Aerosol data.
1604  type(aero_data_t), intent(in) :: aero_data
1605  !> Weight group number.
1606  integer, intent(in) :: i_group
1607  !> Weight class number.
1608  integer, intent(in) :: i_class
1609  !> Ratio of <tt>new_weight / old_weight</tt>.
1610  real(kind=dp), intent(in) :: weight_ratio
1611  !> Whether to allow doubling of the population.
1612  logical, intent(in) :: allow_doubling
1613  !> Whether to allow halving of the population.
1614  logical, intent(in) :: allow_halving
1615 
1616  real(kind=dp) :: ratio
1617  integer :: i_part, i_remove, n_remove, i_entry, n_part
1618  type(aero_info_t) :: aero_info
1619 
1620  ! We could use the ratio < 1 case unconditionally, but that would
1621  ! have higher variance for the ratio > 1 case than the current
1622  ! scheme.
1623 
1624  call aero_state_sort(aero_state, aero_data)
1625  n_part = integer_varray_n_entry( &
1626  aero_state%aero_sorted%group_class%inverse(i_group, i_class))
1627 
1628  if ((weight_ratio > 1d0) .and. (allow_halving .or. (n_part == 0))) then
1629  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), &
1630  weight_ratio)
1631  n_remove = prob_round(real(n_part, kind=dp) &
1632  * (1d0 - 1d0 / weight_ratio))
1633  do i_remove = 1,n_remove
1635  aero_state%aero_sorted%group_class%inverse(i_group, i_class)))
1636  i_part = aero_state%aero_sorted%group_class%inverse(i_group, &
1637  i_class)%entry(i_entry)
1638  aero_info%id = aero_state%apa%particle(i_part)%id
1639  aero_info%action = aero_info_halved
1640  call aero_state_remove_particle(aero_state, i_part, .true., &
1641  aero_info)
1642  end do
1643  elseif ((weight_ratio < 1d0) &
1644  .and. (allow_doubling .or. (n_part == 0))) then
1645  call aero_weight_scale(aero_state%awa%weight(i_group, i_class), &
1646  weight_ratio)
1647  do i_entry = n_part,1,-1
1648  i_part = aero_state%aero_sorted%group_class%inverse(i_group, &
1649  i_class)%entry(i_entry)
1650  call aero_state_dup_particle(aero_state, aero_data, i_part, &
1651  1d0 / weight_ratio)
1652  end do
1653  end if
1654 
1655  end subroutine aero_state_scale_weight
1656 
1657 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1658 
1659  !> Mix the aero_states between all processes. Currently uses a
1660  !> simple all-to-all diffusion.
1661  subroutine aero_state_mix(aero_state, del_t, mix_timescale, &
1662  aero_data, specify_prob_transfer)
1663 
1664  !> Aerosol state.
1665  type(aero_state_t), intent(inout) :: aero_state
1666  !> Timestep (s).
1667  real(kind=dp), intent(in) :: del_t
1668  !> Mixing timescale (s).
1669  real(kind=dp), intent(in) :: mix_timescale
1670  !> Aero data values.
1671  type(aero_data_t), intent(in) :: aero_data
1672  !> Transfer probability of each particle (0 means no mixing, 1
1673  !> means total mixing).
1674  real(kind=dp), optional, intent(in) :: specify_prob_transfer
1675 
1676 #ifdef PMC_USE_MPI
1677  integer :: rank, n_proc, i_proc, ierr
1678  integer :: buffer_size, buffer_size_check
1679  character, allocatable :: buffer(:)
1680  type(aero_state_t), allocatable :: aero_state_sends(:)
1681  type(aero_state_t), allocatable :: aero_state_recvs(:)
1682  real(kind=dp) :: prob_transfer, prob_not_transferred
1683  real(kind=dp) :: prob_transfer_given_not_transferred
1684 
1685  ! process information
1686  rank = pmc_mpi_rank()
1687  n_proc = pmc_mpi_size()
1688  if (n_proc == 1) then
1689  ! buffer allocation below fails if n_proc == 1
1690  ! so bail out early (nothing to mix anyway)
1691  return
1692  end if
1693 
1694  ! allocate aero_state arrays
1695  allocate(aero_state_sends(n_proc))
1696  allocate(aero_state_recvs(n_proc))
1697 
1698  ! compute the transfer probability
1699  if (present(specify_prob_transfer)) then
1700  prob_transfer = specify_prob_transfer / real(n_proc, kind=dp)
1701  else
1702  prob_transfer = (1d0 - exp(- del_t / mix_timescale)) &
1703  / real(n_proc, kind=dp)
1704  end if
1705 
1706  ! extract particles to send
1707  prob_not_transferred = 1d0
1708  do i_proc = 0,(n_proc - 1)
1709  if (i_proc /= rank) then
1710  ! because we are doing sequential sampling from the aero_state
1711  ! we need to scale up the later transfer probabilities, because
1712  ! the later particles are being transferred conditioned on the
1713  ! fact that they did not transfer already
1714  prob_transfer_given_not_transferred = prob_transfer &
1715  / prob_not_transferred
1716  call aero_state_sample(aero_state, &
1717  aero_state_sends(i_proc + 1), aero_data, &
1718  prob_transfer_given_not_transferred, aero_info_none)
1719  prob_not_transferred = prob_not_transferred - prob_transfer
1720  end if
1721  end do
1722 
1723  ! exchange the particles
1724  call aero_state_mpi_alltoall(aero_state_sends, aero_state_recvs)
1725 
1726  ! process the received particles
1727  do i_proc = 0,(n_proc - 1)
1728  if (i_proc /= rank) then
1729  call aero_state_add(aero_state, aero_state_recvs(i_proc + 1), &
1730  aero_data)
1731  end if
1732  end do
1733 
1734  ! cleanup
1735  deallocate(aero_state_sends)
1736  deallocate(aero_state_recvs)
1737 #endif
1738 
1739  end subroutine aero_state_mix
1740 
1741 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1742 
1743  !> Do an MPI all-to-all transfer of aerosol states.
1744  !!
1745  !! States without particles are not sent.
1746  subroutine aero_state_mpi_alltoall(send, recv)
1747 
1748  !> Array of aero_states to send (one per process).
1749  type(aero_state_t), intent(in) :: send(:)
1750  !> Array of aero_states to receives (one per process).
1751  type(aero_state_t), intent(inout) :: recv(size(send))
1752 
1753 #ifdef PMC_USE_MPI
1754  character, allocatable :: sendbuf(:), recvbuf(:)
1755  integer :: sendcounts(size(send)), sdispls(size(send))
1756  integer :: recvcounts(size(send)), rdispls(size(send))
1757  integer :: i_proc, position, old_position, max_sendbuf_size, ierr
1758 
1759  call assert(978709842, size(send) == pmc_mpi_size())
1760 
1761  max_sendbuf_size = 0
1762  do i_proc = 1,pmc_mpi_size()
1763  if (aero_state_total_particles(send(i_proc)) > 0) then
1764  max_sendbuf_size = max_sendbuf_size &
1765  + pmc_mpi_pack_size_aero_state(send(i_proc))
1766  end if
1767  end do
1768 
1769  allocate(sendbuf(max_sendbuf_size))
1770 
1771  position = 0
1772  do i_proc = 1,pmc_mpi_size()
1773  old_position = position
1774  if (aero_state_total_particles(send(i_proc)) > 0) then
1775  call pmc_mpi_pack_aero_state(sendbuf, position, send(i_proc))
1776  end if
1777  sendcounts(i_proc) = position - old_position
1778  end do
1779  call assert(393267406, position <= max_sendbuf_size)
1780 
1781  call pmc_mpi_alltoall_integer(sendcounts, recvcounts)
1782  allocate(recvbuf(sum(recvcounts)))
1783 
1784  sdispls(1) = 0
1785  rdispls(1) = 0
1786  do i_proc = 2,pmc_mpi_size()
1787  sdispls(i_proc) = sdispls(i_proc - 1) + sendcounts(i_proc - 1)
1788  rdispls(i_proc) = rdispls(i_proc - 1) + recvcounts(i_proc - 1)
1789  end do
1790 
1791  call mpi_alltoallv(sendbuf, sendcounts, sdispls, mpi_character, recvbuf, &
1792  recvcounts, rdispls, mpi_character, mpi_comm_world, ierr)
1793  call pmc_mpi_check_ierr(ierr)
1794 
1795  position = 0
1796  do i_proc = 1,pmc_mpi_size()
1797  call assert(189739257, position == rdispls(i_proc))
1798  if (recvcounts(i_proc) > 0) then
1799  call pmc_mpi_unpack_aero_state(recvbuf, position, recv(i_proc))
1800  end if
1801  end do
1802 
1803  deallocate(sendbuf)
1804  deallocate(recvbuf)
1805 #endif
1806 
1807  end subroutine aero_state_mpi_alltoall
1808 
1809 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1810 
1811  !> Set each aerosol particle to have its original total volume, but
1812  !> species volume ratios given by the total species volume ratio
1813  !> within each bin. This preserves the (weighted) total species
1814  !> volume per bin as well as per-particle total volumes.
1815  subroutine aero_state_bin_average_comp(aero_state, bin_grid, aero_data)
1816 
1817  !> Aerosol state to average.
1818  type(aero_state_t), intent(inout) :: aero_state
1819  !> Bin grid to average within.
1820  type(bin_grid_t), intent(in) :: bin_grid
1821  !> Aerosol data.
1822  type(aero_data_t), intent(in) :: aero_data
1823 
1824  real(kind=dp) :: species_volume_conc(aero_data_n_spec(aero_data))
1825  real(kind=dp) :: total_volume_conc, particle_volume, num_conc
1826  integer :: i_bin, i_class, i_entry, i_part, i_spec
1827 
1828  call aero_state_sort(aero_state, aero_data, bin_grid)
1829 
1830  do i_bin = 1,bin_grid_size(bin_grid)
1831  species_volume_conc = 0d0
1832  total_volume_conc = 0d0
1833  do i_class = 1,size(aero_state%awa%weight, 2)
1834  do i_entry = 1,integer_varray_n_entry( &
1835  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1836  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1837  i_class)%entry(i_entry)
1838  num_conc = aero_weight_array_num_conc(aero_state%awa, &
1839  aero_state%apa%particle(i_part), aero_data)
1840  particle_volume = aero_particle_volume( &
1841  aero_state%apa%particle(i_part))
1842  species_volume_conc = species_volume_conc &
1843  + num_conc * aero_state%apa%particle(i_part)%vol
1844  total_volume_conc = total_volume_conc + num_conc * particle_volume
1845  end do
1846  end do
1847  do i_class = 1,size(aero_state%awa%weight, 2)
1848  do i_entry = 1,integer_varray_n_entry( &
1849  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1850  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1851  i_class)%entry(i_entry)
1852  particle_volume = aero_particle_volume( &
1853  aero_state%apa%particle(i_part))
1854  aero_state%apa%particle(i_part)%vol &
1855  = particle_volume * species_volume_conc / total_volume_conc
1856  end do
1857  end do
1858  end do
1859 
1860  end subroutine aero_state_bin_average_comp
1861 
1862 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1863 
1864  !> Set each aerosol particle to have its original species ratios,
1865  !> but total volume given by the average volume of all particles
1866  !> within each bin.
1867  !!
1868  !! This does not preserve the total species volume
1869  !! per bin. If the \c bin_center parameter is \c .true. then the
1870  !! particles in each bin are set to have the bin center volume,
1871  !! rather than the average volume of the particles in that bin.
1872  !!
1873  !! If the weighting function is not constant (AERO_WEIGHT_TYPE_NONE)
1874  !! then the averaging can be performed in either a number-preserving
1875  !! way or in a volume-preserving way. The volume-preserving way does
1876  !! not preserve species volume ratios in gernal, but will do so if
1877  !! the particle population has already been composition-averaged.
1878  subroutine aero_state_bin_average_size(aero_state, bin_grid, aero_data, &
1879  bin_center, preserve_number)
1880 
1881  !> Aerosol state to average.
1882  type(aero_state_t), intent(inout) :: aero_state
1883  !> Bin grid to average within.
1884  type(bin_grid_t), intent(in) :: bin_grid
1885  !> Aerosol data.
1886  type(aero_data_t), intent(in) :: aero_data
1887  !> Whether to assign the bin center volume (rather than the average
1888  !> volume).
1889  logical, intent(in) :: bin_center
1890  !> Whether to use the number-preserving scheme (otherwise will use
1891  !> the volume-preserving scheme). This parameter has no effect if
1892  !> \c bin_center is \c .true.
1893  logical, intent(in) :: preserve_number
1894 
1895  real(kind=dp) :: total_volume_conc, particle_volume
1896  real(kind=dp) :: new_particle_volume, num_conc, total_num_conc
1897  real(kind=dp) :: lower_volume, upper_volume, center_volume
1898  real(kind=dp) :: lower_function, upper_function, center_function
1899  integer :: i_bin, i_class, i_entry, i_part, i_bisect, n_part
1900  logical :: monotone_increasing, monotone_decreasing
1901 
1902  call aero_state_sort(aero_state, aero_data, bin_grid)
1903 
1904  do i_bin = 1,bin_grid_size(bin_grid)
1905  do i_class = 1,size(aero_state%awa%weight, 2)
1906  if (integer_varray_n_entry( &
1907  aero_state%aero_sorted%size_class%inverse(i_bin, i_class)) &
1908  == 0) then
1909  cycle
1910  end if
1911 
1912  n_part = integer_varray_n_entry( &
1913  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1914  total_num_conc = 0d0
1915  total_volume_conc = 0d0
1916  do i_entry = 1,integer_varray_n_entry( &
1917  aero_state%aero_sorted%size_class%inverse(i_bin, i_class))
1918  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1919  i_class)%entry(i_entry)
1920  num_conc = aero_weight_array_num_conc(aero_state%awa, &
1921  aero_state%apa%particle(i_part), aero_data)
1922  total_num_conc = total_num_conc + num_conc
1923  particle_volume = aero_particle_volume( &
1924  aero_state%apa%particle(i_part))
1925  total_volume_conc = total_volume_conc &
1926  + num_conc * particle_volume
1927  end do
1928 
1929  ! determine the new_particle_volume for all particles in this bin
1930  if (bin_center) then
1931  new_particle_volume = aero_data_rad2vol(aero_data, &
1932  bin_grid%centers(i_bin))
1933  elseif (aero_weight_array_check_flat(aero_state%awa)) then
1934  num_conc & ! any radius will have the same num_conc
1935  = aero_weight_array_num_conc_at_radius(aero_state%awa, &
1936  i_class, 1d0)
1937  new_particle_volume = total_volume_conc / num_conc &
1938  / real(integer_varray_n_entry( & aero_state%aero_sorted%size_class%inverse(i_bin, i_class)), &
1939  kind=dp)
1940  elseif (preserve_number) then
1941  ! number-preserving scheme: Solve the implicit equation:
1942  ! n_part * W(new_vol) = total_num_conc
1943  !
1944  ! We assume that the weighting function is strictly monotone
1945  ! so this equation has a unique solution and the solution
1946  ! lies between the min and max particle volumes. We use
1947  ! bisection as this doesn't really need to be fast, just
1948  ! robust.
1949 
1950  call aero_weight_array_check_monotonicity(aero_state%awa, &
1951  monotone_increasing, monotone_decreasing)
1952  call assert_msg(214077200, &
1953  monotone_increasing .or. monotone_decreasing, &
1954  "monotone weight function required for averaging")
1955 
1956  ! initialize to min and max particle volumes
1957  do i_entry = 1,n_part
1958  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
1959  i_class)%entry(i_entry)
1960  particle_volume = aero_particle_volume( &
1961  aero_state%apa%particle(i_part))
1962  if (i_part == 1) then
1963  lower_volume = particle_volume
1964  upper_volume = particle_volume
1965  end if
1966  lower_volume = min(lower_volume, particle_volume)
1967  upper_volume = max(upper_volume, particle_volume)
1968  end do
1969 
1970  lower_function = real(n_part, kind=dp) &
1971  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1972  i_class, aero_data_vol2rad(aero_data, lower_volume)) &
1973  - total_num_conc
1974  upper_function = real(n_part, kind=dp) &
1975  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1976  i_class, aero_data_vol2rad(aero_data, upper_volume)) &
1977  - total_num_conc
1978 
1979  ! do 50 rounds of bisection (2^50 = 10^15)
1980  do i_bisect = 1,50
1981  center_volume = (lower_volume + upper_volume) / 2d0
1982  center_function = real(n_part, kind=dp) &
1983  * aero_weight_array_num_conc_at_radius(aero_state%awa, &
1984  i_class, aero_data_vol2rad(aero_data, center_volume)) &
1985  - total_num_conc
1986  if ((lower_function > 0d0 .and. center_function > 0d0) &
1987  .or. (lower_function < 0d0 .and. center_function < 0d0)) &
1988  then
1989  lower_volume = center_volume
1990  lower_function = center_function
1991  else
1992  upper_volume = center_volume
1993  upper_function = center_function
1994  end if
1995  end do
1996 
1997  new_particle_volume = center_volume
1998  else
1999  ! volume-preserving scheme: Solve the implicit equation:
2000  ! n_part * W(new_vol) * new_vol = total_volume_conc
2001  !
2002  ! We assume that the weighting function is strictly monotone
2003  ! so this equation has a unique solution and the solution
2004  ! lies between the min and max particle volumes. We use
2005  ! bisection as this doesn't really need to be fast, just
2006  ! robust.
2007 
2008  call aero_weight_array_check_monotonicity(aero_state%awa, &
2009  monotone_increasing, monotone_decreasing)
2010  call assert_msg(483078128, &
2011  monotone_increasing .or. monotone_decreasing, &
2012  "monotone weight function required for averaging")
2013 
2014  ! initialize to min and max particle volumes
2015  do i_entry = 1,n_part
2016  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
2017  i_class)%entry(i_entry)
2018  particle_volume = aero_particle_volume( &
2019  aero_state%apa%particle(i_part))
2020  if (i_part == 1) then
2021  lower_volume = particle_volume
2022  upper_volume = particle_volume
2023  end if
2024  lower_volume = min(lower_volume, particle_volume)
2025  upper_volume = max(upper_volume, particle_volume)
2026  end do
2027 
2028  lower_function = real(n_part, kind=dp) &
2029  * aero_weight_array_num_conc_at_radius( &
2030  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2031  lower_volume)) * lower_volume - total_volume_conc
2032  upper_function = real(n_part, kind=dp) &
2033  * aero_weight_array_num_conc_at_radius( &
2034  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2035  upper_volume)) * upper_volume - total_volume_conc
2036 
2037  ! do 50 rounds of bisection (2^50 = 10^15)
2038  do i_bisect = 1,50
2039  center_volume = (lower_volume + upper_volume) / 2d0
2040  center_function = real(n_part, kind=dp) &
2041  * aero_weight_array_num_conc_at_radius( &
2042  aero_state%awa, i_class, aero_data_vol2rad(aero_data, &
2043  center_volume)) * center_volume - total_volume_conc
2044  if ((lower_function > 0d0 .and. center_function > 0d0) &
2045  .or. (lower_function < 0d0 .and. center_function < 0d0)) &
2046  then
2047  lower_volume = center_volume
2048  lower_function = center_function
2049  else
2050  upper_volume = center_volume
2051  upper_function = center_function
2052  end if
2053  end do
2054 
2055  new_particle_volume = center_volume
2056  end if
2057 
2058  do i_entry = 1,n_part
2059  i_part = aero_state%aero_sorted%size_class%inverse(i_bin, &
2060  i_class)%entry(i_entry)
2061  particle_volume = aero_particle_volume( &
2062  aero_state%apa%particle(i_part))
2063  aero_state%apa%particle(i_part)%vol &
2064  = aero_state%apa%particle(i_part)%vol &
2065  / particle_volume * new_particle_volume
2066  end do
2067  end do
2068  end do
2069 
2070  end subroutine aero_state_bin_average_size
2071 
2072 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2073 
2074  !> Make all particles dry (water set to zero).
2075  subroutine aero_state_make_dry(aero_state, aero_data)
2076 
2077  !> Aerosol state to make dry.
2078  type(aero_state_t), intent(inout) :: aero_state
2079  !> Aerosol data.
2080  type(aero_data_t), intent(in) :: aero_data
2081 
2082  integer :: i_part
2083  real(kind=dp) :: reweight_num_conc(aero_state%apa%n_part)
2084 
2085  ! We're modifying particle diameters, so bin sorting is now invalid
2086  aero_state%valid_sort = .false.
2087 
2088  call aero_state_num_conc_for_reweight(aero_state, aero_data, &
2089  reweight_num_conc)
2090  if (aero_data%i_water > 0) then
2091  do i_part = 1,aero_state_n_part(aero_state)
2092  aero_state%apa%particle(i_part)%vol(aero_data%i_water) = 0d0
2093  end do
2094  aero_state%valid_sort = .false.
2095  end if
2096  ! adjust particles to account for weight changes
2097  call aero_state_reweight(aero_state, aero_data, reweight_num_conc)
2098 
2099  end subroutine aero_state_make_dry
2100 
2101 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2102 
2103  !> Determines the number of bytes required to pack the given value.
2104  integer function pmc_mpi_pack_size_aero_state(val)
2105 
2106  !> Value to pack.
2107  type(aero_state_t), intent(in) :: val
2108 
2109  integer :: total_size, i_group
2110 
2111  total_size = 0
2112  total_size = total_size + pmc_mpi_pack_size_apa(val%apa)
2113  total_size = total_size + pmc_mpi_pack_size_aero_weight_array(val%awa)
2114  total_size = total_size + pmc_mpi_pack_size_real_array_2d(val%n_part_ideal)
2115  total_size = total_size + pmc_mpi_pack_size_aia(val%aero_info_array)
2116  pmc_mpi_pack_size_aero_state = total_size
2117 
2118  end function pmc_mpi_pack_size_aero_state
2119 
2120 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2121 
2122  !> Packs the given value into the buffer, advancing position.
2123  subroutine pmc_mpi_pack_aero_state(buffer, position, val)
2124 
2125  !> Memory buffer.
2126  character, intent(inout) :: buffer(:)
2127  !> Current buffer position.
2128  integer, intent(inout) :: position
2129  !> Value to pack.
2130  type(aero_state_t), intent(in) :: val
2131 
2132 #ifdef PMC_USE_MPI
2133  integer :: prev_position, i_group
2134 
2135  prev_position = position
2136  call pmc_mpi_pack_aero_particle_array(buffer, position, val%apa)
2137  call pmc_mpi_pack_aero_weight_array(buffer,position,val%awa)
2138  call pmc_mpi_pack_real_array_2d(buffer, position, val%n_part_ideal)
2139  call pmc_mpi_pack_aero_info_array(buffer, position, val%aero_info_array)
2140  call assert(850997402, &
2141  position - prev_position <= pmc_mpi_pack_size_aero_state(val))
2142 #endif
2143 
2144  end subroutine pmc_mpi_pack_aero_state
2145 
2146 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2147 
2148  !> Unpacks the given value from the buffer, advancing position.
2149  subroutine pmc_mpi_unpack_aero_state(buffer, position, val)
2150 
2151  !> Memory buffer.
2152  character, intent(inout) :: buffer(:)
2153  !> Current buffer position.
2154  integer, intent(inout) :: position
2155  !> Value to pack.
2156  type(aero_state_t), intent(inout) :: val
2157 
2158 #ifdef PMC_USE_MPI
2159  integer :: prev_position, i_group, n_group
2160 
2161  val%valid_sort = .false.
2162  prev_position = position
2163  call pmc_mpi_unpack_aero_particle_array(buffer, position, val%apa)
2164  call pmc_mpi_unpack_aero_weight_array(buffer,position,val%awa)
2165  call pmc_mpi_unpack_real_array_2d(buffer, position, val%n_part_ideal)
2166  call pmc_mpi_unpack_aero_info_array(buffer, position, val%aero_info_array)
2167  call assert(132104747, &
2168  position - prev_position <= pmc_mpi_pack_size_aero_state(val))
2169 #endif
2170 
2171  end subroutine pmc_mpi_unpack_aero_state
2172 
2173 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2174 
2175  !> Gathers data from all processes into one aero_state on process 0.
2176  subroutine aero_state_mpi_gather(aero_state, aero_state_total, aero_data)
2177 
2178  !> Local aero_state.
2179  type(aero_state_t), intent(in) :: aero_state
2180  !> Centralized aero_state (only on process 0).
2181  type(aero_state_t), intent(inout) :: aero_state_total
2182  !> Aero data values.
2183  type(aero_data_t), intent(in) :: aero_data
2184 
2185 #ifdef PMC_USE_MPI
2186  type(aero_state_t) :: aero_state_transfer
2187  integer :: n_proc, ierr, status(mpi_status_size)
2188  integer :: buffer_size, max_buffer_size, i_proc, position
2189  character, allocatable :: buffer(:)
2190 #endif
2191 
2192  if (pmc_mpi_rank() == 0) then
2193  aero_state_total = aero_state
2194  end if
2195 
2196 #ifdef PMC_USE_MPI
2197 
2198  if (pmc_mpi_rank() /= 0) then
2199  ! send data from remote processes
2200  max_buffer_size = 0
2201  max_buffer_size = max_buffer_size &
2202  + pmc_mpi_pack_size_aero_state(aero_state)
2203  allocate(buffer(max_buffer_size))
2204  position = 0
2205  call pmc_mpi_pack_aero_state(buffer, position, aero_state)
2206  call assert(542772170, position <= max_buffer_size)
2207  buffer_size = position
2208  call mpi_send(buffer, buffer_size, mpi_character, 0, &
2209  aero_state_tag_gather, mpi_comm_world, ierr)
2210  call pmc_mpi_check_ierr(ierr)
2211  deallocate(buffer)
2212  else
2213  ! root process receives data from each remote proc
2214  n_proc = pmc_mpi_size()
2215  do i_proc = 1,(n_proc - 1)
2216  ! determine buffer size at root process
2217  call mpi_probe(i_proc, aero_state_tag_gather, mpi_comm_world, &
2218  status, ierr)
2219  call pmc_mpi_check_ierr(ierr)
2220  call mpi_get_count(status, mpi_character, buffer_size, ierr)
2221  call pmc_mpi_check_ierr(ierr)
2222 
2223  ! get buffer at root process
2224  allocate(buffer(buffer_size))
2225  call mpi_recv(buffer, buffer_size, mpi_character, i_proc, &
2226  aero_state_tag_gather, mpi_comm_world, status, ierr)
2227 
2228  ! unpack it
2229  position = 0
2230  call pmc_mpi_unpack_aero_state(buffer, position, &
2231  aero_state_transfer)
2232  call assert(518174881, position == buffer_size)
2233  deallocate(buffer)
2234 
2235  call aero_state_add(aero_state_total, aero_state_transfer, &
2236  aero_data)
2237  end do
2238  end if
2239 
2240 #endif
2241 
2242  end subroutine aero_state_mpi_gather
2243 
2244 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2245 
2246  !> Write the aero particle dimension to the given NetCDF file if it
2247  !> is not already present and in any case return the associated
2248  !> dimid.
2249  subroutine aero_state_netcdf_dim_aero_particle(aero_state, ncid, &
2250  dimid_aero_particle)
2251 
2252  !> aero_state structure.
2253  type(aero_state_t), intent(in) :: aero_state
2254  !> NetCDF file ID, in data mode.
2255  integer, intent(in) :: ncid
2256  !> Dimid of the aero particle dimension.
2257  integer, intent(out) :: dimid_aero_particle
2258 
2259  integer :: status, i_part
2260  integer :: varid_aero_particle
2261  integer :: aero_particle_centers(aero_state_n_part(aero_state))
2262 
2263  ! try to get the dimension ID
2264  status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
2265  if (status == nf90_noerr) return
2266  if (status /= nf90_ebaddim) call pmc_nc_check(status)
2267 
2268  ! dimension not defined, so define now define it
2269  call pmc_nc_check(nf90_redef(ncid))
2270 
2271  call pmc_nc_check(nf90_def_dim(ncid, "aero_particle", &
2272  aero_state_n_part(aero_state), dimid_aero_particle))
2273 
2274  call pmc_nc_check(nf90_enddef(ncid))
2275 
2276  do i_part = 1,aero_state_n_part(aero_state)
2277  aero_particle_centers(i_part) = i_part
2278  end do
2279  call pmc_nc_write_integer_1d(ncid, aero_particle_centers, &
2280  "aero_particle", (/ dimid_aero_particle /), &
2281  description="dummy dimension variable (no useful value)")
2282 
2284 
2285 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2286 
2287  !> Write the aero removed dimension to the given NetCDF file if it
2288  !> is not already present and in any case return the associated
2289  !> dimid.
2290  subroutine aero_state_netcdf_dim_aero_removed(aero_state, ncid, &
2291  dimid_aero_removed)
2292 
2293  !> aero_state structure.
2294  type(aero_state_t), intent(in) :: aero_state
2295  !> NetCDF file ID, in data mode.
2296  integer, intent(in) :: ncid
2297  !> Dimid of the aero removed dimension.
2298  integer, intent(out) :: dimid_aero_removed
2299 
2300  integer :: status, i_remove, dim_size
2301  integer :: varid_aero_removed
2302  integer :: aero_removed_centers( &
2303  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2304 
2305  ! try to get the dimension ID
2306  status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed)
2307  if (status == nf90_noerr) return
2308  if (status /= nf90_ebaddim) call pmc_nc_check(status)
2309 
2310  ! dimension not defined, so define now define it
2311  call pmc_nc_check(nf90_redef(ncid))
2312 
2313  dim_size = max(1, aero_info_array_n_item(aero_state%aero_info_array))
2314  call pmc_nc_check(nf90_def_dim(ncid, "aero_removed", &
2315  dim_size, dimid_aero_removed))
2316 
2317  call pmc_nc_check(nf90_enddef(ncid))
2318 
2319  do i_remove = 1,dim_size
2320  aero_removed_centers(i_remove) = i_remove
2321  end do
2322  call pmc_nc_write_integer_1d(ncid, aero_removed_centers, &
2323  "aero_removed", (/ dimid_aero_removed /), &
2324  description="dummy dimension variable (no useful value)")
2325 
2326  end subroutine aero_state_netcdf_dim_aero_removed
2327 
2328 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2329 
2330  !> Write full state.
2331  subroutine aero_state_output_netcdf(aero_state, ncid, aero_data, &
2332  record_removals, record_optical)
2333 
2334  !> aero_state to write.
2335  type(aero_state_t), intent(in) :: aero_state
2336  !> NetCDF file ID, in data mode.
2337  integer, intent(in) :: ncid
2338  !> aero_data structure.
2339  type(aero_data_t), intent(in) :: aero_data
2340  !> Whether to output particle removal info.
2341  logical, intent(in) :: record_removals
2342  !> Whether to output aerosol optical properties.
2343  logical, intent(in) :: record_optical
2344 
2345  integer :: dimid_aero_particle, dimid_aero_species, dimid_aero_source
2346  integer :: dimid_aero_removed
2347  integer :: i_part, i_remove
2348  real(kind=dp) :: aero_particle_mass(aero_state_n_part(aero_state), &
2349  aero_data_n_spec(aero_data))
2350  integer :: aero_n_orig_part(aero_state_n_part(aero_state), &
2351  aero_data_n_source(aero_data))
2352  integer :: aero_particle_weight_group(aero_state_n_part(aero_state))
2353  integer :: aero_particle_weight_class(aero_state_n_part(aero_state))
2354  real(kind=dp) :: aero_absorb_cross_sect(aero_state_n_part(aero_state))
2355  real(kind=dp) :: aero_scatter_cross_sect(aero_state_n_part(aero_state))
2356  real(kind=dp) :: aero_asymmetry(aero_state_n_part(aero_state))
2357  real(kind=dp) :: aero_refract_shell_real(aero_state_n_part(aero_state))
2358  real(kind=dp) :: aero_refract_shell_imag(aero_state_n_part(aero_state))
2359  real(kind=dp) :: aero_refract_core_real(aero_state_n_part(aero_state))
2360  real(kind=dp) :: aero_refract_core_imag(aero_state_n_part(aero_state))
2361  real(kind=dp) :: aero_core_vol(aero_state_n_part(aero_state))
2362  integer :: aero_water_hyst_leg(aero_state_n_part(aero_state))
2363  real(kind=dp) :: aero_num_conc(aero_state_n_part(aero_state))
2364  integer :: aero_id(aero_state_n_part(aero_state))
2365  real(kind=dp) :: aero_least_create_time(aero_state_n_part(aero_state))
2366  real(kind=dp) :: aero_greatest_create_time(aero_state_n_part(aero_state))
2367  integer :: aero_removed_id( &
2368  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2369  integer :: aero_removed_action( &
2370  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2371  integer :: aero_removed_other_id( &
2372  max(1, aero_info_array_n_item(aero_state%aero_info_array)))
2373 
2374  !> \page output_format_aero_state Output File Format: Aerosol Particle State
2375  !!
2376  !! The aerosol state consists of a set of individual aerosol
2377  !! particles, each with its own individual properties. The
2378  !! properties of all particles are stored in arrays, one per
2379  !! property. For example, <tt>aero_absorb_cross_sect(i)</tt> gives
2380  !! the absorption cross section of particle number \c i, while
2381  !! <tt>aero_particle_mass(i,s)</tt> gives the mass of species \c s
2382  !! in particle \c i. The aerosol species are described in \ref
2383  !! output_format_aero_data.
2384  !!
2385  !! Each aerosol particle \c i represents a number concentration
2386  !! given by <tt>aero_num_conc(i)</tt>. Multiplying a per-particle
2387  !! quantity by the respective number concentration gives the
2388  !! concentration of that quantity contributed by the particle. For
2389  !! example, summing <tt>aero_particle_mass(i,s) *
2390  !! aero_num_conc(i)</tt> over all \c i gives the total mass
2391  !! concentration of species \c s in kg/m^3. Similarly, summing
2392  !! <tt>aero_absorb_cross_sect(i) * aero_num_conc(i)</tt> over all
2393  !! \c i will give the concentration of scattering cross section in
2394  !! m^2/m^3.
2395  !!
2396  !! FIXME: the aero_weight is also output
2397  !!
2398  !! The aerosol state uses the \c aero_species NetCDF dimension as
2399  !! specified in the \ref output_format_aero_data section, as well
2400  !! as the NetCDF dimension:
2401  !! - \b aero_particle: number of aerosol particles
2402  !!
2403  !! The aerosol state NetCDF variables are:
2404  !! - \b aero_particle (dim \c aero_particle): dummy dimension variable
2405  !! (no useful value)
2406  !! - \b aero_particle_mass (unit kg,
2407  !! dim <tt>aero_particle x aero_species</tt>): constituent masses of
2408  !! each aerosol particle - <tt>aero_particle_mass(i,s)</tt> gives the
2409  !! mass of species \c s in particle \c i
2410  !! - \b aero_n_orig_part (dim <tt>aero_particle x
2411  !! aero_source</tt>): number of original particles from each
2412  !! source that formed each aerosol particle -
2413  !! <tt>aero_n_orig_part(i,s)</tt> is the number of particles
2414  !! from source \c s that contributed to particle \c i - when
2415  !! particle \c i first enters the simulation (by emissions,
2416  !! dilution, etc.) it has <tt>aero_n_orig_part(i,s) = 1</tt>
2417  !! for the source number \c s it came from (otherwise zero)
2418  !! and when two particles coagulate, their values of \c
2419  !! aero_n_orig_part are added (the number of coagulation
2420  !! events that formed each particle is thus
2421  !! <tt>sum(aero_n_orig_part(i,:)) - 1</tt>)
2422  !! - \b aero_particle_weight_group (dim <tt>aero_particle</tt>):
2423  !! weight group number (1 to <tt>aero_weight_group</tt>) of
2424  !! each aerosol particle
2425  !! - \b aero_particle_weight_class (dim <tt>aero_particle</tt>):
2426  !! weight class number (1 to <tt>aero_weight_class</tt>) of each
2427  !! aerosol particle
2428  !! - \b aero_absorb_cross_sect (unit m^2, dim \c aero_particle):
2429  !! optical absorption cross sections of each aerosol particle
2430  !! - \b aero_scatter_cross_sect (unit m^2, dim \c aero_particle):
2431  !! optical scattering cross sections of each aerosol particle
2432  !! - \b aero_asymmetry (dimensionless, dim \c aero_particle): optical
2433  !! asymmetry parameters of each aerosol particle
2434  !! - \b aero_refract_shell_real (dimensionless, dim \c aero_particle):
2435  !! real part of the refractive indices of the shell of each
2436  !! aerosol particle
2437  !! - \b aero_refract_shell_imag (dimensionless, dim \c aero_particle):
2438  !! imaginary part of the refractive indices of the shell of each
2439  !! aerosol particle
2440  !! - \b aero_refract_core_real (dimensionless, dim \c aero_particle):
2441  !! real part of the refractive indices of the core of each
2442  !! aerosol particle
2443  !! - \b aero_refract_core_imag (dimensionless, dim \c aero_particle):
2444  !! imaginary part of the refractive indices of the core of each
2445  !! aerosol particle
2446  !! - \b aero_core_vol (unit m^3, dim \c aero_particle): volume of the
2447  !! optical cores of each aerosol particle
2448  !! - \b aero_water_hyst_leg (dim \c aero_particle): integers
2449  !! specifying which leg of the water hysteresis curve each
2450  !! particle is on, using the MOSAIC numbering convention
2451  !! - \b aero_num_conc (unit m^{-3}, dim \c aero_particle): number
2452  !! concentration associated with each particle
2453  !! - \b aero_id (dim \c aero_particle): unique ID number of each
2454  !! aerosol particle
2455  !! - \b aero_least_create_time (unit s, dim \c aero_particle): least
2456  !! (earliest) creation time of any original constituent particles
2457  !! that coagulated to form each particle, measured from the start
2458  !! of the simulation - a particle is said to be created when it
2459  !! first enters the simulation (by emissions, dilution, etc.)
2460  !! - \b aero_greatest_create_time (unit s, dim \c
2461  !! aero_particle): greatest (latest) creation time of any
2462  !! original constituent particles that coagulated to form each
2463  !! particle, measured from the start of the simulation - a
2464  !! particle is said to be created when it first enters the
2465  !! simulation (by emissions, dilution, etc.)
2466 
2467  call aero_weight_array_output_netcdf(aero_state%awa, ncid)
2468 
2469  call aero_data_netcdf_dim_aero_species(aero_data, ncid, &
2470  dimid_aero_species)
2471  call aero_data_netcdf_dim_aero_source(aero_data, ncid, &
2472  dimid_aero_source)
2473 
2474  if (aero_state_n_part(aero_state) > 0) then
2475  call aero_state_netcdf_dim_aero_particle(aero_state, ncid, &
2476  dimid_aero_particle)
2477 
2478  do i_part = 1,aero_state_n_part(aero_state)
2479  aero_particle_mass(i_part, :) &
2480  = aero_state%apa%particle(i_part)%vol * aero_data%density
2481  aero_n_orig_part(i_part, :) &
2482  = aero_state%apa%particle(i_part)%n_orig_part
2483  aero_particle_weight_group(i_part) &
2484  = aero_state%apa%particle(i_part)%weight_group
2485  aero_particle_weight_class(i_part) &
2486  = aero_state%apa%particle(i_part)%weight_class
2487  aero_water_hyst_leg(i_part) &
2488  = aero_state%apa%particle(i_part)%water_hyst_leg
2489  aero_num_conc(i_part) &
2490  = aero_state_particle_num_conc(aero_state, &
2491  aero_state%apa%particle(i_part), aero_data)
2492  aero_id(i_part) = aero_state%apa%particle(i_part)%id
2493  aero_least_create_time(i_part) &
2494  = aero_state%apa%particle(i_part)%least_create_time
2495  aero_greatest_create_time(i_part) &
2496  = aero_state%apa%particle(i_part)%greatest_create_time
2497  if (record_optical) then
2498  aero_absorb_cross_sect(i_part) &
2499  = aero_state%apa%particle(i_part)%absorb_cross_sect
2500  aero_scatter_cross_sect(i_part) &
2501  = aero_state%apa%particle(i_part)%scatter_cross_sect
2502  aero_asymmetry(i_part) = aero_state%apa%particle(i_part)%asymmetry
2503  aero_refract_shell_real(i_part) &
2504  = real(aero_state%apa%particle(i_part)%refract_shell)
2505  aero_refract_shell_imag(i_part) &
2506  = aimag(aero_state%apa%particle(i_part)%refract_shell)
2507  aero_refract_core_real(i_part) &
2508  = real(aero_state%apa%particle(i_part)%refract_core)
2509  aero_refract_core_imag(i_part) &
2510  = aimag(aero_state%apa%particle(i_part)%refract_core)
2511  aero_core_vol(i_part) = aero_state%apa%particle(i_part)%core_vol
2512  end if
2513  end do
2514  call pmc_nc_write_real_2d(ncid, aero_particle_mass, &
2515  "aero_particle_mass", (/ dimid_aero_particle, &
2516  dimid_aero_species /), unit="kg", &
2517  long_name="constituent masses of each aerosol particle")
2518  call pmc_nc_write_integer_2d(ncid, aero_n_orig_part, &
2519  "aero_n_orig_part", (/ dimid_aero_particle, &
2520  dimid_aero_source /), &
2521  long_name="number of original constituent particles from " &
2522  // "each source that coagulated to form each aerosol particle")
2523  call pmc_nc_write_integer_1d(ncid, aero_particle_weight_group, &
2524  "aero_particle_weight_group", (/ dimid_aero_particle /), &
2525  long_name="weight group number of each aerosol particle")
2526  call pmc_nc_write_integer_1d(ncid, aero_particle_weight_class, &
2527  "aero_particle_weight_class", (/ dimid_aero_particle /), &
2528  long_name="weight class number of each aerosol particle")
2529  call pmc_nc_write_integer_1d(ncid, aero_water_hyst_leg, &
2530  "aero_water_hyst_leg", (/ dimid_aero_particle /), &
2531  long_name="leg of the water hysteresis curve leg of each "&
2532  // "aerosol particle")
2533  call pmc_nc_write_real_1d(ncid, aero_num_conc, &
2534  "aero_num_conc", (/ dimid_aero_particle /), unit="m^{-3}", &
2535  long_name="number concentration for each particle")
2536  call pmc_nc_write_integer_1d(ncid, aero_id, &
2537  "aero_id", (/ dimid_aero_particle /), &
2538  long_name="unique ID number of each aerosol particle")
2539  call pmc_nc_write_real_1d(ncid, aero_least_create_time, &
2540  "aero_least_create_time", (/ dimid_aero_particle /), unit="s", &
2541  long_name="least creation time of each aerosol particle", &
2542  description="least (earliest) creation time of any original " &
2543  // "constituent particles that coagulated to form each " &
2544  // "particle, measured from the start of the simulation")
2545  call pmc_nc_write_real_1d(ncid, aero_greatest_create_time, &
2546  "aero_greatest_create_time", (/ dimid_aero_particle /), &
2547  unit="s", &
2548  long_name="greatest creation time of each aerosol particle", &
2549  description="greatest (latest) creation time of any original " &
2550  // "constituent particles that coagulated to form each " &
2551  // "particle, measured from the start of the simulation")
2552  if (record_optical) then
2553  call pmc_nc_write_real_1d(ncid, aero_absorb_cross_sect, &
2554  "aero_absorb_cross_sect", (/ dimid_aero_particle /), &
2555  unit="m^2", &
2556  long_name="optical absorption cross sections of each " &
2557  // "aerosol particle")
2558  call pmc_nc_write_real_1d(ncid, aero_scatter_cross_sect, &
2559  "aero_scatter_cross_sect", (/ dimid_aero_particle /), &
2560  unit="m^2", &
2561  long_name="optical scattering cross sections of each " &
2562  // "aerosol particle")
2563  call pmc_nc_write_real_1d(ncid, aero_asymmetry, &
2564  "aero_asymmetry", (/ dimid_aero_particle /), unit="1", &
2565  long_name="optical asymmetry parameters of each " &
2566  // "aerosol particle")
2567  call pmc_nc_write_real_1d(ncid, aero_refract_shell_real, &
2568  "aero_refract_shell_real", (/ dimid_aero_particle /), &
2569  unit="1", &
2570  long_name="real part of the refractive indices of the " &
2571  // "shell of each aerosol particle")
2572  call pmc_nc_write_real_1d(ncid, aero_refract_shell_imag, &
2573  "aero_refract_shell_imag", (/ dimid_aero_particle /), &
2574  unit="1", &
2575  long_name="imaginary part of the refractive indices of " &
2576  // "the shell of each aerosol particle")
2577  call pmc_nc_write_real_1d(ncid, aero_refract_core_real, &
2578  "aero_refract_core_real", (/ dimid_aero_particle /), &
2579  unit="1", &
2580  long_name="real part of the refractive indices of the core " &
2581  // "of each aerosol particle")
2582  call pmc_nc_write_real_1d(ncid, aero_refract_core_imag, &
2583  "aero_refract_core_imag", (/ dimid_aero_particle /), &
2584  unit="1", &
2585  long_name="imaginary part of the refractive indices of " &
2586  // "the core of each aerosol particle")
2587  call pmc_nc_write_real_1d(ncid, aero_core_vol, &
2588  "aero_core_vol", (/ dimid_aero_particle /), unit="m^3", &
2589  long_name="volume of the optical cores of each " &
2590  // "aerosol particle")
2591  end if
2592  end if
2593 
2594  ! FIXME: move this to aero_info_array.F90, together with
2595  ! aero_state_netcdf_dim_aero_removed() ?
2596  if (record_removals) then
2597  call aero_state_netcdf_dim_aero_removed(aero_state, ncid, &
2598  dimid_aero_removed)
2599  if (aero_info_array_n_item(aero_state%aero_info_array) >= 1) then
2600  do i_remove = 1,aero_info_array_n_item(aero_state%aero_info_array)
2601  aero_removed_id(i_remove) = &
2602  aero_state%aero_info_array%aero_info(i_remove)%id
2603  aero_removed_action(i_remove) = &
2604  aero_state%aero_info_array%aero_info(i_remove)%action
2605  aero_removed_other_id(i_remove) = &
2606  aero_state%aero_info_array%aero_info(i_remove)%other_id
2607  end do
2608  else
2609  aero_removed_id(1) = 0
2610  aero_removed_action(1) = aero_info_none
2611  aero_removed_other_id(1) = 0
2612  end if
2613  call pmc_nc_write_integer_1d(ncid, aero_removed_id, &
2614  "aero_removed_id", (/ dimid_aero_removed /), &
2615  long_name="ID of removed particles")
2616  call pmc_nc_write_integer_1d(ncid, aero_removed_action, &
2617  "aero_removed_action", (/ dimid_aero_removed /), &
2618  long_name="reason for particle removal", &
2619  description="valid is 0 (invalid entry), 1 (removed due to " &
2620  // "dilution), 2 (removed due to coagulation -- combined " &
2621  // "particle ID is in \c aero_removed_other_id), 3 (removed " &
2622  // "due to populating halving), or 4 (removed due to " &
2623  // "weighting changes")
2624  call pmc_nc_write_integer_1d(ncid, aero_removed_other_id, &
2625  "aero_removed_other_id", (/ dimid_aero_removed /), &
2626  long_name="ID of other particle involved in removal", &
2627  description="if <tt>aero_removed_action(i)</tt> is 2 " &
2628  // "(due to coagulation), then " &
2629  // "<tt>aero_removed_other_id(i)</tt> is the ID of the " &
2630  // "resulting combined particle, or 0 if the new particle " &
2631  // "was not created")
2632  end if
2633 
2634  end subroutine aero_state_output_netcdf
2635 
2636  ! this belongs in the subroutine above, but is outside because
2637  ! Doxygen 1.8.7 doesn't resolve references when multiple \page
2638  ! blocks are in one subroutine
2639 
2640  !> \page output_format_aero_removed Output File Format: Aerosol Particle Removal Information
2641  !!
2642  !! When an aerosol particle is introduced into the simulation it
2643  !! is assigned a unique ID number. This ID number will persist
2644  !! over time, allowing tracking of a paticular particle's
2645  !! evolution. If the \c record_removals variable in the input spec
2646  !! file is \c yes, then the every time a particle is removed from
2647  !! the simulation its removal will be recorded in the removal
2648  !! information.
2649  !!
2650  !! The removal information written at timestep \c n contains
2651  !! information about every particle ID that is present at time
2652  !! <tt>(n - 1)</tt> but not present at time \c n.
2653  !!
2654  !! The removal information is always written in the output files,
2655  !! even if no particles were removed in the previous
2656  !! timestep. Unfortunately, NetCDF files cannot contain arrays of
2657  !! length 0. In the case of no particles being removed, the \c
2658  !! aero_removed dimension will be set to 1 and
2659  !! <tt>aero_removed_action(1)</tt> will be 0 (\c AERO_INFO_NONE).
2660  !!
2661  !! When two particles coagulate, the ID number of the combined
2662  !! particle will be the ID particle of the largest constituent, if
2663  !! possible (weighting functions can make this impossible to
2664  !! achieve). A given particle ID may thus be lost due to
2665  !! coagulation (if the resulting combined particle has a different
2666  !! ID), or the ID may be preserved (as the ID of the combined
2667  !! particle). Only if the ID is lost will the particle be recorded
2668  !! in the removal information, and in this case
2669  !! <tt>aero_removed_action(i)</tt> will be 2 (\c AERO_INFO_COAG)
2670  !! and <tt>aero_removed_other_id(i)</tt> will be the ID number of
2671  !! the combined particle.
2672  !!
2673  !! The aerosol removal information NetCDF dimensions are:
2674  !! - \b aero_removed: number of aerosol particles removed from the
2675  !! simulation during the previous timestep (or 1, as described
2676  !! above)
2677  !!
2678  !! The aerosol removal information NetCDF variables are:
2679  !! - \b aero_removed (dim \c aero_removed): dummy dimension variable
2680  !! (no useful value)
2681  !! - \b aero_removed_id (dim \c aero_removed): the ID number of each
2682  !! removed particle
2683  !! - \b aero_removed_action (dim \c aero_removed): the reasons for
2684  !! removal for each particle, with values:
2685  !! - 0 (\c AERO_INFO_NONE): no information (invalid entry)
2686  !! - 1 (\c AERO_INFO_DILUTION): particle was removed due to dilution
2687  !! with outside air
2688  !! - 2 (\c AERO_INFO_COAG): particle was removed due to coagulation
2689  !! - 3 (\c AERO_INFO_HALVED): particle was removed due to halving of
2690  !! the aerosol population
2691  !! - 4 (\c AERO_INFO_WEIGHT): particle was removed due to adjustments
2692  !! in the particle's weighting function
2693  !! - \b aero_removed_other_id (dim \c aero_removed): the ID number of
2694  !! the combined particle formed by coagulation, if the removal reason
2695  !! was coagulation (2, \c AERO_INFO_COAG). May be 0, if the new
2696  !! coagulated particle was not created due to weighting.
2697 
2698 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2699 
2700  !> Read full state.
2701  subroutine aero_state_input_netcdf(aero_state, ncid, aero_data)
2702 
2703  !> aero_state to read.
2704  type(aero_state_t), intent(inout) :: aero_state
2705  !> NetCDF file ID, in data mode.
2706  integer, intent(in) :: ncid
2707  !> aero_data structure.
2708  type(aero_data_t), intent(in) :: aero_data
2709 
2710  integer :: dimid_aero_particle, dimid_aero_removed, n_info_item, n_part
2711  integer :: i_bin, i_part_in_bin, i_part, i_remove, status
2712  type(aero_particle_t) :: aero_particle
2713  character(len=1000) :: name
2714 
2715  real(kind=dp), allocatable :: aero_particle_mass(:,:)
2716  integer, allocatable :: aero_n_orig_part(:,:)
2717  integer, allocatable :: aero_particle_weight_group(:)
2718  integer, allocatable :: aero_particle_weight_class(:)
2719  real(kind=dp), allocatable :: aero_absorb_cross_sect(:)
2720  real(kind=dp), allocatable :: aero_scatter_cross_sect(:)
2721  real(kind=dp), allocatable :: aero_asymmetry(:)
2722  real(kind=dp), allocatable :: aero_refract_shell_real(:)
2723  real(kind=dp), allocatable :: aero_refract_shell_imag(:)
2724  real(kind=dp), allocatable :: aero_refract_core_real(:)
2725  real(kind=dp), allocatable :: aero_refract_core_imag(:)
2726  real(kind=dp), allocatable :: aero_core_vol(:)
2727  integer, allocatable :: aero_water_hyst_leg(:)
2728  real(kind=dp), allocatable :: aero_num_conc(:)
2729  integer, allocatable :: aero_id(:)
2730  real(kind=dp), allocatable :: aero_least_create_time(:)
2731  real(kind=dp), allocatable :: aero_greatest_create_time(:)
2732  integer, allocatable :: aero_removed_id(:)
2733  integer, allocatable :: aero_removed_action(:)
2734  integer, allocatable :: aero_removed_other_id(:)
2735 
2736  call aero_state_zero(aero_state)
2737 
2738  status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
2739  if (status == nf90_ebaddim) then
2740  ! no aero_particle dimension means no particles present
2741  call aero_weight_array_input_netcdf(aero_state%awa, ncid)
2742  return
2743  end if
2744  call pmc_nc_check(status)
2745  call pmc_nc_check(nf90_inquire_dimension(ncid, dimid_aero_particle, &
2746  name, n_part))
2747 
2748  call pmc_nc_read_real_2d(ncid, aero_particle_mass, &
2749  "aero_particle_mass")
2750  call pmc_nc_read_integer_2d(ncid, aero_n_orig_part, &
2751  "aero_n_orig_part")
2752  call pmc_nc_read_integer_1d(ncid, aero_particle_weight_group, &
2753  "aero_particle_weight_group")
2754  call pmc_nc_read_integer_1d(ncid, aero_particle_weight_class, &
2755  "aero_particle_weight_class")
2756  call pmc_nc_read_real_1d(ncid, aero_absorb_cross_sect, &
2757  "aero_absorb_cross_sect", must_be_present=.false.)
2758  call pmc_nc_read_real_1d(ncid, aero_scatter_cross_sect, &
2759  "aero_scatter_cross_sect", must_be_present=.false.)
2760  call pmc_nc_read_real_1d(ncid, aero_asymmetry, &
2761  "aero_asymmetry", must_be_present=.false.)
2762  call pmc_nc_read_real_1d(ncid, aero_refract_shell_real, &
2763  "aero_refract_shell_real", must_be_present=.false.)
2764  call pmc_nc_read_real_1d(ncid, aero_refract_shell_imag, &
2765  "aero_refract_shell_imag", must_be_present=.false.)
2766  call pmc_nc_read_real_1d(ncid, aero_refract_core_real, &
2767  "aero_refract_core_real", must_be_present=.false.)
2768  call pmc_nc_read_real_1d(ncid, aero_refract_core_imag, &
2769  "aero_refract_core_imag", must_be_present=.false.)
2770  call pmc_nc_read_real_1d(ncid, aero_core_vol, &
2771  "aero_core_vol", must_be_present=.false.)
2772  call pmc_nc_read_integer_1d(ncid, aero_water_hyst_leg, &
2773  "aero_water_hyst_leg")
2774  call pmc_nc_read_real_1d(ncid, aero_num_conc, &
2775  "aero_num_conc")
2776  call pmc_nc_read_integer_1d(ncid, aero_id, &
2777  "aero_id")
2778  call pmc_nc_read_real_1d(ncid, aero_least_create_time, &
2779  "aero_least_create_time")
2780  call pmc_nc_read_real_1d(ncid, aero_greatest_create_time, &
2781  "aero_greatest_create_time")
2782 
2783  call aero_weight_array_input_netcdf(aero_state%awa, ncid)
2784  call aero_state_set_n_part_ideal(aero_state, 0d0)
2785 
2786  do i_part = 1,n_part
2787  aero_particle%vol = aero_particle_mass(i_part, :) / aero_data%density
2788  aero_particle%n_orig_part = aero_n_orig_part(i_part, :)
2789  aero_particle%weight_group = aero_particle_weight_group(i_part)
2790  aero_particle%weight_class = aero_particle_weight_class(i_part)
2791  if (size(aero_absorb_cross_sect) == n_part) then
2792  aero_particle%absorb_cross_sect = aero_absorb_cross_sect(i_part)
2793  end if
2794  if (size(aero_scatter_cross_sect) == n_part) then
2795  aero_particle%scatter_cross_sect = aero_scatter_cross_sect(i_part)
2796  end if
2797  if (size(aero_asymmetry) == n_part) then
2798  aero_particle%asymmetry = aero_asymmetry(i_part)
2799  end if
2800  if ((size(aero_refract_shell_real) == n_part) &
2801  .and. (size(aero_refract_shell_imag) == n_part)) then
2802  aero_particle%refract_shell = &
2803  cmplx(aero_refract_shell_real(i_part), &
2804  aero_refract_shell_imag(i_part), kind=dc)
2805  end if
2806  if ((size(aero_refract_core_real) == n_part) &
2807  .and. (size(aero_refract_core_imag) == n_part)) then
2808  aero_particle%refract_core = cmplx(aero_refract_core_real(i_part), &
2809  aero_refract_core_imag(i_part), kind=dc)
2810  end if
2811  if (size(aero_core_vol) == n_part) then
2812  aero_particle%core_vol = aero_core_vol(i_part)
2813  end if
2814  aero_particle%water_hyst_leg = aero_water_hyst_leg(i_part)
2815  aero_particle%id = aero_id(i_part)
2816  aero_particle%least_create_time = aero_least_create_time(i_part)
2817  aero_particle%greatest_create_time = aero_greatest_create_time(i_part)
2818 
2819  call assert(314368871, almost_equal(aero_num_conc(i_part), &
2820  aero_weight_array_num_conc(aero_state%awa, aero_particle, &
2821  aero_data)))
2822 
2823  call aero_state_add_particle(aero_state, aero_particle, aero_data)
2824  end do
2825 
2826  call pmc_nc_read_integer_1d(ncid, aero_removed_id, &
2827  "aero_removed_id", must_be_present=.false.)
2828  call pmc_nc_read_integer_1d(ncid, aero_removed_action, &
2829  "aero_removed_action", must_be_present=.false.)
2830  call pmc_nc_read_integer_1d(ncid, aero_removed_other_id, &
2831  "aero_removed_other_id", must_be_present=.false.)
2832 
2833  n_info_item = size(aero_removed_id)
2834  if (n_info_item >= 1) then
2835  if ((n_info_item > 1) &
2836  .or. ((n_info_item == 1) .and. (aero_removed_id(1) /= 0))) then
2837  call aero_info_array_enlarge_to(aero_state%aero_info_array, &
2838  n_info_item)
2839  do i_remove = 1,n_info_item
2840  aero_state%aero_info_array%aero_info(i_remove)%id &
2841  = aero_removed_id(i_remove)
2842  aero_state%aero_info_array%aero_info(i_remove)%action &
2843  = aero_removed_action(i_remove)
2844  aero_state%aero_info_array%aero_info(i_remove)%other_id &
2845  = aero_removed_other_id(i_remove)
2846  end do
2847  end if
2848  end if
2849 
2850  end subroutine aero_state_input_netcdf
2851 
2852 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2853 
2854  !> Sorts the particles if necessary.
2855  subroutine aero_state_sort(aero_state, aero_data, bin_grid, all_procs_same)
2856 
2857  !> Aerosol state to sort.
2858  type(aero_state_t), intent(inout) :: aero_state
2859  !> Aerosol data.
2860  type(aero_data_t), intent(in) :: aero_data
2861  !> Bin grid.
2862  type(bin_grid_t), optional, intent(in) :: bin_grid
2863  !> Whether all processors should use the same bin grid.
2864  logical, optional, intent(in) :: all_procs_same
2865 
2866  call aero_sorted_remake_if_needed(aero_state%aero_sorted, aero_state%apa, &
2867  aero_data, aero_state%valid_sort, size(aero_state%awa%weight, 1), &
2868  size(aero_state%awa%weight, 2), bin_grid, all_procs_same)
2869  aero_state%valid_sort = .true.
2870 
2871  end subroutine aero_state_sort
2872 
2873 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2874 
2875  !> Check that aerosol state data is consistent.
2876  subroutine aero_state_check(aero_state, aero_data)
2877 
2878  !> Aerosol state to check.
2879  type(aero_state_t), intent(in) :: aero_state
2880  !> Aerosol data.
2881  type(aero_data_t), intent(in) :: aero_data
2882 
2883  logical, parameter :: continue_on_error = .false.
2884 
2885  call aero_particle_array_check(aero_state%apa, aero_data, &
2886  continue_on_error)
2887  if (aero_state%valid_sort) then
2888  call aero_sorted_check(aero_state%aero_sorted, aero_state%apa, &
2889  aero_data, size(aero_state%awa%weight, 1), &
2890  size(aero_state%awa%weight, 2), continue_on_error)
2891  end if
2892 
2893  end subroutine aero_state_check
2894 
2895 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2896 
2897 end module pmc_aero_state
2898 
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.
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
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 real(kind=dp) function aero_particle_diameter(aero_particle, aero_data)
Total diameter of the particle (m).
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
real(kind=dp) function, dimension(aero_state_n_part(aero_state)) aero_state_diameters(aero_state, aero_data)
Returns the diameters of all particles.
Definition: aero_state.F90:989
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.