PartMC  2.8.0
aero_particle.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012, 2016 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_particle module.
7 
8 !> The aero_particle_t structure and associated subroutines.
10 
11  use pmc_util
13  use pmc_aero_data
14  use pmc_spec_file
15  use pmc_env_state
16  use pmc_mpi
17 #ifdef PMC_USE_MPI
18  use mpi
19 #endif
20 
21  !> Maximum size of aero_components for a single particle.
22  integer, parameter :: max_aero_component_size = 20
23 
24  !> Single aerosol particle data structure.
25  !!
26  !! The \c vol array stores the total volumes of the different
27  !! species that make up the particle. This array must have length
28  !! equal to aero_data%%n_spec, so that \c vol(i) is the volume (in
29  !! m^3) of the i'th aerosol species.
31  !> Constituent species volumes [length aero_data_n_spec()] (m^3).
32  real(kind=dp), allocatable :: vol(:)
33  !> Weighting function group number (see \c aero_weight_array_t).
34  integer :: weight_group
35  !> Weighting function class number (see \c aero_weight_array_t).
36  integer :: weight_class
37  !> Absorption cross-section (m^2).
38  real(kind=dp), allocatable :: absorb_cross_sect(:)
39  !> Scattering cross-section (m^2).
40  real(kind=dp), allocatable :: scatter_cross_sect(:)
41  !> Asymmetry parameter (1).
42  real(kind=dp), allocatable :: asymmetry(:)
43  !> Refractive index of the shell (1).
44  complex(kind=dc), allocatable :: refract_shell(:)
45  !> Refractive index of the core (1).
46  complex(kind=dc), allocatable :: refract_core(:)
47  !> Volume of the core (m^3).
48  real(kind=dp) :: core_vol
49  !> Water hysteresis curve section (0 = lower, 1 = upper)
50  integer :: water_hyst_leg
51  !> Unique ID number.
52  integer(kind=8) :: id
53  !> Array of aerosol components.
54  type(aero_component_t), allocatable :: component(:)
55  !> First time a constituent was created (s).
56  real(kind=dp) :: least_create_time
57  !> Last time a constituent was created (s).
58  real(kind=dp) :: greatest_create_time
59  !> True number of primary particles contributing to particle.
60  integer :: n_primary_parts
61  end type aero_particle_t
62 
63  !> Next unique ID number to use for a particle.
64  integer(kind=8), save :: next_id = 1
65 
66 contains
67 
68 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
69 
70  !> Shift data from one aero_particle_t to another and free the first
71  !> one.
72  subroutine aero_particle_shift(aero_particle_from, aero_particle_to)
73 
74  !> Reference particle (will be deallocated on return).
75  type(aero_particle_t), intent(inout) :: aero_particle_from
76  !> Destination particle (not allocated on entry).
77  type(aero_particle_t), intent(inout) :: aero_particle_to
78 
79  call move_alloc(aero_particle_from%vol, aero_particle_to%vol)
80  aero_particle_to%weight_group = aero_particle_from%weight_group
81  aero_particle_to%weight_class = aero_particle_from%weight_class
82  call move_alloc(aero_particle_from%absorb_cross_sect, &
83  aero_particle_to%absorb_cross_sect)
84  call move_alloc(aero_particle_from%scatter_cross_sect, &
85  aero_particle_to%scatter_cross_sect)
86  call move_alloc(aero_particle_from%asymmetry, &
87  aero_particle_to%asymmetry)
88  call move_alloc(aero_particle_from%refract_shell, &
89  aero_particle_to%refract_shell)
90  call move_alloc(aero_particle_from%refract_core, &
91  aero_particle_to%refract_core)
92  aero_particle_to%core_vol = aero_particle_from%core_vol
93  aero_particle_to%water_hyst_leg = aero_particle_from%water_hyst_leg
94  aero_particle_to%id = aero_particle_from%id
95  call move_alloc(aero_particle_from%component, aero_particle_to%component)
96  aero_particle_to%least_create_time = aero_particle_from%least_create_time
97  aero_particle_to%greatest_create_time = &
98  aero_particle_from%greatest_create_time
99  aero_particle_to%n_primary_parts = aero_particle_from%n_primary_parts
100 
101  end subroutine aero_particle_shift
102 
103 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
104 
105  !> Resets an aero_particle to be zero.
106  subroutine aero_particle_zero(aero_particle, aero_data)
107 
108  !> Particle to zero.
109  type(aero_particle_t), intent(inout) :: aero_particle
110  !> Aerosol data.
111  type(aero_data_t), intent(in) :: aero_data
112 
113  call ensure_real_array_size(aero_particle%vol, aero_data_n_spec(aero_data))
114  aero_particle%vol = 0d0
115  aero_particle%weight_group = 0
116  aero_particle%weight_class = 0
117  call ensure_real_array_size(aero_particle%absorb_cross_sect, n_swbands)
118  aero_particle%absorb_cross_sect = 0d0
119  call ensure_real_array_size(aero_particle%scatter_cross_sect, n_swbands)
120  aero_particle%scatter_cross_sect = 0d0
121  call ensure_real_array_size(aero_particle%asymmetry, n_swbands)
122  aero_particle%asymmetry = 0d0
123  call ensure_complex_array_size(aero_particle%refract_shell, n_swbands)
124  aero_particle%refract_shell = (0d0, 0d0)
125  call ensure_complex_array_size(aero_particle%refract_core, n_swbands)
126  aero_particle%refract_core = (0d0, 0d0)
127  aero_particle%core_vol = 0d0
128  aero_particle%water_hyst_leg = 0
129  aero_particle%id = 0
130  if (allocated(aero_particle%component)) deallocate(aero_particle%component)
131  allocate(aero_particle%component(0))
132  aero_particle%least_create_time = 0d0
133  aero_particle%greatest_create_time = 0d0
134  aero_particle%n_primary_parts = 0
135 
136  end subroutine aero_particle_zero
137 
138 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
139 
140  !> Assigns a globally-unique new ID number to the particle.
141  subroutine aero_particle_new_id(aero_particle)
142 
143  !> Particle to set ID for.
144  type(aero_particle_t), intent(inout) :: aero_particle
145 
146  aero_particle%id = (next_id - 1) * pmc_mpi_size() + pmc_mpi_rank() + 1
147  next_id = next_id + 1
148 
149  end subroutine aero_particle_new_id
150 
151 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
152 
153  !> Sets the creation times for the particle.
154  subroutine aero_particle_set_create_time(aero_particle, create_time)
155 
156  !> Particle to set time for.
157  type(aero_particle_t), intent(inout) :: aero_particle
158  !> Creation time.
159  real(kind=dp), intent(in) :: create_time
160 
161  aero_particle%component(1)%create_time = create_time
162  aero_particle%least_create_time = create_time
163  aero_particle%greatest_create_time = create_time
164 
165  end subroutine aero_particle_set_create_time
166 
167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168 
169  !> Set the initial aero_component for a particle.
170  subroutine aero_particle_set_component(aero_particle, i_source, create_time)
171 
172  !> Particle to set time for.
173  type(aero_particle_t), intent(inout) :: aero_particle
174  !> Source number for the particle.
175  integer, intent(in) :: i_source
176  !> Creation time.
177  real(kind=dp), intent(in) :: create_time
178 
179  if (allocated(aero_particle%component)) deallocate(aero_particle%component)
180  allocate(aero_particle%component(1))
181  call aero_particle_set_source(aero_particle, i_source)
182  call aero_particle_set_create_time(aero_particle, create_time)
183 
184  end subroutine aero_particle_set_component
185 
186 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
187 
188  !> Gets the list of sources that a particle consists of.
189  subroutine aero_particle_get_component_sources(aero_particle, source_list)
190 
191  !> Particle.
192  type(aero_particle_t), intent(in) :: aero_particle
193  !> List of sources.
194  integer, intent(inout) :: source_list(:)
195 
196  integer :: i_comp, i_source
197 
198  do i_comp = 1,aero_particle_n_components(aero_particle)
199  i_source = aero_particle%component(i_comp)%source_id
200  source_list(i_source) = source_list(i_source) + 1
201  end do
202 
204 
205 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
206 
207  !> Sets the aerosol particle volumes.
208  subroutine aero_particle_set_vols(aero_particle, vols)
209 
210  !> Particle.
211  type(aero_particle_t), intent(inout) :: aero_particle
212  !> New volumes.
213  real(kind=dp), intent(in) :: vols(:)
214 
215  aero_particle%vol = vols
216 
217  end subroutine aero_particle_set_vols
218 
219 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
220 
221  !> Sets the aerosol particle source.
222  subroutine aero_particle_set_source(aero_particle, i_source)
223 
224  !> Particle.
225  type(aero_particle_t), intent(inout) :: aero_particle
226  !> Source number for the particle.
227  integer, intent(in) :: i_source
228 
229  aero_particle%component(1)%source_id = i_source
230 
231  end subroutine aero_particle_set_source
232 
233 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
234 
235  !> Sets the aerosol particle weight group.
236  subroutine aero_particle_set_weight(aero_particle, i_group, i_class)
237 
238  !> Particle.
239  type(aero_particle_t), intent(inout) :: aero_particle
240  !> Weight group number for the particle.
241  integer, intent(in), optional :: i_group
242  !> Weight class number for the particle.
243  integer, intent(in), optional :: i_class
244 
245  if (present(i_group)) aero_particle%weight_group = i_group
246  if (present(i_class)) aero_particle%weight_class = i_class
247 
248  end subroutine aero_particle_set_weight
249 
250 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
251 
252  !> Total mass of the particle (kg).
253  elemental real(kind=dp) function aero_particle_mass(aero_particle, aero_data)
254 
255  !> Particle.
256  type(aero_particle_t), intent(in) :: aero_particle
257  !> Aerosol data.
258  type(aero_data_t), intent(in) :: aero_data
259 
260  aero_particle_mass = sum(aero_particle%vol * aero_data%density)
261 
262  end function aero_particle_mass
263 
264 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
265 
266  !> Mass of a single species in the particle (kg).
267  elemental real(kind=dp) function aero_particle_species_mass(aero_particle, &
268  i_spec, aero_data)
269 
270  !> Particle.
271  type(aero_particle_t), intent(in) :: aero_particle
272  !> Species number to find mass of.
273  integer, intent(in) :: i_spec
274  !> Aerosol data.
275  type(aero_data_t), intent(in) :: aero_data
276 
277  aero_particle_species_mass = aero_particle%vol(i_spec) &
278  * aero_data%density(i_spec)
279 
280  end function aero_particle_species_mass
281 
282 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
283 
284  !> Mass of all species in the particle (kg).
285  function aero_particle_species_masses(aero_particle, aero_data)
286 
287  !> Particle.
288  type(aero_particle_t), intent(in) :: aero_particle
289  !> Aerosol data.
290  type(aero_data_t), intent(in) :: aero_data
291  !> Return value.
292  real(kind=dp) :: aero_particle_species_masses(aero_data_n_spec(aero_data))
293 
294  aero_particle_species_masses = aero_particle%vol * aero_data%density
295 
296  end function aero_particle_species_masses
297 
298 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
299 
300  !> Total moles in the particle (1).
301  elemental real(kind=dp) function aero_particle_moles(aero_particle, &
302  aero_data)
303 
304  !> Particle.
305  type(aero_particle_t), intent(in) :: aero_particle
306  !> Aerosol data.
307  type(aero_data_t), intent(in) :: aero_data
308 
309  aero_particle_moles = sum(aero_particle%vol * aero_data%density &
310  / aero_data%molec_weight)
311 
312  end function aero_particle_moles
313 
314 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
315 
316  !> Total volume of the particle (m^3).
317  elemental real(kind=dp) function aero_particle_volume(aero_particle)
318 
319  !> Particle.
320  type(aero_particle_t), intent(in) :: aero_particle
321 
322  aero_particle_volume = sum(aero_particle%vol)
323 
324  end function aero_particle_volume
325 
326 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
327 
328  !> Volume of a single species in the particle (m^3).
329  elemental real(kind=dp) function aero_particle_species_volume( &
330  aero_particle, i_spec)
331 
332  !> Particle.
333  type(aero_particle_t), intent(in) :: aero_particle
334  !> Species number to find volume of.
335  integer, intent(in) :: i_spec
336 
337  aero_particle_species_volume = aero_particle%vol(i_spec)
338 
339  end function aero_particle_species_volume
340 
341 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
342 
343  !> Total dry volume of the particle (m^3).
344  elemental real(kind=dp) function aero_particle_dry_volume(aero_particle, &
345  aero_data)
346 
347  !> Particle.
348  type(aero_particle_t), intent(in) :: aero_particle
349  !> Aerosol data.
350  type(aero_data_t), intent(in) :: aero_data
351 
352  integer :: i_spec
353 
355  do i_spec = 1,aero_data_n_spec(aero_data)
356  if (i_spec /= aero_data%i_water) then
358  + aero_particle%vol(i_spec)
359  end if
360  end do
361 
362  end function aero_particle_dry_volume
363 
364 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
365 
366  !> Total volume (dry or wet) of the particle (m^3).
367  elemental real(kind=dp) function aero_particle_volume_maybe_dry( &
368  aero_particle, aero_data, dry_volume)
369 
370  !> Particle.
371  type(aero_particle_t), intent(in) :: aero_particle
372  !> Aerosol data.
373  type(aero_data_t), intent(in) :: aero_data
374  !> Whether the desired volume is dry (otherwise wet).
375  logical, intent(in) :: dry_volume
376 
377  if (dry_volume) then
379  = aero_particle_dry_volume(aero_particle, aero_data)
380  else
382  end if
383 
384  end function aero_particle_volume_maybe_dry
385 
386 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
387 
388  !> Total radius of the particle (m).
389  elemental real(kind=dp) function aero_particle_radius(aero_particle, &
390  aero_data)
391 
392  !> Particle.
393  type(aero_particle_t), intent(in) :: aero_particle
394  !> Aerosol data.
395  type(aero_data_t), intent(in) :: aero_data
396 
398  aero_particle_volume(aero_particle))
399 
400  end function aero_particle_radius
401 
402 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
403 
404  !> Total dry radius of the particle (m).
405  elemental real(kind=dp) function aero_particle_dry_radius(aero_particle, &
406  aero_data)
407 
408  !> Particle.
409  type(aero_particle_t), intent(in) :: aero_particle
410  !> Aerosol data.
411  type(aero_data_t), intent(in) :: aero_data
412 
414  aero_particle_dry_volume(aero_particle, aero_data))
415 
416  end function aero_particle_dry_radius
417 
418 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
419 
420  !> Total diameter of the particle (m).
421  elemental real(kind=dp) function aero_particle_diameter(aero_particle, &
422  aero_data)
423 
424  !> Particle.
425  type(aero_particle_t), intent(in) :: aero_particle
426  !> Aerosol data.
427  type(aero_data_t), intent(in) :: aero_data
428 
429  aero_particle_diameter = 2d0 * aero_particle_radius(aero_particle, &
430  aero_data)
431 
432  end function aero_particle_diameter
433 
434 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
435 
436  !> Total dry diameter of the particle (m).
437  elemental real(kind=dp) function aero_particle_dry_diameter(aero_particle, &
438  aero_data)
439 
440  !> Particle.
441  type(aero_particle_t), intent(in) :: aero_particle
442  !> Aerosol data.
443  type(aero_data_t), intent(in) :: aero_data
444 
446  = 2d0 * aero_particle_dry_radius(aero_particle, aero_data)
447 
448  end function aero_particle_dry_diameter
449 
450 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
451 
452  !> Mobility diameter of the particle (m).
453  real(kind=dp) function aero_particle_mobility_diameter(aero_particle, &
454  aero_data, env_state)
455 
456  !> Particle.
457  type(aero_particle_t), intent(in) :: aero_particle
458  !> Aerosol data.
459  type(aero_data_t), intent(in) :: aero_data
460  !> Environment state.
461  type(env_state_t), intent(in) :: env_state
462 
463  real(kind=dp) :: volume, mobility_radius
464 
465  volume = aero_particle_volume(aero_particle)
466  mobility_radius = fractal_vol_to_mobility_rad(aero_data%fractal, &
467  volume, env_state%temp, env_state%pressure)
468  aero_particle_mobility_diameter = rad2diam(mobility_radius)
469 
471 
472 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
473 
474  !> Average density of the particle (kg/m^3).
475  real(kind=dp) function aero_particle_density(aero_particle, aero_data)
476 
477  !> Particle.
478  type(aero_particle_t), intent(in) :: aero_particle
479  !> Aerosol data.
480  type(aero_data_t), intent(in) :: aero_data
481 
482  aero_particle_density = aero_particle_mass(aero_particle, aero_data) &
483  / aero_particle_volume(aero_particle)
484 
485  end function aero_particle_density
486 
487 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
488 
489  !> Returns the volume-average of the non-water elements of quantity.
491  aero_particle, aero_data, quantity)
492 
493  !> Aerosol particle.
494  type(aero_particle_t), intent(in) :: aero_particle
495  !> Aerosol data.
496  type(aero_data_t), intent(in) :: aero_data
497  !> Quantity to average.
498  real(kind=dp), intent(in) :: quantity(:)
499 
500  real(kind=dp) :: ones(aero_data_n_spec(aero_data))
501 
502  ones = 1d0
504  aero_particle_total_solute_quantity(aero_particle, &
505  aero_data, quantity) &
506  / aero_particle_total_solute_quantity(aero_particle, aero_data, ones)
507 
509 
510 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
511 
512  !> Returns the volume-total of the non-water elements of quantity.
513  real(kind=dp) function aero_particle_total_solute_quantity(aero_particle, &
514  aero_data, quantity)
515 
516  !> Aerosol particle.
517  type(aero_particle_t), intent(in) :: aero_particle
518  !> Aerosol data.
519  type(aero_data_t), intent(in) :: aero_data
520  !> Quantity to total.
521  real(kind=dp), intent(in) :: quantity(:)
522 
523  real(kind=dp) total
524  integer i
525 
526  total = 0d0
527  do i = 1,aero_data_n_spec(aero_data)
528  if (i /= aero_data%i_water) then
529  total = total + aero_particle%vol(i) * quantity(i)
530  end if
531  end do
533 
535 
536 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
537 
538  !> Returns the water element of quantity.
539  real(kind=dp) function aero_particle_average_water_quantity(aero_particle, &
540  aero_data, quantity)
541 
542  !> Aerosol particle.
543  type(aero_particle_t), intent(in) :: aero_particle
544  !> Aerosol data.
545  type(aero_data_t), intent(in) :: aero_data
546  !> Quantity to average.
547  real(kind=dp), intent(in) :: quantity(:)
548 
549  call assert(420016623, aero_data%i_water > 0)
550  aero_particle_average_water_quantity = quantity(aero_data%i_water)
551 
553 
554 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
555 
556  !> Returns the volume-total of the water element of quantity.
557  real(kind=dp) function aero_particle_total_water_quantity(aero_particle, &
558  aero_data, quantity)
559 
560  !> Aerosol particle.
561  type(aero_particle_t), intent(in) :: aero_particle
562  !> Aerosol data.
563  type(aero_data_t), intent(in) :: aero_data
564  !> Quantity to total.
565  real(kind=dp), intent(in) :: quantity(:)
566 
567  call assert(223343210, aero_data%i_water > 0)
569  = aero_particle%vol(aero_data%i_water) &
570  * quantity(aero_data%i_water)
571 
573 
574 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
575 
576  !> Returns the water molecular weight.
577  !> (kg/mole)
578  real(kind=dp) function aero_particle_water_molec_weight(aero_data)
579 
580  !> Aerosol data.
581  type(aero_data_t), intent(in) :: aero_data
582 
583  call assert(772012490, aero_data%i_water > 0)
585  = aero_data%molec_weight(aero_data%i_water)
586 
588 
589 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
590 
591  !> Returns the average of the solute molecular weight (kg/mole).
592  real(kind=dp) function aero_particle_solute_molec_weight(aero_particle, &
593  aero_data)
594 
595  !> Aerosol particle.
596  type(aero_particle_t), intent(in) :: aero_particle
597  !> Aerosol data.
598  type(aero_data_t), intent(in) :: aero_data
599 
601  = aero_particle_average_solute_quantity(aero_particle, &
602  aero_data, aero_data%molec_weight)
603 
605 
606 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
607 
608  !> Returns the average of the solute ion number (1).
609  real(kind=dp) function aero_particle_solute_num_ions(aero_particle, &
610  aero_data)
611 
612  !> Aerosol data.
613  type(aero_data_t), intent(in) :: aero_data
614  !> Aerosol particle.
615  type(aero_particle_t), intent(in) :: aero_particle
616 
618  = aero_particle_average_solute_quantity(aero_particle, &
619  aero_data, real(aero_data%num_ions, kind=dp))
620 
621  end function aero_particle_solute_num_ions
622 
623 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
624 
625  !> Returns the water density (kg/m^3).
626  real(kind=dp) function aero_particle_water_density(aero_data)
627 
628  !> Aerosol data.
629  type(aero_data_t), intent(in) :: aero_data
630 
631  call assert(235482108, aero_data%i_water > 0)
632  aero_particle_water_density = aero_data%density(aero_data%i_water)
633 
634  end function aero_particle_water_density
635 
636 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
637 
638  !> Returns the average of the solute densities (kg/m^3).
639  real(kind=dp) function aero_particle_solute_density(aero_particle, &
640  aero_data)
641 
642  !> Aerosol data.
643  type(aero_data_t), intent(in) :: aero_data
644  !> Aerosol particle.
645  type(aero_particle_t), intent(in) :: aero_particle
646 
648  = aero_particle_average_solute_quantity(aero_particle, &
649  aero_data, aero_data%density)
650 
651  end function aero_particle_solute_density
652 
653 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
654 
655  !> Returns the water mass (kg).
656  real(kind=dp) function aero_particle_water_mass(aero_particle, aero_data)
657 
658  !> Aerosol data.
659  type(aero_data_t), intent(in) :: aero_data
660  !> Aerosol particle.
661  type(aero_particle_t), intent(in) :: aero_particle
662 
663  call assert(888636139, aero_data%i_water > 0)
664  aero_particle_water_mass = aero_particle%vol(aero_data%i_water) &
665  * aero_data%density(aero_data%i_water)
666 
667  end function aero_particle_water_mass
668 
669 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
670 
671  !> Returns the total solute mass (kg).
672  real(kind=dp) function aero_particle_solute_mass(aero_particle, aero_data)
673 
674  !> Aerosol data.
675  type(aero_data_t), intent(in) :: aero_data
676  !> Aerosol particle.
677  type(aero_particle_t), intent(in) :: aero_particle
678 
680  = aero_particle_total_solute_quantity(aero_particle, &
681  aero_data, aero_data%density)
682 
683  end function aero_particle_solute_mass
684 
685 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
686 
687  !> Returns the total solute volume (m^3).
688  real(kind=dp) function aero_particle_solute_volume(aero_particle, &
689  aero_data)
690 
691  !> Aerosol data.
692  type(aero_data_t), intent(in) :: aero_data
693  !> Aerosol particle.
694  type(aero_particle_t), intent(in) :: aero_particle
695 
696  real(kind=dp) :: ones(aero_data_n_spec(aero_data))
697 
698  ones = 1d0
700  = aero_particle_total_solute_quantity(aero_particle, &
701  aero_data, ones)
702 
703  end function aero_particle_solute_volume
704 
705 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
706 
707  !> Returns the total solute radius (m).
708  real(kind=dp) function aero_particle_solute_radius(aero_particle, &
709  aero_data)
710 
711  !> Aerosol particle.
712  type(aero_particle_t), intent(in) :: aero_particle
713  !> Aerosol data.
714  type(aero_data_t), intent(in) :: aero_data
715 
717  = aero_data_vol2rad(aero_data, &
718  aero_particle_solute_volume(aero_particle, aero_data))
719 
720  end function aero_particle_solute_radius
721 
722 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
723 
724  !> Returns the average of the solute kappas (1).
725  real(kind=dp) function aero_particle_solute_kappa(aero_particle, aero_data)
726 
727  !> Aerosol data.
728  type(aero_data_t), intent(in) :: aero_data
729  !> Aerosol particle.
730  type(aero_particle_t), intent(in) :: aero_particle
731 
732  real(kind=dp) :: kappa(aero_data_n_spec(aero_data)), m_w, rho_w, m_a, rho_a
733  integer :: i_spec
734 
735  m_w = aero_particle_water_molec_weight(aero_data)
736  rho_w = aero_particle_water_density(aero_data)
737  do i_spec = 1,aero_data_n_spec(aero_data)
738  if (i_spec == aero_data%i_water) then
739  kappa(i_spec) = 0d0
740  elseif (aero_data%num_ions(i_spec) > 0) then
741  call assert_msg(123681459, aero_data%kappa(i_spec) == 0d0, &
742  "species has nonzero num_ions and kappa: " &
743  // trim(aero_data%name(i_spec)))
744  m_a = aero_data%molec_weight(i_spec)
745  rho_a = aero_data%density(i_spec)
746  kappa(i_spec) = m_w * rho_a / (m_a * rho_w) &
747  * real(aero_data%num_ions(i_spec), kind=dp)
748  else
749  kappa(i_spec) = aero_data%kappa(i_spec)
750  end if
751  end do
753  = aero_particle_average_solute_quantity(aero_particle, &
754  aero_data, kappa)
755 
756  end function aero_particle_solute_kappa
757 
758 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
759 
760  !> Returns the approximate critical relative humidity (1).
761  real(kind=dp) function aero_particle_approx_crit_rel_humid(aero_particle, &
762  aero_data, env_state)
763 
764  !> Aerosol particle.
765  type(aero_particle_t), intent(in) :: aero_particle
766  !> Aerosol data.
767  type(aero_data_t), intent(in) :: aero_data
768  !> Environment state.
769  type(env_state_t), intent(in) :: env_state
770 
771  real(kind=dp) :: kappa, diam, c, a
772 
773  kappa = aero_particle_solute_kappa(aero_particle, aero_data)
774  c = sqrt(4d0 * env_state_a(env_state)**3 / 27d0)
775  diam = aero_particle_diameter(aero_particle, aero_data)
776  aero_particle_approx_crit_rel_humid = c / sqrt(kappa * diam**3) + 1d0
777 
779 
780 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
781 
782  !> Returns the critical relative humidity (1).
783  real(kind=dp) function aero_particle_crit_rel_humid(aero_particle, &
784  aero_data, env_state)
785 
786  !> Aerosol particle.
787  type(aero_particle_t), intent(in) :: aero_particle
788  !> Aerosol data.
789  type(aero_data_t), intent(in) :: aero_data
790  !> Environment state.
791  type(env_state_t), intent(in) :: env_state
792 
793  real(kind=dp) :: kappa, crit_diam, dry_diam, a
794 
795  a = env_state_a(env_state)
796  dry_diam = aero_particle_dry_diameter(aero_particle, aero_data)
797  crit_diam = aero_particle_crit_diameter(aero_particle, aero_data, &
798  env_state)
799  kappa = aero_particle_solute_kappa(aero_particle, aero_data)
800  if (kappa < 1d-30) then
801  aero_particle_crit_rel_humid = exp(a / crit_diam)
802  else
803  aero_particle_crit_rel_humid = (crit_diam**3 - dry_diam**3) &
804  / (crit_diam**3 - dry_diam**3 * (1 - kappa)) * exp(a / crit_diam)
805  end if
806 
807  end function aero_particle_crit_rel_humid
808 
809 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
810 
811  !> Returns the critical diameter (m).
812  !!
813  !! The method is as follows. We need to solve the polynomial
814  !! equation \f$f(D)\f$ given in Eq. (7) of
815  !!
816  !! N. Riemer, M. West, R. Zaveri, D. Easter, "Estimating black
817  !! carbon aging time-scales with a particle-resolved aerosol
818  !! model", Aerosol Science 41, 143-158, 2010.
819  !!
820  !! This is the equation:
821  !! \f[
822  !! f(D) = D^6 + c_4 D^4 + c_3 D^3 + c_0,
823  !! \f]
824  !! where \f$c_4 < 0\f$, \f$c_3 < 0\f$, and \f$c_0 > 0\f$. There is
825  !! unique solution for \f$D > D_{\rm dry}\f$, as shown in the above
826  !! paper.
827  !!
828  !! The polynomial has first derivative which factors like
829  !! \f[
830  !! f'(D) = 6 D^5 + 4 c_4 D^3 + 3 c_3 D^2
831  !! = (3 D^2 + 4 c_4) D^3 + (D^3 + c_3) 3 D^2.
832  !! \f]
833  !! The first term is positive for \f$D > (-4 c_4 / 3)^{1/2}\f$ and
834  !! the second is positive for \f$D > (-c_3)^{1/3}\f$. If we take
835  !! \f[
836  !! D_0 = max((-4 c_4 / 3)^{1/2}, (-c_3)^{1/3})
837  !! \f]
838  !! then we have that \f$f'(D) > 0\f$ for \f$D > D_0\f$. Similarly,
839  !! \f[
840  !! f''(D) = 30 D^4 + 12 c_4 D^2 + 6 c_3 D
841  !! = (5 D^2 + 4 c_4) 3 D^2 + (5 D^3 + 2 c_3) 3 D,
842  !! \f]
843  !! so \f$f''(D) > 0\f$ for \f$D > D_0\f$ (as this ensures that \f$D
844  !! > (-4 c_4 / 5)^{1/2}\f$ and \f$D > (-2 c_3 / 5)^{1/3}\f$).
845  !!
846  !! Thus for $D > D_0$ we have that the first and second derivatives
847  !! of $f(D)$ are positive, so Newton's method starting from
848  !! \f$D_0\f$ will converge quadratically. This is the scheme used
849  !! here.
850  real(kind=dp) function aero_particle_crit_diameter(aero_particle, &
851  aero_data, env_state)
852 
853  !> Aerosol particle.
854  type(aero_particle_t), intent(in) :: aero_particle
855  !> Aerosol data.
856  type(aero_data_t), intent(in) :: aero_data
857  !> Environment state.
858  type(env_state_t), intent(in) :: env_state
859 
860  integer, parameter :: crit_diam_max_iter = 100
861 
862  real(kind=dp) :: kappa, dry_diam, a, c4, c3, c0, d, f, df, dd
863  integer :: i_newton
864 
865  a = env_state_a(env_state)
866  dry_diam = aero_particle_dry_diameter(aero_particle, aero_data)
867  kappa = aero_particle_solute_kappa(aero_particle, aero_data)
868  if (kappa < 1d-30) then
869  ! bail out early for hydrophobic particles
870  aero_particle_crit_diameter = dry_diam
871  return
872  end if
873 
874  c4 = - 3d0 * dry_diam**3 * kappa / a
875  c3 = - dry_diam**3 * (2d0 - kappa)
876  c0 = dry_diam**6 * (1d0 - kappa)
877 
878  ! Newton's method for f(d) = d^6 + c4 d^4 + c3 d^3 + c0
879  d = max(sqrt(-4d0 / 3d0 * c4), (-c3)**(1d0/3d0))
880  do i_newton = 1,crit_diam_max_iter
881  f = d**6 + c4 * d**4 + c3 * d**3 + c0
882  df = 6 * d**5 + 4 * c4 * d**3 + 3 * c3 * d**2
883  dd = f / df
884  d = d - dd
885  if (abs(dd / d) < 1d-14) then
886  exit
887  end if
888  end do
889  call warn_assert_msg(408545686, i_newton < crit_diam_max_iter, &
890  "critical diameter Newton loop failed to converge")
891  call warn_assert_msg(353290871, d >= dry_diam, &
892  "critical diameter Newton loop converged to invalid solution")
894 
895  end function aero_particle_crit_diameter
896 
897 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
898 
899  !> Coagulate two particles together to make a new one. The new
900  !> particle will not have its ID set.
901  subroutine aero_particle_coagulate(aero_particle_1, &
902  aero_particle_2, aero_particle_new)
903 
904  !> First particle.
905  type(aero_particle_t), intent(in) :: aero_particle_1
906  !> Second particle.
907  type(aero_particle_t), intent(in) :: aero_particle_2
908  !> Result particle.
909  type(aero_particle_t), intent(inout) :: aero_particle_new
910 
911  integer :: n_comp_1, n_comp_2, n_comp_1_new, n_comp_2_new, i
912  type(aero_component_t), allocatable :: new_aero_component(:)
913  integer :: n_swbands
914  integer, allocatable :: sample(:)
915 
916  call assert(203741686, size(aero_particle_1%vol) &
917  == size(aero_particle_2%vol))
918  aero_particle_new%vol = aero_particle_1%vol + aero_particle_2%vol
919  aero_particle_new%weight_group = 0
920  aero_particle_new%weight_class = 0
921  n_swbands = size(aero_particle_1%absorb_cross_sect)
922  call ensure_real_array_size(aero_particle_new%absorb_cross_sect, n_swbands)
923  call ensure_real_array_size(aero_particle_new%scatter_cross_sect, &
924  n_swbands)
925  call ensure_real_array_size(aero_particle_new%asymmetry, n_swbands)
926  call ensure_complex_array_size(aero_particle_new%refract_shell, n_swbands)
927  call ensure_complex_array_size(aero_particle_new%refract_core, n_swbands)
928  aero_particle_new%absorb_cross_sect = 0d0
929  aero_particle_new%scatter_cross_sect = 0d0
930  aero_particle_new%asymmetry = 0d0
931  aero_particle_new%refract_shell = (0d0, 0d0)
932  aero_particle_new%refract_core = (0d0, 0d0)
933  aero_particle_new%core_vol = 0d0
934  if ((aero_particle_1%water_hyst_leg == 1) &
935  .and. (aero_particle_2%water_hyst_leg == 1)) then
936  aero_particle_new%water_hyst_leg = 1
937  else
938  aero_particle_new%water_hyst_leg = 0
939  end if
940  aero_particle_new%id = 0
941  n_comp_1 = aero_particle_n_components(aero_particle_1)
942  call assert_msg(465791384, aero_particle_1%n_primary_parts >= &
943  n_comp_1, 'n_primary_parts = ' &
944  // trim(integer_to_string(aero_particle_1%n_primary_parts)) &
945  // ' is less than n_components = ' &
946  // trim(integer_to_string(n_comp_1)))
947  n_comp_2 = aero_particle_n_components(aero_particle_2)
948  call assert_msg(465791385, aero_particle_2%n_primary_parts >= &
949  n_comp_2, 'n_primary_parts = ' &
950  // trim(integer_to_string(aero_particle_2%n_primary_parts)) &
951  // ' is less than n_components = ' &
952  // trim(integer_to_string(n_comp_2)))
953  if (n_comp_1 + n_comp_2 > max_aero_component_size) then
954  n_comp_1_new = prob_round(real(aero_particle_1%n_primary_parts, &
955  kind=dp) / (aero_particle_1%n_primary_parts &
956  + aero_particle_2%n_primary_parts) * max_aero_component_size)
957  n_comp_2_new = max_aero_component_size - n_comp_1_new
958  allocate(new_aero_component(max_aero_component_size))
959  call sample_without_replacement(n_comp_1_new, n_comp_1, sample)
960  do i = 1,n_comp_1_new
961  new_aero_component(i) = aero_particle_1%component(sample(i))
962  end do
963  call sample_without_replacement(n_comp_2_new, n_comp_2, sample)
964  do i = 1,n_comp_2_new
965  new_aero_component(i+n_comp_1_new) = aero_particle_2%component( &
966  sample(i))
967  end do
968  aero_particle_new%component = new_aero_component
969  else
970  new_aero_component = [aero_particle_1%component, &
971  aero_particle_2%component]
972  call move_alloc(new_aero_component, aero_particle_new%component)
973  end if
974  aero_particle_new%least_create_time = &
975  min(aero_particle_1%least_create_time, &
976  aero_particle_2%least_create_time)
977  aero_particle_new%greatest_create_time = &
978  max(aero_particle_1%greatest_create_time, &
979  aero_particle_2%greatest_create_time)
980 
981  aero_particle_new%n_primary_parts = aero_particle_1%n_primary_parts &
982  + aero_particle_2%n_primary_parts
983 
984  end subroutine aero_particle_coagulate
985 
986 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
987 
988  !> Return the number of aerosol components, or -1 if uninitialized.
989  elemental integer function aero_particle_n_components(particle)
990 
991  !> Particle.
992  type(aero_particle_t), intent(in) :: particle
993 
994  if (allocated(particle%component)) then
995  aero_particle_n_components = size(particle%component)
996  else
998  end if
999 
1000  end function aero_particle_n_components
1001 
1002 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1003 
1004  !> Determines the number of bytes required to pack the given value.
1005  integer function pmc_mpi_pack_size_aero_particle(val)
1006 
1007  !> Value to pack.
1008  type(aero_particle_t), intent(in) :: val
1009 
1010  integer :: i
1011 
1013  pmc_mpi_pack_size_real_array(val%vol) &
1014  + pmc_mpi_pack_size_integer(val%weight_group) &
1015  + pmc_mpi_pack_size_integer(val%weight_class) &
1016  + pmc_mpi_pack_size_real_array(val%absorb_cross_sect) &
1017  + pmc_mpi_pack_size_real_array(val%scatter_cross_sect) &
1018  + pmc_mpi_pack_size_real_array(val%asymmetry) &
1019  + pmc_mpi_pack_size_complex_array(val%refract_shell) &
1020  + pmc_mpi_pack_size_complex_array(val%refract_core) &
1021  + pmc_mpi_pack_size_real(val%core_vol) &
1022  + pmc_mpi_pack_size_integer(val%water_hyst_leg) &
1023  + pmc_mpi_pack_size_integer64(val%id) &
1025  + pmc_mpi_pack_size_real(val%least_create_time) &
1026  + pmc_mpi_pack_size_real(val%greatest_create_time) &
1027  + pmc_mpi_pack_size_integer(val%n_primary_parts)
1028 
1029  do i = 1,aero_particle_n_components(val)
1031  + pmc_mpi_pack_size_aero_component(val%component(i))
1032  end do
1033 
1034  end function pmc_mpi_pack_size_aero_particle
1035 
1036 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1037 
1038  !> Packs the given value into the buffer, advancing position.
1039  subroutine pmc_mpi_pack_aero_particle(buffer, position, val)
1040 
1041  !> Memory buffer.
1042  character, intent(inout) :: buffer(:)
1043  !> Current buffer position.
1044  integer, intent(inout) :: position
1045  !> Value to pack.
1046  type(aero_particle_t), intent(in) :: val
1047 
1048 #ifdef PMC_USE_MPI
1049  integer :: prev_position, i
1050 
1051  prev_position = position
1052  call pmc_mpi_pack_real_array(buffer, position, val%vol)
1053  call pmc_mpi_pack_integer(buffer, position, val%weight_group)
1054  call pmc_mpi_pack_integer(buffer, position, val%weight_class)
1055  call pmc_mpi_pack_real_array(buffer, position, val%absorb_cross_sect)
1056  call pmc_mpi_pack_real_array(buffer, position, val%scatter_cross_sect)
1057  call pmc_mpi_pack_real_array(buffer, position, val%asymmetry)
1058  call pmc_mpi_pack_complex_array(buffer, position, val%refract_shell)
1059  call pmc_mpi_pack_complex_array(buffer, position, val%refract_core)
1060  call pmc_mpi_pack_real(buffer, position, val%core_vol)
1061  call pmc_mpi_pack_integer(buffer, position, val%water_hyst_leg)
1062  call pmc_mpi_pack_integer64(buffer, position, val%id)
1063  call pmc_mpi_pack_integer(buffer, position, &
1065  do i = 1,aero_particle_n_components(val)
1066  call pmc_mpi_pack_aero_component(buffer, position, val%component(i))
1067  end do
1068  call pmc_mpi_pack_real(buffer, position, val%least_create_time)
1069  call pmc_mpi_pack_real(buffer, position, val%greatest_create_time)
1070  call pmc_mpi_pack_integer(buffer, position, val%n_primary_parts)
1071  call assert(810223998, position - prev_position &
1073 #endif
1074 
1075  end subroutine pmc_mpi_pack_aero_particle
1076 
1077 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1078 
1079  !> Unpacks the given value from the buffer, advancing position.
1080  subroutine pmc_mpi_unpack_aero_particle(buffer, position, val)
1081 
1082  !> Memory buffer.
1083  character, intent(inout) :: buffer(:)
1084  !> Current buffer position.
1085  integer, intent(inout) :: position
1086  !> Value to pack.
1087  type(aero_particle_t), intent(inout) :: val
1088 
1089 #ifdef PMC_USE_MPI
1090  integer :: prev_position, n_components, i
1091 
1092  prev_position = position
1093  call pmc_mpi_unpack_real_array(buffer, position, val%vol)
1094  call pmc_mpi_unpack_integer(buffer, position, val%weight_group)
1095  call pmc_mpi_unpack_integer(buffer, position, val%weight_class)
1096  call pmc_mpi_unpack_real_array(buffer, position, val%absorb_cross_sect)
1097  call pmc_mpi_unpack_real_array(buffer, position, val%scatter_cross_sect)
1098  call pmc_mpi_unpack_real_array(buffer, position, val%asymmetry)
1099  call pmc_mpi_unpack_complex_array(buffer, position, val%refract_shell)
1100  call pmc_mpi_unpack_complex_array(buffer, position, val%refract_core)
1101  call pmc_mpi_unpack_real(buffer, position, val%core_vol)
1102  call pmc_mpi_unpack_integer(buffer, position, val%water_hyst_leg)
1103  call pmc_mpi_unpack_integer64(buffer, position, val%id)
1104  call pmc_mpi_unpack_integer(buffer, position, n_components)
1105  if (n_components > -1) then
1106  allocate(val%component(n_components))
1107  end if
1108  do i = 1,n_components
1109  call pmc_mpi_unpack_aero_component(buffer, position, val%component(i))
1110  end do
1111  call pmc_mpi_unpack_real(buffer, position, val%least_create_time)
1112  call pmc_mpi_unpack_real(buffer, position, val%greatest_create_time)
1113  call pmc_mpi_unpack_integer(buffer, position, val%n_primary_parts)
1114  call assert(287447241, position - prev_position &
1116 #endif
1117 
1118  end subroutine pmc_mpi_unpack_aero_particle
1119 
1120 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1121 
1122  !> Check that the particle data is consistent.
1123  subroutine aero_particle_check(aero_particle, aero_data, &
1124  continue_on_error)
1125 
1126  !> Aerosol particle to check.
1127  type(aero_particle_t), intent(in) :: aero_particle
1128  !> Aerosol data.
1129  type(aero_data_t), intent(in) :: aero_data
1130  !> Whether to continue despite error.
1131  logical, intent(in) :: continue_on_error
1132 
1133  if (allocated(aero_particle%vol)) then
1134  if (size(aero_particle%vol) /= aero_data_n_spec(aero_data)) then
1135  write(0, *) 'ERROR aero_particle A:'
1136  write(0, *) 'size(aero_particle%vol)', size(aero_particle%vol)
1137  write(0, *) 'aero_data_n_spec(aero_data)', &
1138  aero_data_n_spec(aero_data)
1139  call assert(185878626, continue_on_error)
1140  end if
1141  end if
1142 
1143  end subroutine aero_particle_check
1144 
1145 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1146 
1147 end module pmc_aero_particle
pmc_aero_particle::aero_particle_solute_molec_weight
real(kind=dp) function aero_particle_solute_molec_weight(aero_particle, aero_data)
Returns the average of the solute molecular weight (kg/mole).
Definition: aero_particle.F90:594
pmc_aero_particle::aero_particle_t
Single aerosol particle data structure.
Definition: aero_particle.F90:30
pmc_rand::prob_round
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
pmc_aero_particle::aero_particle_set_weight
subroutine aero_particle_set_weight(aero_particle, i_group, i_class)
Sets the aerosol particle weight group.
Definition: aero_particle.F90:237
pmc_util::rad2diam
real(kind=dp) elemental function rad2diam(r)
Convert radius (m) to diameter (m).
Definition: util.F90:253
pmc_aero_data::aero_data_n_spec
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Definition: aero_data.F90:246
pmc_mpi::pmc_mpi_size
integer function pmc_mpi_size()
Returns the total number of processes.
Definition: mpi.F90:134
pmc_aero_particle::aero_particle_average_solute_quantity
real(kind=dp) function aero_particle_average_solute_quantity(aero_particle, aero_data, quantity)
Returns the volume-average of the non-water elements of quantity.
Definition: aero_particle.F90:492
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_aero_particle::aero_particle_solute_kappa
real(kind=dp) function aero_particle_solute_kappa(aero_particle, aero_data)
Returns the average of the solute kappas (1).
Definition: aero_particle.F90:726
pmc_aero_particle::aero_particle_crit_rel_humid
real(kind=dp) function aero_particle_crit_rel_humid(aero_particle, aero_data, env_state)
Returns the critical relative humidity (1).
Definition: aero_particle.F90:785
pmc_aero_particle::aero_particle_solute_num_ions
real(kind=dp) function aero_particle_solute_num_ions(aero_particle, aero_data)
Returns the average of the solute ion number (1).
Definition: aero_particle.F90:611
pmc_mpi::pmc_mpi_pack_size_real
integer function pmc_mpi_pack_size_real(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:385
pmc_aero_particle::aero_particle_total_solute_quantity
real(kind=dp) function aero_particle_total_solute_quantity(aero_particle, aero_data, quantity)
Returns the volume-total of the non-water elements of quantity.
Definition: aero_particle.F90:515
pmc_aero_particle::aero_particle_set_create_time
subroutine aero_particle_set_create_time(aero_particle, create_time)
Sets the creation times for the particle.
Definition: aero_particle.F90:155
pmc_aero_particle
The aero_particle_t structure and associated subroutines.
Definition: aero_particle.F90:9
pmc_aero_particle::aero_particle_moles
elemental real(kind=dp) function aero_particle_moles(aero_particle, aero_data)
Total moles in the particle (1).
Definition: aero_particle.F90:303
pmc_aero_particle::aero_particle_new_id
subroutine aero_particle_new_id(aero_particle)
Assigns a globally-unique new ID number to the particle.
Definition: aero_particle.F90:142
pmc_aero_component
The aero_particle_t structure and associated subroutines.
Definition: aero_component.F90:9
pmc_aero_particle::aero_particle_species_volume
elemental real(kind=dp) function aero_particle_species_volume(aero_particle, i_spec)
Volume of a single species in the particle (m^3).
Definition: aero_particle.F90:331
pmc_aero_particle::aero_particle_solute_density
real(kind=dp) function aero_particle_solute_density(aero_particle, aero_data)
Returns the average of the solute densities (kg/m^3).
Definition: aero_particle.F90:641
pmc_aero_particle::aero_particle_crit_diameter
real(kind=dp) function aero_particle_crit_diameter(aero_particle, aero_data, env_state)
Returns the critical diameter (m).
Definition: aero_particle.F90:852
pmc_mpi::pmc_mpi_unpack_integer64
subroutine pmc_mpi_unpack_integer64(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1164
pmc_aero_particle::aero_particle_diameter
elemental real(kind=dp) function aero_particle_diameter(aero_particle, aero_data)
Total diameter of the particle (m).
Definition: aero_particle.F90:423
pmc_aero_particle::aero_particle_water_density
real(kind=dp) function aero_particle_water_density(aero_data)
Returns the water density (kg/m^3).
Definition: aero_particle.F90:627
pmc_aero_particle::aero_particle_n_components
elemental integer function aero_particle_n_components(particle)
Return the number of aerosol components, or -1 if uninitialized.
Definition: aero_particle.F90:990
pmc_mpi::pmc_mpi_pack_integer64
subroutine pmc_mpi_pack_integer64(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:736
pmc_aero_particle::aero_particle_solute_radius
real(kind=dp) function aero_particle_solute_radius(aero_particle, aero_data)
Returns the total solute radius (m).
Definition: aero_particle.F90:710
pmc_util::warn_assert_msg
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition: util.F90:60
pmc_env_state::env_state_a
real(kind=dp) function env_state_a(env_state)
Condensation parameter.
Definition: env_state.F90:228
pmc_mpi::pmc_mpi_rank
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
pmc_aero_component::pmc_mpi_pack_size_aero_component
integer function pmc_mpi_pack_size_aero_component(val)
Determines the number of bytes required to pack the given value.
Definition: aero_component.F90:34
pmc_aero_particle::aero_particle_solute_volume
real(kind=dp) function aero_particle_solute_volume(aero_particle, aero_data)
Returns the total solute volume (m^3).
Definition: aero_particle.F90:690
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
pmc_util::ensure_real_array_size
subroutine ensure_real_array_size(x, n, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1058
pmc_aero_component::aero_component_t
Aerosol particle component data structure.
Definition: aero_component.F90:21
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_mpi::pmc_mpi_pack_real
subroutine pmc_mpi_pack_real(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:761
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:104
pmc_aero_particle::aero_particle_density
real(kind=dp) function aero_particle_density(aero_particle, aero_data)
Average density of the particle (kg/m^3).
Definition: aero_particle.F90:476
pmc_mpi::pmc_mpi_pack_size_real_array
integer function pmc_mpi_pack_size_real_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:527
pmc_aero_particle::aero_particle_solute_mass
real(kind=dp) function aero_particle_solute_mass(aero_particle, aero_data)
Returns the total solute mass (kg).
Definition: aero_particle.F90:673
pmc_aero_particle::aero_particle_mobility_diameter
real(kind=dp) function aero_particle_mobility_diameter(aero_particle, aero_data, env_state)
Mobility diameter of the particle (m).
Definition: aero_particle.F90:455
pmc_aero_data::n_swbands
integer, parameter n_swbands
Definition: aero_data.F90:34
pmc_aero_particle::aero_particle_set_vols
subroutine aero_particle_set_vols(aero_particle, vols)
Sets the aerosol particle volumes.
Definition: aero_particle.F90:209
pmc_mpi::pmc_mpi_pack_complex_array
subroutine pmc_mpi_pack_complex_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:863
pmc_aero_particle::aero_particle_shift
subroutine aero_particle_shift(aero_particle_from, aero_particle_to)
Shift data from one aero_particle_t to another and free the first one.
Definition: aero_particle.F90:73
pmc_mpi::pmc_mpi_unpack_complex_array
subroutine pmc_mpi_unpack_complex_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1292
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:767
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:78
pmc_aero_particle::aero_particle_approx_crit_rel_humid
real(kind=dp) function aero_particle_approx_crit_rel_humid(aero_particle, aero_data, env_state)
Returns the approximate critical relative humidity (1).
Definition: aero_particle.F90:763
pmc_aero_particle::aero_particle_species_masses
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).
Definition: aero_particle.F90:286
pmc_aero_component::pmc_mpi_pack_aero_component
subroutine pmc_mpi_pack_aero_component(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_component.F90:48
pmc_aero_particle::aero_particle_water_mass
real(kind=dp) function aero_particle_water_mass(aero_particle, aero_data)
Returns the water mass (kg).
Definition: aero_particle.F90:657
pmc_aero_particle::pmc_mpi_unpack_aero_particle
subroutine pmc_mpi_unpack_aero_particle(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_particle.F90:1081
pmc_aero_particle::max_aero_component_size
integer, parameter max_aero_component_size
Maximum size of aero_components for a single particle.
Definition: aero_particle.F90:22
pmc_aero_particle::aero_particle_zero
subroutine aero_particle_zero(aero_particle, aero_data)
Resets an aero_particle to be zero.
Definition: aero_particle.F90:107
pmc_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_aero_data::aero_data_vol2rad
real(kind=dp) elemental function aero_data_vol2rad(aero_data, v)
Convert mass-equivalent volume (m^3) to geometric radius (m).
Definition: aero_data.F90:102
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
pmc_mpi::pmc_mpi_unpack_integer
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1139
pmc_rand::sample_without_replacement
subroutine sample_without_replacement(sample_size, population_size, sample)
Generates a vector of indices for a sample set from a population without replacement.
Definition: rand.F90:516
pmc_aero_particle::aero_particle_get_component_sources
subroutine aero_particle_get_component_sources(aero_particle, source_list)
Gets the list of sources that a particle consists of.
Definition: aero_particle.F90:190
pmc_aero_particle::pmc_mpi_pack_size_aero_particle
integer function pmc_mpi_pack_size_aero_particle(val)
Determines the number of bytes required to pack the given value.
Definition: aero_particle.F90:1006
pmc_aero_particle::aero_particle_volume_maybe_dry
elemental real(kind=dp) function aero_particle_volume_maybe_dry(aero_particle, aero_data, dry_volume)
Total volume (dry or wet) of the particle (m^3).
Definition: aero_particle.F90:369
pmc_aero_particle::aero_particle_dry_diameter
elemental real(kind=dp) function aero_particle_dry_diameter(aero_particle, aero_data)
Total dry diameter of the particle (m).
Definition: aero_particle.F90:439
pmc_aero_particle::aero_particle_check
subroutine aero_particle_check(aero_particle, aero_data, continue_on_error)
Check that the particle data is consistent.
Definition: aero_particle.F90:1125
pmc_aero_particle::aero_particle_dry_radius
elemental real(kind=dp) function aero_particle_dry_radius(aero_particle, aero_data)
Total dry radius of the particle (m).
Definition: aero_particle.F90:407
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_aero_particle::aero_particle_dry_volume
elemental real(kind=dp) function aero_particle_dry_volume(aero_particle, aero_data)
Total dry volume of the particle (m^3).
Definition: aero_particle.F90:346
pmc_aero_particle::aero_particle_radius
elemental real(kind=dp) function aero_particle_radius(aero_particle, aero_data)
Total radius of the particle (m).
Definition: aero_particle.F90:391
pmc_mpi::pmc_mpi_pack_size_integer64
integer function pmc_mpi_pack_size_integer64(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:365
pmc_aero_particle::next_id
integer(kind=8), save next_id
Next unique ID number to use for a particle.
Definition: aero_particle.F90:64
pmc_mpi::pmc_mpi_pack_size_integer
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:345
pmc_aero_component::pmc_mpi_unpack_aero_component
subroutine pmc_mpi_unpack_aero_component(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_component.F90:72
pmc_aero_particle::aero_particle_species_mass
elemental real(kind=dp) function aero_particle_species_mass(aero_particle, i_spec, aero_data)
Mass of a single species in the particle (kg).
Definition: aero_particle.F90:269
pmc_mpi::pmc_mpi_pack_integer
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:711
pmc_aero_particle::aero_particle_water_molec_weight
real(kind=dp) function aero_particle_water_molec_weight(aero_data)
Returns the water molecular weight. (kg/mole)
Definition: aero_particle.F90:579
pmc_mpi::pmc_mpi_pack_size_complex_array
integer function pmc_mpi_pack_size_complex_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:467
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_mpi::pmc_mpi_pack_real_array
subroutine pmc_mpi_pack_real_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:927
pmc_aero_particle::aero_particle_set_source
subroutine aero_particle_set_source(aero_particle, i_source)
Sets the aerosol particle source.
Definition: aero_particle.F90:223
pmc_util::ensure_complex_array_size
subroutine ensure_complex_array_size(x, n, only_grow)
Allocate or reallocate the given array to ensure it is of the given size, preserving any data and/or ...
Definition: util.F90:1298
pmc_mpi::pmc_mpi_unpack_real
subroutine pmc_mpi_unpack_real(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1189
pmc_mpi::pmc_mpi_unpack_real_array
subroutine pmc_mpi_unpack_real_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1356
pmc_aero_particle::aero_particle_coagulate
subroutine aero_particle_coagulate(aero_particle_1, aero_particle_2, aero_particle_new)
Coagulate two particles together to make a new one. The new particle will not have its ID set.
Definition: aero_particle.F90:903
pmc_aero_particle::aero_particle_mass
elemental real(kind=dp) function aero_particle_mass(aero_particle, aero_data)
Total mass of the particle (kg).
Definition: aero_particle.F90:254
pmc_aero_particle::pmc_mpi_pack_aero_particle
subroutine pmc_mpi_pack_aero_particle(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_particle.F90:1040
pmc_aero_particle::aero_particle_volume
elemental real(kind=dp) function aero_particle_volume(aero_particle)
Total volume of the particle (m^3).
Definition: aero_particle.F90:318
pmc_aero_particle::aero_particle_total_water_quantity
real(kind=dp) function aero_particle_total_water_quantity(aero_particle, aero_data, quantity)
Returns the volume-total of the water element of quantity.
Definition: aero_particle.F90:559
pmc_aero_particle::aero_particle_set_component
subroutine aero_particle_set_component(aero_particle, i_source, create_time)
Set the initial aero_component for a particle.
Definition: aero_particle.F90:171
pmc_aero_particle::aero_particle_average_water_quantity
real(kind=dp) function aero_particle_average_water_quantity(aero_particle, aero_data, quantity)
Returns the water element of quantity.
Definition: aero_particle.F90:541