PartMC  2.3.0
aero_weight.F90
Go to the documentation of this file.
1 ! Copyright (C) 2010-2015 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_weight module.
7 
8 !> The aero_weight_t structure and associated subroutines.
10 
11  use pmc_util
12  use pmc_constants
13  use pmc_rand
14  use pmc_spec_file
16  use pmc_netcdf
17  use pmc_mpi
18 #ifdef PMC_USE_MPI
19  use mpi
20 #endif
21 
22  !> Type code for an undefined or invalid weighting.
23  integer, parameter :: AERO_WEIGHT_TYPE_INVALID = 0
24  !> Type code for no (or flat) weighting.
25  integer, parameter :: AERO_WEIGHT_TYPE_NONE = 1
26  !> Type code for power function weighting.
27  integer, parameter :: AERO_WEIGHT_TYPE_POWER = 2
28  !> Type code for MFA weighting.
29  integer, parameter :: AERO_WEIGHT_TYPE_MFA = 3
30 
31  !> An aerosol size distribution weighting function.
33  !> Weight type (given by module constants).
34  integer :: type
35  !> Overall weight magnitude (m^{-3}).
36  real(kind=dp) :: magnitude
37  !> Exponent for "power" weight.
38  real(kind=dp) :: exponent
39  end type aero_weight_t
40 
41 contains
42 
43 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
44 
45  !> Allocates an aero_weight.
46  elemental subroutine aero_weight_allocate(aero_weight)
47 
48  !> Aerosol weight.
49  type(aero_weight_t), intent(out) :: aero_weight
50 
51  call aero_weight_zero(aero_weight)
52 
53  end subroutine aero_weight_allocate
54 
55 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
56 
57  !> Allocates an aero_weight of the given size.
58  elemental subroutine aero_weight_allocate_size(aero_weight)
59 
60  !> Aerosol weight.
61  type(aero_weight_t), intent(out) :: aero_weight
62 
63  call aero_weight_zero(aero_weight)
64 
65  end subroutine aero_weight_allocate_size
66 
67 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
68 
69  !> Free all storage.
70  elemental subroutine aero_weight_deallocate(aero_weight)
71 
72  !> Aerosol weight.
73  type(aero_weight_t), intent(in) :: aero_weight
74 
75  end subroutine aero_weight_deallocate
76 
77 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
78 
79  !> Zeros the contents of the \c aero_weight.
80  elemental subroutine aero_weight_zero(aero_weight)
81 
82  !> Aerosol weight.
83  type(aero_weight_t), intent(inout) :: aero_weight
84 
85  aero_weight%type = aero_weight_type_invalid
86  aero_weight%magnitude = 1d0
87  aero_weight%exponent = 0d0
88 
89  end subroutine aero_weight_zero
90 
91 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
92 
93  !> Sets the \c aero_weight to a non-zero normalized value.
94  elemental subroutine aero_weight_normalize(aero_weight)
95 
96  !> Aerosol weight.
97  type(aero_weight_t), intent(inout) :: aero_weight
98 
99  aero_weight%magnitude = 1d0
100 
101  end subroutine aero_weight_normalize
102 
103 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
104 
105  !> Copy an aero_weight.
106  elemental subroutine aero_weight_copy(aero_weight_from, aero_weight_to)
107 
108  !> Aerosol weight original.
109  type(aero_weight_t), intent(in) :: aero_weight_from
110  !> Aerosol weight copy.
111  type(aero_weight_t), intent(inout) :: aero_weight_to
112 
113  aero_weight_to%type = aero_weight_from%type
114  aero_weight_to%magnitude = aero_weight_from%magnitude
115  aero_weight_to%exponent = aero_weight_from%exponent
116 
117  end subroutine aero_weight_copy
118 
119 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
120 
121  !> Scale the weight by the given fraction, so <tt>new_weight =
122  !> old_weight * factor</tt>.
123  elemental subroutine aero_weight_scale(aero_weight, factor)
124 
125  !> Aerosol weight to scale.
126  type(aero_weight_t), intent(inout) :: aero_weight
127  !> Factor to scale by.
128  real(kind=dp), intent(in) :: factor
129 
130  aero_weight%magnitude = aero_weight%magnitude * factor
131 
132  end subroutine aero_weight_scale
133 
134 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
135 
136  !> Combine \c aero_weight_delta into \c aero_weight with a harmonic mean.
137  elemental subroutine aero_weight_combine(aero_weight, aero_weight_delta)
138 
139  !> Aerosol weight to add into.
140  type(aero_weight_t), intent(inout) :: aero_weight
141  !> Aerosol weight to add from.
142  type(aero_weight_t), intent(in) :: aero_weight_delta
143 
144  aero_weight%magnitude = harmonic_mean(aero_weight%magnitude, &
145  aero_weight_delta%magnitude)
146 
147  end subroutine aero_weight_combine
148 
149 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
150 
151  !> Adjust source and destination weights to reflect moving \c
152  !> sample_prop proportion of particles from \c aero_weight_from to
153  !> \c aero_weight_to.
154  elemental subroutine aero_weight_shift(aero_weight_from, &
155  aero_weight_to, sample_prop, overwrite_to)
156 
157  !> Aerosol weight to shift from.
158  type(aero_weight_t), intent(inout) :: aero_weight_from
159  !> Aerosol weight to shift to.
160  type(aero_weight_t), intent(inout) :: aero_weight_to
161  !> Proportion of particles being transfered.
162  real(kind=dp), intent(in) :: sample_prop
163  !> Whether to overwrite the destination weight (default: no).
164  logical, intent(in), optional :: overwrite_to
165 
166  real(kind=dp) :: magnitude_transfer
167  logical :: use_overwrite_to
168 
169  magnitude_transfer = aero_weight_from%magnitude / sample_prop
170  aero_weight_from%magnitude = harmonic_mean(aero_weight_from%magnitude, &
171  - magnitude_transfer)
172  use_overwrite_to = .false.
173  if (present(overwrite_to)) then
174  if (overwrite_to) then
175  use_overwrite_to = .true.
176  end if
177  end if
178  if (use_overwrite_to) then
179  aero_weight_to%magnitude = magnitude_transfer
180  else
181  aero_weight_to%magnitude = harmonic_mean(aero_weight_to%magnitude, &
182  magnitude_transfer)
183  end if
184 
185  end subroutine aero_weight_shift
186 
187 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
188 
189  !> Compute the number concentration at a given radius (m^{-3}).
190  real(kind=dp) function aero_weight_num_conc_at_radius(aero_weight, radius)
191 
192  !> Aerosol weight.
193  type(aero_weight_t), intent(in) :: aero_weight
194  !> Radius to compute number concentration at (m).
195  real(kind=dp), intent(in) :: radius
196 
198  if (aero_weight%type == aero_weight_type_none) then
199  aero_weight_num_conc_at_radius = aero_weight%magnitude
200  elseif ((aero_weight%type == aero_weight_type_power) &
201  .or. (aero_weight%type == aero_weight_type_mfa)) then
202  aero_weight_num_conc_at_radius = aero_weight%magnitude &
203  * radius**aero_weight%exponent
204  else
205  call die_msg(700421478, "unknown aero_weight type: " &
206  // trim(integer_to_string(aero_weight%type)))
207  end if
208 
209  end function aero_weight_num_conc_at_radius
210 
211 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
212 
213  !> Compute the number concentration for a particle (m^{-3}).
214  real(kind=dp) function aero_weight_num_conc(aero_weight, aero_particle)
215 
216  !> Aerosol weight.
217  type(aero_weight_t), intent(in) :: aero_weight
218  !> Aerosol particle to compute number concentration for.
219  type(aero_particle_t), intent(in) :: aero_particle
220 
222  aero_particle_radius(aero_particle))
223 
224  end function aero_weight_num_conc
225 
226 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
227 
228  !> Compute the radius at a given number concentration (m). Inverse
229  !> of aero_weight_num_conc_at_radius().
230  real(kind=dp) function aero_weight_radius_at_num_conc(aero_weight, num_conc)
231 
232  !> Aerosol weight.
233  type(aero_weight_t), intent(in) :: aero_weight
234  !> Number concentration to compute the radius at (m^{-3}).
235  real(kind=dp), intent(in) :: num_conc
236 
238  if (aero_weight%type == aero_weight_type_none) then
239  call die_msg(545688537, "cannot invert weight type 'none'")
240  elseif ((aero_weight%type == aero_weight_type_power) &
241  .or. (aero_weight%type == aero_weight_type_mfa)) then
242  call assert_msg(902242996, aero_weight%exponent /= 0d0, &
243  "cannot invert weight with zero exponent")
245  = exp(log(num_conc / aero_weight%magnitude) / aero_weight%exponent)
246  else
247  call die_msg(307766845, "unknown aero_weight type: " &
248  // trim(integer_to_string(aero_weight%type)))
249  end if
250 
251  end function aero_weight_radius_at_num_conc
252 
253 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
254 
255  !> Ensures that a weight function exponent is valid.
256  subroutine aero_weight_check_valid_exponent(aero_weight)
257 
258  !> Aerosol weight array.
259  type(aero_weight_t), intent(in) :: aero_weight
260 
261  ! check we know about all the weight type
262  call assert(269952052, (aero_weight%type == aero_weight_type_none) &
263  .or. (aero_weight%type == aero_weight_type_power) &
264  .or. (aero_weight%type == aero_weight_type_mfa))
265  if (aero_weight%type == aero_weight_type_none) then
266  call assert(853998284, aero_weight%exponent == 0d0)
267  end if
268  if (aero_weight%type == aero_weight_type_mfa) then
269  call assert(829651126, aero_weight%exponent == -3d0)
270  end if
271 
272  end subroutine aero_weight_check_valid_exponent
273 
274 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
275 
276  !> Determine whether a weight function is monotone increasing,
277  !> monotone decreasing, or neither.
278  subroutine aero_weight_check_monotonicity(aero_weight, &
279  monotone_increasing, monotone_decreasing)
280 
281  !> Aerosol weight array.
282  type(aero_weight_t), intent(in) :: aero_weight
283  !> Whether all weights are monotone increasing.
284  logical, intent(out) :: monotone_increasing
285  !> Whether all weights are monotone decreasing.
286  logical, intent(out) :: monotone_decreasing
287 
288  ! check we know about all the weight types
289  call assert(610698264, (aero_weight%type == aero_weight_type_none) &
290  .or. (aero_weight%type == aero_weight_type_power) &
291  .or. (aero_weight%type == aero_weight_type_mfa))
292  monotone_increasing = .true.
293  monotone_decreasing = .true.
294  if (aero_weight%type == aero_weight_type_power) then
295  if (aero_weight%exponent < 0d0) then
296  monotone_increasing = .false.
297  end if
298  if (aero_weight%exponent > 0d0) then
299  monotone_decreasing = .false.
300  end if
301  end if
302  if (aero_weight%type == aero_weight_type_mfa) then
303  monotone_increasing = .false.
304  end if
305 
306  end subroutine aero_weight_check_monotonicity
307 
308 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
309 
310  !> Read an aero_weight from a spec file.
311  subroutine spec_file_read_aero_weight(file, aero_weight)
312 
313  !> Spec file.
314  type(spec_file_t), intent(inout) :: file
315  !> Aerosol weight.
316  type(aero_weight_t), intent(inout) :: aero_weight
317 
318  character(len=SPEC_LINE_MAX_VAR_LEN) :: weight_type
319 
320  !> \page input_format_aero_weight Input File Format: Aerosol Weighting Function
321  !!
322  !! For efficiency the aerosol population can be weighted so that
323  !! the true number distribution \f$n(D)\f$ is given by
324  !! \f[ n(D) = w(D) c(D) \f]
325  !! where \f$w(D)\f$ is a fixed weighting function, \f$c(D)\f$ is
326  !! the computational (simulated) number distribution, and \f$D\f$
327  !! is the diameter. Thus a large value of \f$w(D)\f$ means that
328  !! relatively few computational particles are used at diameter
329  !! \f$D\f$, while a small value of \f$w(D)\f$ means that
330  !! relatively many computational particles will be used at that
331  !! diameter.
332  !!
333  !! The aerosol weighting function is specified by the parameters:
334  !! - \b weight (string): the type of weighting function --- must
335  !! be one of: "none" for no weighting (\f$w(D) = 1\f$);
336  !! "power" for a power-law weighting (\f$w(D) = D^\alpha\f$),
337  !! or "mfa" for the mass flow algorithm weighting (\f$w(D) =
338  !! D^{-3}\f$ with dependent coagulation particle removal)
339  !! - if the \c weight is \c power then the next parameter is:
340  !! - \b exponent (real, dimensionless): the exponent
341  !! \f$\alpha\f$ in the power law relationship --- setting
342  !! the \c exponent to 0 is equivalent to no weighting, while
343  !! setting the exponent negative uses more computational
344  !! particles at larger diameters and setting the exponent
345  !! positive uses more computational particles at smaller
346  !! diameters; in practice exponents between 0 and -3 are
347  !! most useful
348  !!
349  !! See also:
350  !! - \ref spec_file_format --- the input file text format
351 
352  call spec_file_read_string(file, 'weight', weight_type)
353  aero_weight%magnitude = 0d0
354  if (trim(weight_type) == 'none') then
355  aero_weight%type = aero_weight_type_none
356  elseif (trim(weight_type) == 'power') then
357  aero_weight%type = aero_weight_type_power
358  call spec_file_read_real(file, 'exponent', aero_weight%exponent)
359  elseif (trim(weight_type) == 'mfa') then
360  aero_weight%type = aero_weight_type_mfa
361  aero_weight%exponent = -3d0
362  else
363  call spec_file_die_msg(456342050, file, "unknown weight_type: " &
364  // trim(weight_type))
365  end if
366 
367  end subroutine spec_file_read_aero_weight
368 
369 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
370 
371  !> Determines the number of bytes required to pack the given value.
372  integer function pmc_mpi_pack_size_aero_weight(val)
373 
374  !> Value to pack.
375  type(aero_weight_t), intent(in) :: val
376 
378  pmc_mpi_pack_size_integer(val%type) &
379  + pmc_mpi_pack_size_real(val%magnitude) &
380  + pmc_mpi_pack_size_real(val%exponent)
381 
382  end function pmc_mpi_pack_size_aero_weight
383 
384 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
385 
386  !> Packs the given value into the buffer, advancing position.
387  subroutine pmc_mpi_pack_aero_weight(buffer, position, val)
388 
389  !> Memory buffer.
390  character, intent(inout) :: buffer(:)
391  !> Current buffer position.
392  integer, intent(inout) :: position
393  !> Value to pack.
394  type(aero_weight_t), intent(in) :: val
395 
396 #ifdef PMC_USE_MPI
397  integer :: prev_position
398 
399  prev_position = position
400  call pmc_mpi_pack_integer(buffer, position, val%type)
401  call pmc_mpi_pack_real(buffer, position, val%magnitude)
402  call pmc_mpi_pack_real(buffer, position, val%exponent)
403  call assert(579699255, &
404  position - prev_position <= pmc_mpi_pack_size_aero_weight(val))
405 #endif
406 
407  end subroutine pmc_mpi_pack_aero_weight
408 
409 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
410 
411  !> Unpacks the given value from the buffer, advancing position.
412  subroutine pmc_mpi_unpack_aero_weight(buffer, position, val)
413 
414  !> Memory buffer.
415  character, intent(inout) :: buffer(:)
416  !> Current buffer position.
417  integer, intent(inout) :: position
418  !> Value to pack.
419  type(aero_weight_t), intent(inout) :: val
420 
421 #ifdef PMC_USE_MPI
422  integer :: prev_position
423 
424  prev_position = position
425  call pmc_mpi_unpack_integer(buffer, position, val%type)
426  call pmc_mpi_unpack_real(buffer, position, val%magnitude)
427  call pmc_mpi_unpack_real(buffer, position, val%exponent)
428  call assert(639056899, &
429  position - prev_position <= pmc_mpi_pack_size_aero_weight(val))
430 #endif
431 
432  end subroutine pmc_mpi_unpack_aero_weight
433 
434 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
435 
436 end module pmc_aero_weight
real(kind=dp) function aero_weight_num_conc_at_radius(aero_weight, radius)
Compute the number concentration at a given radius (m^{-3}).
integer function pmc_mpi_pack_size_real(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:367
subroutine pmc_mpi_unpack_aero_weight(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
An input file with extra data for printing messages.
Definition: spec_file.F90:59
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:347
elemental subroutine aero_weight_normalize(aero_weight)
Sets the aero_weight to a non-zero normalized value.
Definition: aero_weight.F90:94
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
subroutine pmc_mpi_unpack_real(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:796
elemental subroutine aero_weight_copy(aero_weight_from, aero_weight_to)
Copy an aero_weight.
The aero_weight_t structure and associated subroutines.
Definition: aero_weight.F90:9
Physical constants.
Definition: constants.F90:9
The aero_particle_t structure and associated subroutines.
elemental subroutine aero_weight_allocate(aero_weight)
Allocates an aero_weight.
Definition: aero_weight.F90:46
elemental subroutine aero_weight_shift(aero_weight_from, aero_weight_to, sample_prop, overwrite_to)
Adjust source and destination weights to reflect moving sample_prop proportion of particles from aero...
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:76
subroutine aero_weight_check_valid_exponent(aero_weight)
Ensures that a weight function exponent is valid.
elemental subroutine aero_weight_scale(aero_weight, factor)
Scale the weight by the given fraction, so new_weight = old_weight * factor.
elemental real(kind=dp) function aero_particle_radius(aero_particle)
Total radius of the particle (m).
An aerosol size distribution weighting function.
Definition: aero_weight.F90:32
Random number generators.
Definition: rand.F90:9
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:534
integer function pmc_mpi_pack_size_aero_weight(val)
Determines the number of bytes required to pack the given value.
elemental subroutine aero_weight_combine(aero_weight, aero_weight_delta)
Combine aero_weight_delta into aero_weight with a harmonic mean.
subroutine pmc_mpi_pack_aero_weight(buffer, position, val)
Packs the given value into the buffer, advancing position.
Common utility subroutines.
Definition: util.F90:9
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:771
Wrapper functions for MPI.
Definition: mpi.F90:13
subroutine spec_file_read_string(file, name, var)
Read a string from a spec file that must have a given name.
Definition: spec_file.F90:625
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:743
Single aerosol particle data structure.
Reading formatted text input.
Definition: spec_file.F90:43
real(kind=dp) function aero_weight_radius_at_num_conc(aero_weight, num_conc)
Compute the radius at a given number concentration (m). Inverse of aero_weight_num_conc_at_radius().
Wrapper functions for NetCDF. These all take a NetCDF ncid in data mode and return with it again in d...
Definition: netcdf.F90:11
subroutine pmc_mpi_pack_real(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:559
subroutine aero_weight_check_monotonicity(aero_weight, monotone_increasing, monotone_decreasing)
Determine whether a weight function is monotone increasing, monotone decreasing, or neither...
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:73
elemental real(kind=dp) function harmonic_mean(x1, x2)
Compute the harmonic mean of two numbers.
Definition: util.F90:1552
elemental subroutine aero_weight_allocate_size(aero_weight)
Allocates an aero_weight of the given size.
Definition: aero_weight.F90:58
real(kind=dp) function aero_weight_num_conc(aero_weight, aero_particle)
Compute the number concentration for a particle (m^{-3}).
elemental subroutine aero_weight_zero(aero_weight)
Zeros the contents of the aero_weight.
Definition: aero_weight.F90:80
subroutine spec_file_read_real(file, name, var)
Read a real number from a spec file that must have the given name.
Definition: spec_file.F90:579
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:102
elemental subroutine aero_weight_deallocate(aero_weight)
Free all storage.
Definition: aero_weight.F90:70