PartMC  2.6.1
aero_data.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012, 2016, 2021 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_data module.
7 
8 !> The aero_data_t structure and associated subroutines.
10 
11  use pmc_spec_file
12  use pmc_mpi
13  use pmc_util
14  use pmc_fractal
15  use pmc_netcdf
16 #ifdef PMC_USE_CAMP
17  use camp_util, only: string_t
18  use camp_camp_core
19  use camp_chem_spec_data
20  use camp_aero_rep_data
21  use camp_aero_rep_single_particle
22  use camp_property
23 #endif
24 #ifdef PMC_USE_MPI
25  use mpi
26 #endif
27 
28  integer, parameter :: aero_name_len = 50
29  integer, parameter :: aero_source_name_len = 100
30 
31  !> Aerosol material properties and associated data.
32  !!
33  !! The data in this structure is constant, as it represents physical
34  !! quantities that cannot change over time.
35  !!
36  !! Each aerosol species is identified by an index <tt>i =
37  !! 1,...,aero_data_n_spec(aero_data)</tt>. Then \c name(i) is the name of
38  !! that species, \c density(i) is its density, etc. The ordering of the
39  !! species is arbitrary and should not be relied upon (currently it is the
40  !! order in the species data file). The only exception is that it is
41  !! possible to find out which species is water from the \c i_water
42  !! variable.
43  !!
44  !! The names of the aerosol species are not important to PartMC, as
45  !! only the material properties are used. The names are used for
46  !! input and output, and also for communication with MOSAIC. For the
47  !! MOSAIC interface to work correctly the species must be named the
48  !! same, but without the \c _a suffix.
50  !> Water species number (0 if water is not a species).
51  integer :: i_water
52  !> Length \c aero_data_n_spec(aero_data), species names.
53  character(len=AERO_NAME_LEN), allocatable :: name(:)
54  !> Length \c aero_data_n_spec(aero_data), mosaic_index(i) a positive
55  !> integer giving the mosaic index of species i, or 0 if there is no match.
56  integer, allocatable :: mosaic_index(:)
57  !> Length \c aero_data_n_spec(aero_data), densities (kg m^{-3}).
58  real(kind=dp), allocatable :: density(:)
59  !> Length \c aero_data_n_spec(aero_data), number of ions in solute.
60  integer, allocatable :: num_ions(:)
61  !> Length \c aero_data_n_spec(aero_data), molecular weights (kg mol^{-1}).
62  real(kind=dp), allocatable :: molec_weight(:)
63  !> Length \c aero_data_n_spec(aero_data), kappas (1).
64  real(kind=dp), allocatable :: kappa(:)
65  !> Length \c aero_data_n_source(aero_data), source names.
66  character(len=AERO_SOURCE_NAME_LEN), allocatable :: source_name(:)
67  !> Fractal particle parameters.
68  type(fractal_t) :: fractal
69 #ifdef PMC_USE_CAMP
70  !> CAMP aerosol representation pointer
71  class(aero_rep_data_t), pointer :: aero_rep_ptr
72  !> Aerosol species ids on the camp chem state array for the first
73  !! computational particle
74  integer, allocatable :: camp_particle_spec_id(:)
75  !> Number of elements on the camp chem state array per computational
76  !! particle
77  integer :: camp_particle_state_size = -1
78  contains
79  !> Get the index on the CAMP state array for a specified species and
80  !! computation particle
81  procedure :: camp_spec_id
82 #endif
83  end type aero_data_t
84 
85 contains
86 
87 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
88 
89  !> Convert mass-equivalent volume \f$V\f$ (m^3) to geometric radius
90  !> \f$R_{\rm geo}\f$ (m).
91  real(kind=dp) elemental function aero_data_vol2rad(aero_data, v)
92 
93  !> Aero data structure.
94  type(aero_data_t), intent(in) :: aero_data
95  !> Particle mass-equivalent volume (m^3).
96  real(kind=dp), intent(in) :: v
97 
98  aero_data_vol2rad = fractal_vol2rad(aero_data%fractal, v)
99 
100  end function aero_data_vol2rad
101 
102 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
103 
104  !> Convert mass-equivalent volume \f$V\f$ (m^3) to geometric diameter
105  !> \f$D_{\rm geo}\f$ (m).
106  real(kind=dp) elemental function aero_data_vol2diam(aero_data, v)
107 
108  !> Aero data structure.
109  type(aero_data_t), intent(in) :: aero_data
110  !> Particle mass-equivalent volume (m^3).
111  real(kind=dp), intent(in) :: v
112 
113  aero_data_vol2diam = fractal_vol2diam(aero_data%fractal, v)
114 
115  end function aero_data_vol2diam
116 
117 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
118 
119  !> Convert geometric radius \f$R_{\rm geo}\f$ (m) to mass-equivalent volume
120  !> \f$V\f$ (m^3).
121  real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
122 
123  !> Aero data structure.
124  type(aero_data_t), intent(in) :: aero_data
125  !> Geometric radius (m).
126  real(kind=dp), intent(in) :: r
127 
128  aero_data_rad2vol = fractal_rad2vol(aero_data%fractal, r)
129 
130  end function aero_data_rad2vol
131 
132 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
133 
134  !> Convert geometric diameter \f$D_{\rm geo}\f$ (m) to
135  !> mass-equivalent volume \f$V\f$ (m^3).
136  real(kind=dp) elemental function aero_data_diam2vol(aero_data, d)
137 
138  !> Aero data structure.
139  type(aero_data_t), intent(in) :: aero_data
140  !> Geometric diameter (m).
141  real(kind=dp), intent(in) :: d
142 
143  aero_data_diam2vol = fractal_rad2vol(aero_data%fractal, diam2rad(d))
144 
145  end function aero_data_diam2vol
146 
147 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
148 
149  !> Convert mass-equivalent volume \f$V\f$ (m^3) to number of
150  !> monomers \f$N\f$ in a fractal particle cluster.
151  !!
152  !! Based on Eq. 5 in Naumann [2003].
153  real(kind=dp) elemental function aero_data_vol_to_num_of_monomers( &
154  aero_data, v)
155 
156  !> Aero data structure.
157  type(aero_data_t), intent(in) :: aero_data
158  !> Particle mass-equivalent volume (m^3).
159  real(kind=dp), intent(in) :: v
160 
162  aero_data%fractal, v)
163 
165 
166 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
167 
168  !> Convert mass-equivalent volume \f$V\f$ (m^3) to mobility equivalent
169  !> radius \f$R_{\rm me}\f$ (m).
170  !!
171  !! Based on Eq. 5, 21 and 30 in Naumann [2003].
172  real(kind=dp) function aero_data_vol_to_mobility_rad(aero_data, v, temp, &
173  pressure)
174 
175  !> Aero data structure.
176  type(aero_data_t), intent(in) :: aero_data
177  !> Particle mass-equivalent volume (m^3).
178  real(kind=dp), intent(in) :: v
179  !> Temperature (K).
180  real(kind=dp), intent(in) :: temp
181  !> Pressure (Pa).
182  real(kind=dp), intent(in) :: pressure
183 
185  aero_data%fractal, v, temp, pressure)
186 
187  end function aero_data_vol_to_mobility_rad
188 
189 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
190 
191  !> Convert mobility equivalent radius \f$R_{\rm me}\f$ (m) to
192  !> mass-equivalent volume \f$V\f$ (m^3).
193  real(kind=dp) function aero_data_mobility_rad_to_vol(aero_data, &
194  mobility_rad, temp, pressure)
195 
196  !> Aero data structure.
197  type(aero_data_t), intent(in) :: aero_data
198  !> Mobility equivalent radius (m).
199  real(kind=dp), intent(in) :: mobility_rad
200  !> Temperature (K).
201  real(kind=dp), intent(in) :: temp
202  !> Pressure (Pa).
203  real(kind=dp), intent(in) :: pressure
204 
206  aero_data%fractal, mobility_rad, temp, pressure)
207 
208  end function aero_data_mobility_rad_to_vol
209 
210 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
211 
212  !> Convert mobility equivalent radius \f$R_{\rm me}\f$ (m) to
213  !> geometric radius \f$R_{\rm geo}\f$ (m^3).
214  real(kind=dp) function aero_data_mobility_rad_to_geometric_rad(aero_data, &
215  mobility_rad, temp, pressure)
216 
217  !> Aero data structure.
218  type(aero_data_t), intent(in) :: aero_data
219  !> Mobility equivalent radius (m).
220  real(kind=dp), intent(in) :: mobility_rad
221  !> Temperature (K).
222  real(kind=dp), intent(in) :: temp
223  !> Pressure (Pa).
224  real(kind=dp), intent(in) :: pressure
225 
227  = fractal_mobility_rad_to_geometric_rad(aero_data%fractal, &
228  mobility_rad, temp, pressure)
229 
231 
232 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
233 
234  !> Return the number of aerosol species, or -1 if uninitialized.
235  elemental integer function aero_data_n_spec(aero_data)
236 
237  !> Aero data structure.
238  type(aero_data_t), intent(in) :: aero_data
239 
240  if (allocated(aero_data%name)) then
241  aero_data_n_spec = size(aero_data%name)
242  else
243  aero_data_n_spec = -1
244  end if
245 
246  end function aero_data_n_spec
247 
248 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
249 
250  !> Return the number of aerosol sources, or -1 if uninitialized.
251  elemental integer function aero_data_n_source(aero_data)
252 
253  !> Aero data structure.
254  type(aero_data_t), intent(in) :: aero_data
255 
256  if (allocated(aero_data%source_name)) then
257  aero_data_n_source = size(aero_data%source_name)
258  else
259  aero_data_n_source = -1
260  end if
261 
262  end function aero_data_n_source
263 
264 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
265 
266  !> Returns the number of the species in aero_data with the given name, or
267  !> returns 0 if there is no such species.
268  integer function aero_data_spec_by_name(aero_data, name)
269 
270  !> Aero_data data.
271  type(aero_data_t), intent(in) :: aero_data
272  !> Name of species to find.
273  character(len=*), intent(in) :: name
274 
275  integer i
276  logical found
277 
278  found = .false.
279  do i = 1,aero_data_n_spec(aero_data)
280  if (name == aero_data%name(i)) then
281  found = .true.
282  exit
283  end if
284  end do
285  if (found) then
287  else
289  end if
290 
291  end function aero_data_spec_by_name
292 
293 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
294 
295  !> Returns the number of the source in aero_data with the given name, or
296  !> adds the source if it doesn't exist yet.
297  integer function aero_data_source_by_name(aero_data, name)
298 
299  !> Aero_data data.
300  type(aero_data_t), intent(inout) :: aero_data
301  !> Name of source to find.
302  character(len=*), intent(in) :: name
303 
304  if (.not. allocated(aero_data%source_name)) then
305  aero_data%source_name = [name(1:aero_source_name_len)]
307  return
308  end if
309  aero_data_source_by_name = string_array_find(aero_data%source_name, name)
310  if (aero_data_source_by_name > 0) return
311  aero_data%source_name = [aero_data%source_name, &
312  name(1:aero_source_name_len)]
313  aero_data_source_by_name = size(aero_data%source_name)
314 
315  end function aero_data_source_by_name
316 
317 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
318 
319  !> Fills in aero_data%%i_water.
320  subroutine aero_data_set_water_index(aero_data)
321 
322  !> Aero_data data.
323  type(aero_data_t), intent(inout) :: aero_data
324 
325  integer :: i
326 
327  aero_data%i_water = string_array_find(aero_data%name, "H2O")
328 
329  end subroutine aero_data_set_water_index
330 
331 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
332 
333  !> Fills in aero_data%%mosaic_index.
334  subroutine aero_data_set_mosaic_map(aero_data)
335 
336  !> Aero_data data.
337  type(aero_data_t), intent(inout) :: aero_data
338 
339  integer, parameter :: n_mosaic_spec = 19
340  character(AERO_NAME_LEN), parameter, dimension(n_mosaic_spec) :: &
341  mosaic_spec_name = [ &
342  "SO4 ", "NO3 ", "Cl ", "NH4 ", "MSA ", "ARO1 ", &
343  "ARO2 ", "ALK1 ", "OLE1 ", "API1 ", "API2 ", "LIM1 ", &
344  "LIM2 ", "CO3 ", "Na ", "Ca ", "OIN ", "OC ", &
345  "BC "]
346 
347  integer :: i_spec, i_mosaic_spec, i
348 
349  aero_data%mosaic_index = 0
350  do i_spec = 1,aero_data_n_spec(aero_data)
351  aero_data%mosaic_index(i_spec) = string_array_find(mosaic_spec_name, &
352  aero_data%name(i_spec))
353  end do
354 
355  end subroutine aero_data_set_mosaic_map
356 
357 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
358 
359  !> Get the index on the CAMP state array for a specified species and
360  !! computational particle
361  integer function camp_spec_id(aero_data, i_part, i_spec)
362 
363  !> Aerosol data.
364  class(aero_data_t), intent(in) :: aero_data
365  !> Computational particle index (1...aero_state_t%n_part).
366  integer, intent(in) :: i_part
367  !> Aerosol species index in aero_particle_t%vol(:) array.
368  integer, intent(in) :: i_spec
369 
370 #ifdef PMC_USE_CAMP
371  call assert(106669451, allocated(aero_data%camp_particle_spec_id))
372  call assert(278731889, aero_data%camp_particle_state_size .ge. 0)
373  camp_spec_id = (i_part - 1) * aero_data%camp_particle_state_size + &
374  aero_data%camp_particle_spec_id(i_spec)
375 #else
376  camp_spec_id = 0
377 #endif
378 
379  end function camp_spec_id
380 
381 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
382 
383  !> Read aero_data specification from a spec file.
384  subroutine spec_file_read_aero_data(file, aero_data)
385 
386  !> Spec file to read data from.
387  type(spec_file_t), intent(inout) :: file
388  !> Aero_data data.
389  type(aero_data_t), intent(inout) :: aero_data
390 
391  integer :: n_species, species, i
392  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: species_name(:)
393  real(kind=dp), allocatable :: species_data(:,:)
394 
395  !> \page input_format_aero_data Input File Format: Aerosol Material Data
396  !!
397  !! A aerosol material data file must consist of one line per
398  !! aerosol species, with each line having:
399  !! - species name (string)
400  !! - density (real, unit kg/m^3)
401  !! - ions per fully dissociated molecule (integer) - used to
402  !! compute kappa value if the corresponding kappa value is
403  !! zero
404  !! - molecular weight (real, unit kg/mol)
405  !! - kappa hygroscopicity parameter (real, dimensionless) - if
406  !! zero, then inferred from the ions value
407  !!
408  !! This specifies both which species are to be recognized as
409  !! aerosol consituents, as well as their physical properties. For
410  !! example, an aerosol material data file could contain:
411  !! <pre>
412  !! # species dens (kg/m^3) ions (1) molec wght (kg/mole) kappa (1)
413  !! SO4 1800 0 96e-3 0.65
414  !! NO3 1800 0 62e-3 0.65
415  !! Cl 2200 1 35.5e-3 0
416  !! NH4 1800 0 18e-3 0.65
417  !! </pre>
418  !!
419  !! Note that it is an error to specify a non-zero number of ions
420  !! and a non-zero kappa value for a species. If both values are
421  !! zero then that species has zero hygroscopicity parameter. If
422  !! exactly one of kappa or ions is non-zero then the non-zero
423  !! value is used and the zero value is ignored.
424  !!
425  !! See also:
426  !! - \ref spec_file_format --- the input file text format
427  !! - \ref output_format_aero_data --- the corresponding output format
428 
429  call spec_file_read_real_named_array(file, 0, species_name, species_data)
430 
431  ! check the data size
432  n_species = size(species_data, 1)
433  if (.not. ((size(species_data, 2) == 4) .or. (n_species == 0))) then
434  call die_msg(428926381, 'each line in ' // trim(file%name) &
435  // ' should contain exactly 5 values')
436  end if
437 
438  ! allocate and copy over the data
439  call ensure_string_array_size(aero_data%name, n_species)
440  call ensure_integer_array_size(aero_data%mosaic_index, n_species)
441  call ensure_real_array_size(aero_data%density, n_species)
442  call ensure_integer_array_size(aero_data%num_ions, n_species)
443  call ensure_real_array_size(aero_data%molec_weight, n_species)
444  call ensure_real_array_size(aero_data%kappa, n_species)
445  do i = 1,n_species
446  aero_data%name(i) = species_name(i)(1:aero_name_len)
447  aero_data%density(i) = species_data(i,1)
448  aero_data%num_ions(i) = nint(species_data(i,2))
449  aero_data%molec_weight(i) = species_data(i,3)
450  aero_data%kappa(i) = species_data(i,4)
451  call assert_msg(232362742, &
452  (aero_data%num_ions(i) == 0) .or. (aero_data%kappa(i) == 0d0), &
453  "ions and kappa both non-zero for species " &
454  // trim(aero_data%name(i)) // " in " // trim(file%name))
455  if (species_name(i) == "H2O") then
456  aero_data%i_water = i
457  call warn_assert_msg(945800387, &
458  aero_data%density(i) == const%water_density, &
459  "input H2O density not equal to const%water_density (" &
460  // trim(real_to_string(aero_data%density(i))) // " /= " &
461  // trim(real_to_string(const%water_density)) // ")")
462  call warn_assert_msg(233517437, &
463  aero_data%molec_weight(i) == const%water_molec_weight, &
464  "input H2O molec_weight not equal " &
465  // "to const%water_molec_weight (" &
466  // trim(real_to_string(aero_data%molec_weight(i))) // " /= " &
467  // trim(real_to_string(const%water_molec_weight)) // ")")
468  end if
469  end do
470  call aero_data_set_water_index(aero_data)
471  call aero_data_set_mosaic_map(aero_data)
472 
473  end subroutine spec_file_read_aero_data
474 
475 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
476 
477  !> Read a list of species from the given file with the given name.
478  subroutine spec_file_read_species_list(file, name, aero_data, species_list)
479 
480  !> Spec file.
481  type(spec_file_t), intent(inout) :: file
482  !> Name of line.
483  character(len=*), intent(in) :: name
484  !> Aero_data data.
485  type(aero_data_t), intent(in) :: aero_data
486  !> List of species numbers.
487  integer, allocatable :: species_list(:)
488 
489  type(spec_line_t) :: line
490  integer :: i, spec
491 
492  call spec_file_read_line_no_eof(file, line)
493  call spec_file_check_line_name(file, line, name)
494  call ensure_integer_array_size(species_list, size(line%data))
495  do i = 1,size(line%data)
496  spec = aero_data_spec_by_name(aero_data, line%data(i))
497  if (spec == 0) then
498  call spec_file_die_msg(964771833, file, &
499  'unknown species: ' // trim(line%data(i)))
500  end if
501  species_list(i) = spec
502  end do
503 
504  end subroutine spec_file_read_species_list
505 
506 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
507 
508  !> Determines the number of bytes required to pack the given value.
509  integer function pmc_mpi_pack_size_aero_data(val)
510 
511  !> Value to pack.
512  type(aero_data_t), intent(in) :: val
513 
515  pmc_mpi_pack_size_integer(val%i_water) &
516  + pmc_mpi_pack_size_string_array(val%name) &
517  + pmc_mpi_pack_size_integer_array(val%mosaic_index) &
518  + pmc_mpi_pack_size_real_array(val%density) &
519  + pmc_mpi_pack_size_integer_array(val%num_ions) &
520  + pmc_mpi_pack_size_real_array(val%molec_weight) &
521  + pmc_mpi_pack_size_real_array(val%kappa) &
522  + pmc_mpi_pack_size_string_array(val%source_name) &
523  + pmc_mpi_pack_size_fractal(val%fractal)
524 
525  end function pmc_mpi_pack_size_aero_data
526 
527 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
528 
529  !> Packs the given value into the buffer, advancing position.
530  subroutine pmc_mpi_pack_aero_data(buffer, position, val)
531 
532  !> Memory buffer.
533  character, intent(inout) :: buffer(:)
534  !> Current buffer position.
535  integer, intent(inout) :: position
536  !> Value to pack.
537  type(aero_data_t), intent(in) :: val
538 
539 #ifdef PMC_USE_MPI
540  integer :: prev_position
541 
542  prev_position = position
543  call pmc_mpi_pack_integer(buffer, position, val%i_water)
544  call pmc_mpi_pack_string_array(buffer, position, val%name)
545  call pmc_mpi_pack_integer_array(buffer, position, val%mosaic_index)
546  call pmc_mpi_pack_real_array(buffer, position, val%density)
547  call pmc_mpi_pack_integer_array(buffer, position, val%num_ions)
548  call pmc_mpi_pack_real_array(buffer, position, val%molec_weight)
549  call pmc_mpi_pack_real_array(buffer, position, val%kappa)
550  call pmc_mpi_pack_string_array(buffer, position, val%source_name)
551  call pmc_mpi_pack_fractal(buffer, position, val%fractal)
552  call assert(183834856, &
553  position - prev_position <= pmc_mpi_pack_size_aero_data(val))
554 #endif
555 
556  end subroutine pmc_mpi_pack_aero_data
557 
558 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
559 
560  !> Unpacks the given value from the buffer, advancing position.
561  subroutine pmc_mpi_unpack_aero_data(buffer, position, val)
562 
563  !> Memory buffer.
564  character, intent(inout) :: buffer(:)
565  !> Current buffer position.
566  integer, intent(inout) :: position
567  !> Value to pack.
568  type(aero_data_t), intent(inout) :: val
569 
570 #ifdef PMC_USE_MPI
571  integer :: prev_position
572 
573  prev_position = position
574  call pmc_mpi_unpack_integer(buffer, position, val%i_water)
575  call pmc_mpi_unpack_string_array(buffer, position, val%name)
576  call pmc_mpi_unpack_integer_array(buffer, position, val%mosaic_index)
577  call pmc_mpi_unpack_real_array(buffer, position, val%density)
578  call pmc_mpi_unpack_integer_array(buffer, position, val%num_ions)
579  call pmc_mpi_unpack_real_array(buffer, position, val%molec_weight)
580  call pmc_mpi_unpack_real_array(buffer, position, val%kappa)
581  call pmc_mpi_unpack_string_array(buffer, position, val%source_name)
582  call pmc_mpi_unpack_fractal(buffer, position, val%fractal)
583  call assert(188522823, &
584  position - prev_position <= pmc_mpi_pack_size_aero_data(val))
585 #endif
586 
587  end subroutine pmc_mpi_unpack_aero_data
588 
589 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
590 
591  !> Write the aero species dimension to the given NetCDF file if it
592  !> is not already present and in any case return the associated
593  !> dimid.
594  subroutine aero_data_netcdf_dim_aero_species(aero_data, ncid, &
595  dimid_aero_species)
596 
597  !> Aero_data structure.
598  type(aero_data_t), intent(in) :: aero_data
599  !> NetCDF file ID, in data mode.
600  integer, intent(in) :: ncid
601  !> Dimid of the species dimension.
602  integer, intent(out) :: dimid_aero_species
603 
604  integer :: status, i_spec
605  integer :: varid_aero_species
606  integer :: aero_species_centers(aero_data_n_spec(aero_data))
607  character(len=(AERO_NAME_LEN * aero_data_n_spec(aero_data))) :: &
608  aero_species_names
609 
610  ! try to get the dimension ID
611  status = nf90_inq_dimid(ncid, "aero_species", dimid_aero_species)
612  if (status == nf90_noerr) return
613  if (status /= nf90_ebaddim) call pmc_nc_check(status)
614 
615  ! dimension not defined, so define now define it
616  call pmc_nc_check(nf90_redef(ncid))
617 
618  call pmc_nc_check(nf90_def_dim(ncid, "aero_species", &
619  aero_data_n_spec(aero_data), dimid_aero_species))
620  aero_species_names = ""
621  do i_spec = 1,aero_data_n_spec(aero_data)
622  aero_species_names((len_trim(aero_species_names) + 1):) &
623  = trim(aero_data%name(i_spec))
624  if (i_spec < aero_data_n_spec(aero_data)) then
625  aero_species_names((len_trim(aero_species_names) + 1):) = ","
626  end if
627  end do
628  call pmc_nc_check(nf90_def_var(ncid, "aero_species", nf90_int, &
629  dimid_aero_species, varid_aero_species))
630  call pmc_nc_check(nf90_put_att(ncid, varid_aero_species, "names", &
631  aero_species_names))
632  call pmc_nc_check(nf90_put_att(ncid, varid_aero_species, "description", &
633  "dummy dimension variable (no useful value) - read species names " &
634  // "as comma-separated values from the 'names' attribute"))
635 
636  call pmc_nc_check(nf90_enddef(ncid))
637 
638  do i_spec = 1,aero_data_n_spec(aero_data)
639  aero_species_centers(i_spec) = i_spec
640  end do
641  call pmc_nc_check(nf90_put_var(ncid, varid_aero_species, &
642  aero_species_centers))
643 
644  end subroutine aero_data_netcdf_dim_aero_species
645 
646 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
647 
648  !> Write the aero source dimension to the given NetCDF file if it
649  !> is not already present and in any case return the associated
650  !> dimid.
651  subroutine aero_data_netcdf_dim_aero_source(aero_data, ncid, &
652  dimid_aero_source)
653 
654  !> Aero_data structure.
655  type(aero_data_t), intent(in) :: aero_data
656  !> NetCDF file ID, in data mode.
657  integer, intent(in) :: ncid
658  !> Dimid of the source dimension.
659  integer, intent(out) :: dimid_aero_source
660 
661  integer :: status, i_source
662  integer :: varid_aero_source
663  integer :: aero_source_centers(aero_data_n_source(aero_data))
664  character(len=(AERO_SOURCE_NAME_LEN * aero_data_n_source(aero_data))) &
665  :: aero_source_names
666 
667  ! try to get the dimension ID
668  status = nf90_inq_dimid(ncid, "aero_source", dimid_aero_source)
669  if (status == nf90_noerr) return
670  if (status /= nf90_ebaddim) call pmc_nc_check(status)
671 
672  ! dimension not defined, so define now define it
673  call pmc_nc_check(nf90_redef(ncid))
674 
675  call pmc_nc_check(nf90_def_dim(ncid, "aero_source", &
676  aero_data_n_source(aero_data), dimid_aero_source))
677  aero_source_names = ""
678  do i_source = 1,aero_data_n_source(aero_data)
679  aero_source_names((len_trim(aero_source_names) + 1):) &
680  = trim(aero_data%source_name(i_source))
681  if (i_source < aero_data_n_source(aero_data)) then
682  aero_source_names((len_trim(aero_source_names) + 1):) = ","
683  end if
684  end do
685  call pmc_nc_check(nf90_def_var(ncid, "aero_source", nf90_int, &
686  dimid_aero_source, varid_aero_source))
687  call pmc_nc_check(nf90_put_att(ncid, varid_aero_source, "names", &
688  aero_source_names))
689  call pmc_nc_check(nf90_put_att(ncid, varid_aero_source, "description", &
690  "dummy dimension variable (no useful value) - read source names " &
691  // "as comma-separated values from the 'names' attribute"))
692 
693  call pmc_nc_check(nf90_enddef(ncid))
694 
695  do i_source = 1,aero_data_n_source(aero_data)
696  aero_source_centers(i_source) = i_source
697  end do
698  call pmc_nc_check(nf90_put_var(ncid, varid_aero_source, &
699  aero_source_centers))
700 
701  end subroutine aero_data_netcdf_dim_aero_source
702 
703 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
704 
705  !> Write full state.
706  subroutine aero_data_output_netcdf(aero_data, ncid)
707 
708  !> Aero_data to write.
709  type(aero_data_t), intent(in) :: aero_data
710  !> NetCDF file ID, in data mode.
711  integer, intent(in) :: ncid
712 
713  integer :: dimid_aero_species, dimid_aero_source
714 
715  !> \page output_format_aero_data Output File Format: Aerosol Material Data
716  !!
717  !! The aerosol material data NetCDF dimensions are:
718  !! - \b aero_species: number of aerosol species
719  !! - \b aero_source: number of aerosol sources
720  !!
721  !! The aerosol material data NetCDF variables are:
722  !! - \b aero_species (dim \c aero_species): dummy dimension variable
723  !! (no useful value) - read species names as comma-separated values
724  !! from the 'names' attribute
725  !! - \b aero_source (dim \c aero_source): dummy dimension variable
726  !! (no useful value) - read source names as comma-separated values
727  !! from the 'names' attribute
728  !! - \b aero_mosaic_index (dim \c aero_species): indices of species
729  !! in MOSAIC
730  !! - \b aero_density (unit kg/m^3, dim \c aero_species): densities
731  !! of aerosol species
732  !! - \b aero_num_ions (dim \c aero_species): number of ions produced
733  !! when one molecule of each species fully dissociates in water
734  !! - \b aero_molec_weight (unit kg/mol, dim \c aero_species): molecular
735  !! weights of aerosol species
736  !! - \b aero_kappa (unit kg/mol, dim \c aero_species): hygroscopicity
737  !! parameters of aerosol species
738  !! - \b fractal parameters, see \ref output_format_fractal
739  !!
740  !! See also:
741  !! - \ref input_format_aero_data --- the corresponding input format
742 
743  call aero_data_netcdf_dim_aero_species(aero_data, ncid, &
744  dimid_aero_species)
745  call aero_data_netcdf_dim_aero_source(aero_data, ncid, &
746  dimid_aero_source)
747 
748  call pmc_nc_write_integer_1d(ncid, aero_data%mosaic_index, &
749  "aero_mosaic_index", (/ dimid_aero_species /), &
750  long_name="MOSAIC indices of aerosol species")
751  call pmc_nc_write_real_1d(ncid, aero_data%density, &
752  "aero_density", (/ dimid_aero_species /), unit="kg/m^3", &
753  long_name="densities of aerosol species")
754  call pmc_nc_write_integer_1d(ncid, aero_data%num_ions, &
755  "aero_num_ions", (/ dimid_aero_species /), &
756  long_name="number of ions after dissociation of aerosol species")
757  call pmc_nc_write_real_1d(ncid, aero_data%molec_weight, &
758  "aero_molec_weight", (/ dimid_aero_species /), unit="kg/mol", &
759  long_name="molecular weights of aerosol species")
760  call pmc_nc_write_real_1d(ncid, aero_data%kappa, &
761  "aero_kappa", (/ dimid_aero_species /), unit="1", &
762  long_name="hygroscopicity parameters (kappas) of aerosol species")
763  call fractal_output_netcdf(aero_data%fractal, ncid)
764 
765  end subroutine aero_data_output_netcdf
766 
767 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
768 
769  !> Read full state.
770  subroutine aero_data_input_netcdf(aero_data, ncid)
771 
772  !> Aero_data to read.
773  type(aero_data_t), intent(inout) :: aero_data
774  !> NetCDF file ID, in data mode.
775  integer, intent(in) :: ncid
776 
777  integer, parameter :: MAX_SPECIES = 1000
778  integer, parameter :: MAX_SOURCES = 1000
779 
780  character(len=1000) :: name
781  integer :: dimid_aero_species, n_spec, varid_aero_species, i_spec, i
782  integer :: dimid_aero_source, n_source, varid_aero_source, i_source
783  character(len=((AERO_NAME_LEN + 2) * MAX_SPECIES)) :: aero_species_names
784  character(len=:), allocatable :: aero_source_names
785 
786  call pmc_nc_check(nf90_inq_dimid(ncid, "aero_species", &
787  dimid_aero_species))
788  call pmc_nc_check(nf90_inquire_dimension(ncid, &
789  dimid_aero_species, name, n_spec))
790  call assert(141013948, n_spec < max_species)
791 
792  call pmc_nc_check(nf90_inq_dimid(ncid, "aero_source", &
793  dimid_aero_source))
794  call pmc_nc_check(nf90_inquire_dimension(ncid, &
795  dimid_aero_source, name, n_source))
796  call assert(739238793, n_source < max_sources)
797 
798  call pmc_nc_read_integer_1d(ncid, aero_data%mosaic_index, &
799  "aero_mosaic_index")
800  call pmc_nc_read_real_1d(ncid, aero_data%density, "aero_density")
801  call pmc_nc_read_integer_1d(ncid, aero_data%num_ions, "aero_num_ions")
802  call pmc_nc_read_real_1d(ncid, aero_data%molec_weight, "aero_molec_weight")
803  call pmc_nc_read_real_1d(ncid, aero_data%kappa, "aero_kappa")
804 
805  call pmc_nc_check(nf90_inq_varid(ncid, "aero_species", &
806  varid_aero_species))
807  call pmc_nc_check(nf90_get_att(ncid, varid_aero_species, "names", &
808  aero_species_names))
809  ! aero_species_names are comma-separated, so unpack them
810  call ensure_string_array_size(aero_data%name, n_spec)
811  do i_spec = 1,aero_data_n_spec(aero_data)
812  i = 1
813  do while ((aero_species_names(i:i) /= " ") &
814  .and. (aero_species_names(i:i) /= ","))
815  i = i + 1
816  end do
817  call assert(852937292, i > 1)
818  aero_data%name(i_spec) = aero_species_names(1:(i-1))
819  aero_species_names = aero_species_names((i+1):)
820  end do
821  call assert(729138192, aero_species_names == "")
822 
823  call pmc_nc_check(nf90_inq_varid(ncid, "aero_source", &
824  varid_aero_source))
825  allocate(character(len=((AERO_SOURCE_NAME_LEN + 2) * MAX_SPECIES)) &
826  :: aero_source_names)
827  call pmc_nc_check(nf90_get_att(ncid, varid_aero_source, "names", &
828  aero_source_names))
829  ! aero_source_names are comma-separated, so unpack them
830  call ensure_string_array_size(aero_data%source_name, n_source)
831  do i_source = 1,aero_data_n_source(aero_data)
832  i = 1
833  do while ((aero_source_names(i:i) /= " ") &
834  .and. (aero_source_names(i:i) /= ","))
835  i = i + 1
836  end do
837  call assert(840982478, i > 1)
838  aero_data%source_name(i_source) = aero_source_names(1:(i-1))
839  aero_source_names = aero_source_names((i+1):)
840  end do
841  call assert(377166446, aero_source_names == "")
842 
843  call aero_data_set_water_index(aero_data)
844 
845  call fractal_input_netcdf(aero_data%fractal, ncid)
846 
847  end subroutine aero_data_input_netcdf
848 
849 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
850 
851 
852 #ifdef PMC_USE_CAMP
853  !> Initialize the aero_data_t variable with camp chem data
854  subroutine aero_data_initialize(aero_data, camp_core)
855 
856  !> Aerosol data.
857  class(aero_data_t), intent(inout) :: aero_data
858  !> CAMP core.
859  type(camp_core_t), intent(in) :: camp_core
860 
861  character(len=:), allocatable :: rep_name, prop_name, str_val
862  type(string_t), allocatable :: spec_names(:), tmp_spec_names(:)
863  integer :: num_spec, i_spec, spec_type
864  type(chem_spec_data_t), pointer :: chem_spec_data
865  type(property_t), pointer :: property_set
866 
867  rep_name = "PartMC single particle"
868  if (.not. camp_core%get_aero_rep(rep_name, aero_data%aero_rep_ptr)) then
869  call die_msg(418509983, "Missing 'PartMC single particle' aerosol "// &
870  "representation.")
871  end if
872 
873  call assert_msg(935419266, camp_core%get_chem_spec_data(chem_spec_data), &
874  "No chemical species data in camp_core.")
875 
876  ! Only include real aerosol species (no activity coefficients)
877  spec_names = aero_data%aero_rep_ptr%unique_names()
878  allocate(tmp_spec_names(size(spec_names)))
879  num_spec = 0
880  do i_spec = 1, size(spec_names)
881  call assert(496388827, chem_spec_data%get_type( &
882  aero_data%aero_rep_ptr%spec_name(spec_names(i_spec)%string), &
883  spec_type))
884  if (spec_type /= chem_spec_variable .and. &
885  spec_type /= chem_spec_constant .and. &
886  spec_type /= chem_spec_pssa) cycle
887  if (spec_names(i_spec)%string(1:3) /= "P1.") exit
888  num_spec = num_spec + 1
889  tmp_spec_names(num_spec)%string = spec_names(i_spec)%string(4:)
890  end do
891  deallocate(spec_names)
892  allocate(spec_names(num_spec))
893  spec_names(:) = tmp_spec_names(1:num_spec)
894  deallocate(tmp_spec_names)
895 
896  allocate(aero_data%name(num_spec))
897  allocate(aero_data%mosaic_index(num_spec))
898  allocate(aero_data%density(num_spec))
899  allocate(aero_data%num_ions(num_spec))
900  allocate(aero_data%molec_weight(num_spec))
901  allocate(aero_data%kappa(num_spec))
902  allocate(aero_data%camp_particle_spec_id(num_spec))
903 
904  ! Assume no aerosol water
905  aero_data%i_water = 0
906 
907  do i_spec = 1, num_spec
908  aero_data%name(i_spec) = spec_names(i_spec)%string
909  if (.not. chem_spec_data%get_property_set( &
910  aero_data%aero_rep_ptr%spec_name("P1." &
911  // spec_names(i_spec)%string), property_set)) then
912  call die_msg(934844845, "Missing property set for aerosol species " &
913  // spec_names(i_spec)%string)
914  end if
915  prop_name = "density [kg m-3]"
916  if (.not. property_set%get_real(prop_name, &
917  aero_data%density(i_spec))) then
918  call die_msg(547508215, "Missing density for aerosol species " &
919  // spec_names(i_spec)%string)
920  end if
921  prop_name = "num_ions"
922  if (.not. property_set%get_int(prop_name, &
923  aero_data%num_ions(i_spec))) then
924  call die_msg(324777059, "Missing num_ions for aerosol species " &
925  // spec_names(i_spec)%string)
926  end if
927  prop_name = "molecular weight [kg mol-1]"
928  if (.not. property_set%get_real(prop_name, &
929  aero_data%molec_weight(i_spec))) then
930  call die_msg(549413749, "Missing molec_weight for aerosol species " &
931  // spec_names(i_spec)%string)
932  end if
933  prop_name = "kappa"
934  if (.not. property_set%get_real(prop_name, &
935  aero_data%kappa(i_spec))) then
936  call die_msg(944207343, "Missing kappa for aerosol species " &
937  // spec_names(i_spec)%string)
938  end if
939  prop_name = "PartMC name"
940  if (property_set%get_string(prop_name, str_val)) then
941  if (str_val == "H2O") then
942  call assert_msg(227489086, aero_data%i_water == 0, &
943  "Multiple aerosol water species")
944  aero_data%i_water = i_spec
945  end if
946  end if
947  aero_data%camp_particle_spec_id(i_spec) = &
948  aero_data%aero_rep_ptr%spec_state_id("P1." &
949  // spec_names(i_spec)%string)
950  end do
951 
952  select type( aero_rep => aero_data%aero_rep_ptr)
953  type is(aero_rep_single_particle_t)
954 
955  ! Get the number of elements per-particle on the CAMP state array
956  aero_data%camp_particle_state_size = aero_rep%per_particle_size()
957 
958  class default
959  call die_msg(281737350, "Wrong aerosol representation type")
960  end select
961 
962  end subroutine aero_data_initialize
963 #endif
964 
965 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
966 
967 end module pmc_aero_data
pmc_mpi::pmc_mpi_pack_string_array
subroutine pmc_mpi_pack_string_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:752
pmc_fractal::fractal_vol2rad
real(kind=dp) elemental function fractal_vol2rad(fractal, v)
Convert mass-equivalent volume (m^3) to geometric radius (m).
Definition: fractal.F90:107
pmc_netcdf::pmc_nc_read_integer_1d
subroutine pmc_nc_read_integer_1d(ncid, var, name, must_be_present)
Read a simple integer array from a NetCDF file.
Definition: netcdf.F90:260
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:236
pmc_netcdf::pmc_nc_write_real_1d
subroutine pmc_nc_write_real_1d(ncid, var, name, dimids, dim_name, unit, long_name, standard_name, description)
Write a simple real array to a NetCDF file.
Definition: netcdf.F90:536
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_mpi::pmc_mpi_unpack_string_array
subroutine pmc_mpi_unpack_string_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1010
pmc_aero_data::aero_data_vol_to_num_of_monomers
real(kind=dp) elemental function aero_data_vol_to_num_of_monomers(aero_data, v)
Convert mass-equivalent volume (m^3) to number of monomers in a fractal particle cluster.
Definition: aero_data.F90:155
pmc_aero_data::pmc_mpi_unpack_aero_data
subroutine pmc_mpi_unpack_aero_data(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_data.F90:562
pmc_aero_data::aero_data_vol_to_mobility_rad
real(kind=dp) function aero_data_vol_to_mobility_rad(aero_data, v, temp, pressure)
Convert mass-equivalent volume (m^3) to mobility equivalent radius (m).
Definition: aero_data.F90:174
pmc_mpi::pmc_mpi_unpack_integer_array
subroutine pmc_mpi_unpack_integer_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:946
pmc_fractal::pmc_mpi_pack_size_fractal
integer function pmc_mpi_pack_size_fractal(val)
Determines the number of bytes required to pack the given value.
Definition: fractal.F90:490
pmc_mpi::pmc_mpi_pack_integer_array
subroutine pmc_mpi_pack_integer_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:688
pmc_aero_data::aero_data_source_by_name
integer function aero_data_source_by_name(aero_data, name)
Returns the number of the source in aero_data with the given name, or adds the source if it doesn't e...
Definition: aero_data.F90:298
pmc_netcdf::pmc_nc_read_real_1d
subroutine pmc_nc_read_real_1d(ncid, var, name, must_be_present)
Read a simple real array from a NetCDF file.
Definition: netcdf.F90:219
pmc_netcdf::pmc_nc_write_integer_1d
subroutine pmc_nc_write_integer_1d(ncid, var, name, dimids, dim_name, unit, long_name, standard_name, description)
Write a simple integer array to a NetCDF file.
Definition: netcdf.F90:585
pmc_aero_data::camp_spec_id
integer function camp_spec_id(aero_data, i_part, i_spec)
Get the index on the CAMP state array for a specified species and computational particle.
Definition: aero_data.F90:362
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:59
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:134
pmc_netcdf
Wrapper functions for NetCDF. These all take a NetCDF ncid in data mode and return with it again in d...
Definition: netcdf.F90:11
pmc_aero_data::aero_name_len
integer, parameter aero_name_len
Definition: aero_data.F90:28
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:1041
pmc_aero_data::aero_data_n_source
elemental integer function aero_data_n_source(aero_data)
Return the number of aerosol sources, or -1 if uninitialized.
Definition: aero_data.F90:252
pmc_fractal::pmc_mpi_unpack_fractal
subroutine pmc_mpi_unpack_fractal(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: fractal.F90:528
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
pmc_aero_data::aero_data_initialize
subroutine aero_data_initialize(aero_data, camp_core)
Initialize the aero_data_t variable with camp chem data.
Definition: aero_data.F90:855
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:477
pmc_fractal::fractal_input_netcdf
subroutine fractal_input_netcdf(fractal, ncid)
Read full state.
Definition: fractal.F90:633
pmc_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_fractal::fractal_mobility_rad_to_vol
real(kind=dp) function fractal_mobility_rad_to_vol(fractal, mobility_rad, temp, pressure)
Convert mobility equivalent radius (m) to mass-equivalent volume (m^3).
Definition: fractal.F90:468
pmc_aero_data::aero_data_input_netcdf
subroutine aero_data_input_netcdf(aero_data, ncid)
Read full state.
Definition: aero_data.F90:771
pmc_aero_data::aero_data_set_water_index
subroutine aero_data_set_water_index(aero_data)
Fills in aero_data%i_water.
Definition: aero_data.F90:321
pmc_util::real_to_string
character(len=pmc_util_convert_string_len) function real_to_string(val)
Convert a real to a string format.
Definition: util.F90:782
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:77
pmc_util::diam2rad
real(kind=dp) elemental function diam2rad(d)
Convert diameter (m) to radius (m).
Definition: util.F90:277
pmc_fractal::fractal_vol2diam
real(kind=dp) elemental function fractal_vol2diam(fractal, v)
Convert mass-equivalent volume (m^3) to geometric diameter (m).
Definition: fractal.F90:128
pmc_aero_data::aero_data_mobility_rad_to_geometric_rad
real(kind=dp) function aero_data_mobility_rad_to_geometric_rad(aero_data, mobility_rad, temp, pressure)
Convert mobility equivalent radius (m) to geometric radius (m^3).
Definition: aero_data.F90:216
pmc_fractal::fractal_mobility_rad_to_geometric_rad
real(kind=dp) function fractal_mobility_rad_to_geometric_rad(fractal, mobility_rad, temp, pressure)
Convert mobility equivalent radius (m) to geometric radius (m^3).
Definition: fractal.F90:436
pmc_aero_data::aero_data_mobility_rad_to_vol
real(kind=dp) function aero_data_mobility_rad_to_vol(aero_data, mobility_rad, temp, pressure)
Convert mobility equivalent radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:195
pmc_spec_file::spec_file_read_line_no_eof
subroutine spec_file_read_line_no_eof(file, line)
Read a spec_line from the spec_file. This will always succeed or error out, so should only be called ...
Definition: spec_file.F90:298
fractal
type(fractal_t), intent(in) fractal
Definition: fractal.F90:614
pmc_spec_file::spec_file_read_real_named_array
subroutine spec_file_read_real_named_array(file, max_lines, names, vals)
Read an array of named lines with real data. All lines must have the same number of data elements.
Definition: spec_file.F90:650
pmc_util::ensure_string_array_size
subroutine ensure_string_array_size(x, n)
Allocate or reallocate the given array to ensure it is of the given size, without preserving data.
Definition: util.F90:1189
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:92
pmc_fractal::fractal_vol_to_mobility_rad
real(kind=dp) function fractal_vol_to_mobility_rad(fractal, v, temp, pressure)
Convert mass-equivalent volume (m^3) to mobility equivalent radius (m).
Definition: fractal.F90:317
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:49
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:818
pmc_aero_data::aero_data_netcdf_dim_aero_species
subroutine aero_data_netcdf_dim_aero_species(aero_data, ncid, dimid_aero_species)
Write the aero species dimension to the given NetCDF file if it is not already present and in any cas...
Definition: aero_data.F90:596
pmc_fractal::fractal_vol_to_num_of_monomers
real(kind=dp) elemental function fractal_vol_to_num_of_monomers(fractal, v)
Convert mass-equivalent volume (m^3) to number of monomers in a fractal particle cluster.
Definition: fractal.F90:91
pmc_mpi::pmc_mpi_pack_size_integer_array
integer function pmc_mpi_pack_size_integer_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:447
pmc_fractal
The fractal_t structure and associated subroutines.
Definition: fractal.F90:35
pmc_netcdf::pmc_nc_check
subroutine pmc_nc_check(status)
Check the status of a NetCDF function call.
Definition: netcdf.F90:23
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_aero_data::aero_source_name_len
integer, parameter aero_source_name_len
Definition: aero_data.F90:29
pmc_spec_file::spec_file_die_msg
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:74
pmc_util::string_array_find
integer function string_array_find(array, val)
Return the index of the first occurance of the given value in the array, or 0 if it is not present.
Definition: util.F90:1023
pmc_mpi::pmc_mpi_pack_size_string_array
integer function pmc_mpi_pack_size_string_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:507
pmc_aero_data::aero_data_vol2diam
real(kind=dp) elemental function aero_data_vol2diam(aero_data, v)
Convert mass-equivalent volume (m^3) to geometric diameter (m).
Definition: aero_data.F90:107
pmc_fractal::fractal_t
Fractal data.
Definition: fractal.F90:59
pmc_aero_data::aero_data_set_mosaic_map
subroutine aero_data_set_mosaic_map(aero_data)
Fills in aero_data%mosaic_index.
Definition: aero_data.F90:335
pmc_fractal::fractal_rad2vol
real(kind=dp) elemental function fractal_rad2vol(fractal, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: fractal.F90:143
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_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:561
pmc_aero_data::aero_data_rad2vol
real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:122
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_aero_data::pmc_mpi_pack_size_aero_data
integer function pmc_mpi_pack_size_aero_data(val)
Determines the number of bytes required to pack the given value.
Definition: aero_data.F90:510
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:720
pmc_spec_file::spec_file_check_line_name
subroutine spec_file_check_line_name(file, line, name)
Check that the name of the line data is as given.
Definition: spec_file.F90:390
pmc_fractal::pmc_mpi_pack_fractal
subroutine pmc_mpi_pack_fractal(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: fractal.F90:505
pmc_aero_data::spec_file_read_species_list
subroutine spec_file_read_species_list(file, name, aero_data, species_list)
Read a list of species from the given file with the given name.
Definition: aero_data.F90:479
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:978
pmc_aero_data::aero_data_spec_by_name
integer function aero_data_spec_by_name(aero_data, name)
Returns the number of the species in aero_data with the given name, or returns 0 if there is no such ...
Definition: aero_data.F90:269
pmc_aero_data::aero_data_diam2vol
real(kind=dp) elemental function aero_data_diam2vol(aero_data, d)
Convert geometric diameter (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:137
pmc_aero_data::aero_data_netcdf_dim_aero_source
subroutine aero_data_netcdf_dim_aero_source(aero_data, ncid, dimid_aero_source)
Write the aero source dimension to the given NetCDF file if it is not already present and in any case...
Definition: aero_data.F90:653
pmc_aero_data::pmc_mpi_pack_aero_data
subroutine pmc_mpi_pack_aero_data(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_data.F90:531
pmc_util::ensure_integer_array_size
subroutine ensure_integer_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:1115