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