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