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