PartMC  2.2.1
aero_mode.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West
00002 ! Licensed under the GNU General Public License version 2 or (at your
00003 ! option) any later version. See the file COPYING for details.
00004 
00005 !> \file
00006 !> The pmc_aero_mode module.
00007 
00008 !> The aero_mode_t structure and associated subroutines.
00009 module pmc_aero_mode
00010 
00011   use pmc_bin_grid
00012   use pmc_util
00013   use pmc_constants
00014   use pmc_spec_file
00015   use pmc_aero_data
00016   use pmc_aero_weight
00017   use pmc_mpi
00018   use pmc_rand
00019 #ifdef PMC_USE_MPI
00020   use mpi
00021 #endif
00022 
00023   !> Maximum length of a mode name.
00024   integer, parameter :: AERO_MODE_NAME_LEN = 300
00025   !> Maximum length of a mode type.
00026   integer, parameter :: AERO_MODE_TYPE_LEN = 20
00027 
00028   !> Type code for an undefined or invalid mode.
00029   integer, parameter :: AERO_MODE_TYPE_INVALID    = 0
00030   !> Type code for a log-normal mode.
00031   integer, parameter :: AERO_MODE_TYPE_LOG_NORMAL = 1
00032   !> Type code for an exponential mode.
00033   integer, parameter :: AERO_MODE_TYPE_EXP        = 2
00034   !> Type code for a mono-disperse mode.
00035   integer, parameter :: AERO_MODE_TYPE_MONO       = 3
00036   !> Type code for a sampled mode.
00037   integer, parameter :: AERO_MODE_TYPE_SAMPLED    = 4
00038 
00039   !> An aerosol size distribution mode.
00040   !!
00041   !! Each mode is assumed to be fully internally mixed so that every
00042   !! particle has the same composition. The composition is stored in
00043   !! \c vol_frac, while the other parameters define the size
00044   !! distribution (with \c type defining the type of size distribution
00045   !! function). See \ref input_format_mass_frac for descriptions of
00046   !! the parameters relevant to each mode type.
00047   type aero_mode_t
00048      !> Mode name, used to track particle sources.
00049      character(len=AERO_MODE_NAME_LEN) :: name
00050      !> Mode type (given by module constants).
00051      integer :: type
00052      !> Characteristic radius, with meaning dependent on mode type (m).
00053      real(kind=dp) :: char_radius
00054      !> Log base 10 of geometric standard deviation of radius, (m).
00055      real(kind=dp) :: log10_std_dev_radius
00056      !> Sample bin radii [length <tt>(N + 1)</tt>] (m).
00057      real(kind=dp), pointer :: sample_radius(:)
00058      !> Sample bin number concentrations [length <tt>N</tt>] (m^{-3}).
00059      real(kind=dp), pointer :: sample_num_conc(:)
00060      !> Total number concentration of mode (#/m^3).
00061      real(kind=dp) :: num_conc
00062      !> Species fractions by volume [length \c aero_data%%n_spec] (1).
00063      real(kind=dp), pointer :: vol_frac(:)
00064      !> Source number.
00065      integer :: source
00066   end type aero_mode_t
00067 
00068 contains
00069 
00070 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00071 
00072   !> Return a string representation of a kernel type.
00073   character(len=AERO_MODE_TYPE_LEN) function aero_mode_type_to_string(type)
00074 
00075     !> Aero mode type.
00076     integer, intent(in) :: type
00077    
00078     if (type == AERO_MODE_TYPE_INVALID) then
00079        aero_mode_type_to_string = "invalid"
00080     elseif (type == AERO_MODE_TYPE_LOG_NORMAL) then
00081        aero_mode_type_to_string = "log_normal"
00082     elseif (type == AERO_MODE_TYPE_EXP) then
00083        aero_mode_type_to_string = "exp"
00084     elseif (type == AERO_MODE_TYPE_MONO) then
00085        aero_mode_type_to_string = "mono"
00086     elseif (type == AERO_MODE_TYPE_SAMPLED) then
00087        aero_mode_type_to_string = "sampled"
00088     else
00089        aero_mode_type_to_string = "unknown"
00090     end if
00091 
00092   end function aero_mode_type_to_string
00093 
00094 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00095 
00096   !> Allocates an aero_mode.
00097   subroutine aero_mode_allocate(aero_mode)
00098 
00099     !> Aerosol mode.
00100     type(aero_mode_t), intent(out) :: aero_mode
00101 
00102     allocate(aero_mode%vol_frac(0))
00103     allocate(aero_mode%sample_radius(0))
00104     allocate(aero_mode%sample_num_conc(0))
00105 
00106   end subroutine aero_mode_allocate
00107 
00108 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00109 
00110   !> Allocates an aero_mode of the given size.
00111   subroutine aero_mode_allocate_size(aero_mode, n_spec)
00112 
00113     !> Aerosol mode.
00114     type(aero_mode_t), intent(out) :: aero_mode
00115     !> Number of species.
00116     integer, intent(in) :: n_spec
00117 
00118     allocate(aero_mode%vol_frac(n_spec))
00119     allocate(aero_mode%sample_radius(0))
00120     allocate(aero_mode%sample_num_conc(0))
00121 
00122   end subroutine aero_mode_allocate_size
00123 
00124 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00125 
00126   !> Free all storage.
00127   subroutine aero_mode_deallocate(aero_mode)
00128 
00129     !> Aerosol mode.
00130     type(aero_mode_t), intent(inout) :: aero_mode
00131 
00132     deallocate(aero_mode%vol_frac)
00133     deallocate(aero_mode%sample_radius)
00134     deallocate(aero_mode%sample_num_conc)
00135 
00136   end subroutine aero_mode_deallocate
00137 
00138 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00139 
00140   !> Copy an aero_mode.
00141   subroutine aero_mode_copy(aero_mode_from, aero_mode_to)
00142 
00143     !> Aerosol mode original.
00144     type(aero_mode_t), intent(in) :: aero_mode_from
00145     !> Aerosol mode copy.
00146     type(aero_mode_t), intent(inout) :: aero_mode_to
00147 
00148     call aero_mode_deallocate(aero_mode_to)
00149     call aero_mode_allocate_size(aero_mode_to, size(aero_mode_from%vol_frac))
00150     aero_mode_to%name = aero_mode_from%name
00151     aero_mode_to%type = aero_mode_from%type
00152     aero_mode_to%char_radius = aero_mode_from%char_radius
00153     aero_mode_to%log10_std_dev_radius = aero_mode_from%log10_std_dev_radius
00154     aero_mode_to%num_conc = aero_mode_from%num_conc
00155     aero_mode_to%vol_frac = aero_mode_from%vol_frac
00156     aero_mode_to%source = aero_mode_from%source
00157     deallocate(aero_mode_to%sample_radius)
00158     deallocate(aero_mode_to%sample_num_conc)
00159     allocate(aero_mode_to%sample_radius(size(aero_mode_from%sample_radius)))
00160     allocate(aero_mode_to%sample_num_conc( &
00161          size(aero_mode_from%sample_num_conc)))
00162     aero_mode_to%sample_radius = aero_mode_from%sample_radius
00163     aero_mode_to%sample_num_conc = aero_mode_from%sample_num_conc
00164 
00165   end subroutine aero_mode_copy
00166 
00167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00168 
00169   !> Returns the total number concentration of a mode. (#/m^3)
00170   real(kind=dp) function aero_mode_total_num_conc(aero_mode)
00171 
00172     !> Aerosol mode.
00173     type(aero_mode_t), intent(in) :: aero_mode
00174 
00175     if ((aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) &
00176          .or. (aero_mode%type == AERO_MODE_TYPE_EXP) &
00177        .or. (aero_mode%type == AERO_MODE_TYPE_MONO)) then
00178        aero_mode_total_num_conc = aero_mode%num_conc
00179     elseif (aero_mode%type == AERO_MODE_TYPE_SAMPLED) then
00180        aero_mode_total_num_conc = sum(aero_mode%sample_num_conc)
00181     else
00182        call die_msg(719625922, "unknown aero_mode type: " &
00183             // trim(integer_to_string(aero_mode%type)))
00184     end if
00185 
00186   end function aero_mode_total_num_conc
00187 
00188 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00189 
00190   !> Compute a log-normal distribution.
00191   subroutine num_conc_log_normal(total_num_conc, geom_mean_radius, &
00192        log10_sigma_g, bin_grid, num_conc)
00193 
00194     !> Total number concentration of the mode (m^{-3}).
00195     real(kind=dp), intent(in) :: total_num_conc
00196     !> Geometric mean radius (m).
00197     real(kind=dp), intent(in) :: geom_mean_radius
00198     !> log_10(geom. std. dev.) (1).
00199     real(kind=dp), intent(in) :: log10_sigma_g
00200     !> Bin grid.
00201     type(bin_grid_t), intent(in) :: bin_grid
00202     !> Number concentration (#(ln(r))d(ln(r))).
00203     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00204     
00205     integer :: k
00206     
00207     do k = 1,bin_grid%n_bin
00208        num_conc(k) = total_num_conc / (sqrt(2d0 * const%pi) &
00209             * log10_sigma_g) * dexp(-(dlog10(bin_grid%center_radius(k)) &
00210             - dlog10(geom_mean_radius))**2d0 &
00211             / (2d0 * log10_sigma_g**2d0)) / dlog(10d0)
00212     end do
00213     
00214     ! The formula above was originally for a distribution in
00215     ! log_10(r), while we are using log_e(r) for our bin grid. The
00216     ! division by dlog(10) at the end corrects for this.
00217 
00218     ! Remember that log_e(r) = log_10(r) * log_e(10).
00219     
00220   end subroutine num_conc_log_normal
00221   
00222 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00223 
00224   !> Compute a log-normal distribution in volume.
00225   subroutine vol_conc_log_normal(total_num_conc, geom_mean_radius, &
00226        log10_sigma_g, bin_grid, vol_conc)
00227     
00228     !> Total number concentration of the mode (m^{-3}).
00229     real(kind=dp), intent(in) :: total_num_conc
00230     !> Geometric mean radius (m).
00231     real(kind=dp), intent(in) :: geom_mean_radius
00232     !> log_10(geom. std. dev.) (1).
00233     real(kind=dp), intent(in) :: log10_sigma_g
00234     !> Bin grid.
00235     type(bin_grid_t), intent(in) :: bin_grid
00236     !> Volume concentration (V(ln(r))d(ln(r))).
00237     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00238     
00239     real(kind=dp) :: num_conc(bin_grid%n_bin)
00240 
00241     call num_conc_log_normal(total_num_conc, geom_mean_radius, &
00242          log10_sigma_g, bin_grid, num_conc)
00243     vol_conc = num_conc * rad2vol(bin_grid%center_radius)
00244     
00245   end subroutine vol_conc_log_normal
00246   
00247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00248 
00249   !> Exponential distribution in volume.
00250   !!
00251   !! \f[ n(v) = \frac{1}{\rm mean-vol} \exp(- v / {\rm mean-vol}) \f]
00252   !! Normalized so that sum(num_conc(k) * log_width) = 1
00253   subroutine num_conc_exp(total_num_conc, radius_at_mean_vol, bin_grid, &
00254        num_conc)
00255     
00256     !> Total number concentration of the mode (m^{-3}).
00257     real(kind=dp), intent(in) :: total_num_conc
00258     !> Radius at mean volume (m).
00259     real(kind=dp), intent(in) :: radius_at_mean_vol
00260     !> Bin grid.
00261     type(bin_grid_t), intent(in) :: bin_grid
00262     !> Number concentration (#(ln(r))d(ln(r))).
00263     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00264     
00265     integer :: k
00266     real(kind=dp) :: mean_vol, num_conc_vol
00267     
00268     mean_vol = rad2vol(radius_at_mean_vol)
00269     do k = 1,bin_grid%n_bin
00270        num_conc_vol = total_num_conc / mean_vol &
00271             * exp(-(rad2vol(bin_grid%center_radius(k)) / mean_vol))
00272        call vol_to_lnr(bin_grid%center_radius(k), num_conc_vol, num_conc(k))
00273     end do
00274     
00275   end subroutine num_conc_exp
00276   
00277 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00278 
00279   !> Exponential distribution in volume.
00280   subroutine vol_conc_exp(total_num_conc, radius_at_mean_vol, bin_grid, &
00281        vol_conc)
00282     
00283     !> Total number concentration of the mode (m^{-3}).
00284     real(kind=dp), intent(in) :: total_num_conc
00285     !> Radius at mean volume (m).
00286     real(kind=dp), intent(in) :: radius_at_mean_vol
00287     !> Bin grid.
00288     type(bin_grid_t), intent(in) :: bin_grid
00289     !> Volume concentration (V(ln(r))d(ln(r))).
00290     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00291     
00292     real(kind=dp) :: num_conc(bin_grid%n_bin)
00293 
00294     call num_conc_exp(total_num_conc, radius_at_mean_vol, bin_grid, num_conc)
00295     vol_conc = num_conc * rad2vol(bin_grid%center_radius)
00296     
00297   end subroutine vol_conc_exp
00298   
00299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00300 
00301   !> Mono-disperse distribution.
00302   !> Normalized so that sum(num_conc(k) * log_width) = 1
00303   subroutine num_conc_mono(total_num_conc, radius, bin_grid, num_conc)
00304     
00305     !> Total number concentration of the mode (m^{-3}).
00306     real(kind=dp), intent(in) :: total_num_conc
00307     !> Radius of each particle (m^3).
00308     real(kind=dp), intent(in) :: radius
00309     !> Bin grid.
00310     type(bin_grid_t), intent(in) :: bin_grid
00311     !> Number concentration (#(ln(r))d(ln(r))).
00312     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00313     
00314     integer :: k
00315 
00316     num_conc = 0d0
00317     k = bin_grid_particle_in_bin(bin_grid, radius)
00318     if ((k < 1) .or. (k > bin_grid%n_bin)) then
00319        call warn_msg(825666877, "monodisperse radius outside of bin_grid")
00320     else
00321        num_conc(k) = total_num_conc / bin_grid%log_width
00322     end if
00323     
00324   end subroutine num_conc_mono
00325   
00326 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00327 
00328   !> Mono-disperse distribution in volume.
00329   subroutine vol_conc_mono(total_num_conc, radius, bin_grid, vol_conc)
00330     
00331     !> Total number concentration of the mode (m^{-3}).
00332     real(kind=dp), intent(in) :: total_num_conc
00333     !> Radius of each particle (m^3).
00334     real(kind=dp), intent(in) :: radius
00335     !> Bin grid.
00336     type(bin_grid_t), intent(in) :: bin_grid
00337     !> Volume concentration (V(ln(r))d(ln(r))).
00338     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00339     
00340     integer :: k
00341 
00342     vol_conc = 0d0
00343     k = bin_grid_particle_in_bin(bin_grid, radius)
00344     if ((k < 1) .or. (k > bin_grid%n_bin)) then
00345        call warn_msg(420930707, "monodisperse radius outside of bin_grid")
00346     else
00347        vol_conc(k) = total_num_conc / bin_grid%log_width * rad2vol(radius)
00348     end if
00349     
00350   end subroutine vol_conc_mono
00351   
00352 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00353 
00354   !> Sampled distribution, not normalized.
00355   subroutine num_conc_sampled(sample_radius, sample_num_conc, bin_grid, &
00356        num_conc)
00357     
00358     !> Sampled radius bin edges (m).
00359     real(kind=dp), intent(in) :: sample_radius(:)
00360     !> Sampled number concentrations (m^{-3}).
00361     real(kind=dp), intent(in) :: sample_num_conc(:)
00362     !> Bin grid.
00363     type(bin_grid_t), intent(in) :: bin_grid
00364     !> Number concentration (#(ln(r))d(ln(r))).
00365     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00366 
00367     integer :: i_sample, n_sample, i_lower, i_upper, i_bin
00368     real(kind=dp) :: r_lower, r_upper
00369     real(kind=dp) :: r_bin_lower, r_bin_upper, r1, r2, ratio
00370     
00371     n_sample = size(sample_num_conc)
00372     call assert(188766208, size(sample_radius) == n_sample + 1)
00373     call assert(295384037, n_sample >= 1)
00374 
00375     num_conc = 0d0
00376     do i_sample = 1,n_sample
00377        r_lower = sample_radius(i_sample)
00378        r_upper = sample_radius(i_sample + 1)
00379        i_lower = bin_grid_particle_in_bin(bin_grid, r_lower)
00380        i_upper = bin_grid_particle_in_bin(bin_grid, r_upper)
00381        if (i_upper < 1) cycle
00382        if (i_lower > bin_grid%n_bin) cycle
00383        i_lower = max(1, i_lower)
00384        i_upper = min(bin_grid%n_bin, i_upper)
00385        do i_bin = i_lower,i_upper
00386           r_bin_lower = bin_grid%edge_radius(i_bin)
00387           r_bin_upper = bin_grid%edge_radius(i_bin + 1)
00388           r1 = max(r_lower, r_bin_lower)
00389           r2 = min(r_upper, r_bin_upper)
00390           ratio = (log(r2) - log(r1)) / (log(r_upper) - log(r_lower))
00391           num_conc(i_bin) = num_conc(i_bin) &
00392                + ratio * sample_num_conc(i_sample) / bin_grid%log_width
00393        end do
00394     end do
00395     
00396   end subroutine num_conc_sampled
00397   
00398 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00399 
00400   !> Sampled distribution in volume.
00401   subroutine vol_conc_sampled(sample_radius, sample_num_conc, bin_grid, &
00402        vol_conc)
00403     
00404     !> Sampled radius bin edges (m).
00405     real(kind=dp), intent(in) :: sample_radius(:)
00406     !> Sampled number concentrations (m^{-3}).
00407     real(kind=dp), intent(in) :: sample_num_conc(:)
00408     !> Bin grid.
00409     type(bin_grid_t), intent(in) :: bin_grid
00410     !> Volume concentration (V(ln(r))d(ln(r))).
00411     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00412     
00413     real(kind=dp) :: num_conc(bin_grid%n_bin)
00414 
00415     call num_conc_sampled(sample_radius, sample_num_conc, bin_grid, num_conc)
00416     vol_conc = num_conc * rad2vol(bin_grid%center_radius)
00417     
00418   end subroutine vol_conc_sampled
00419   
00420 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00421 
00422   !> Return the binned number concentration for an aero_mode.
00423   subroutine aero_mode_num_conc(aero_mode, bin_grid, aero_data, &
00424        num_conc)
00425 
00426     !> Aero mode for which to compute number concentration.
00427     type(aero_mode_t), intent(in) :: aero_mode
00428     !> Bin grid.
00429     type(bin_grid_t), intent(in) :: bin_grid
00430     !> Aerosol data.
00431     type(aero_data_t), intent(in) :: aero_data
00432     !> Number concentration (#(ln(r))d(ln(r))).
00433     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00434 
00435     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00436        call num_conc_log_normal(aero_mode%num_conc, aero_mode%char_radius, &
00437             aero_mode%log10_std_dev_radius, bin_grid, num_conc)
00438     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00439        call num_conc_exp(aero_mode%num_conc, aero_mode%char_radius, &
00440             bin_grid, num_conc)
00441     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00442        call num_conc_mono(aero_mode%num_conc, aero_mode%char_radius, &
00443             bin_grid, num_conc)
00444     elseif (aero_mode%type == AERO_MODE_TYPE_SAMPLED) then
00445        call num_conc_sampled(aero_mode%sample_radius, &
00446             aero_mode%sample_num_conc, bin_grid, num_conc)
00447     else
00448        call die_msg(223903246, "unknown aero_mode type: " &
00449             // trim(integer_to_string(aero_mode%type)))
00450     end if
00451 
00452   end subroutine aero_mode_num_conc
00453 
00454 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00455 
00456   !> Return the binned per-species volume concentration for an
00457   !> aero_mode.
00458   subroutine aero_mode_vol_conc(aero_mode, bin_grid, aero_data, &
00459        vol_conc)
00460 
00461     !> Aero mode for which to compute volume concentration.
00462     type(aero_mode_t), intent(in) :: aero_mode
00463     !> Bin grid.
00464     type(bin_grid_t), intent(in) :: bin_grid
00465     !> Aerosol data.
00466     type(aero_data_t), intent(in) :: aero_data
00467     !> Volume concentration (V(ln(r))d(ln(r))).
00468     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin, aero_data%n_spec)
00469 
00470     integer :: i_spec
00471     real(kind=dp) :: vol_conc_total(bin_grid%n_bin)
00472 
00473     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00474        call vol_conc_log_normal(aero_mode%num_conc, aero_mode%char_radius, &
00475             aero_mode%log10_std_dev_radius, bin_grid, vol_conc_total)
00476     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00477        call vol_conc_exp(aero_mode%num_conc, aero_mode%char_radius, &
00478             bin_grid, vol_conc_total)
00479     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00480        call vol_conc_mono(aero_mode%num_conc, aero_mode%char_radius, &
00481             bin_grid, vol_conc_total)
00482     elseif (aero_mode%type == AERO_MODE_TYPE_SAMPLED) then
00483        call vol_conc_sampled(aero_mode%sample_radius, &
00484             aero_mode%sample_num_conc, bin_grid, vol_conc_total)
00485     else
00486        call die_msg(314169653, "Unknown aero_mode type: " &
00487             // trim(integer_to_string(aero_mode%type)))
00488     end if
00489     do i_spec = 1,aero_data%n_spec
00490        vol_conc(:,i_spec) = vol_conc_total &
00491             * aero_mode%vol_frac(i_spec)
00492     end do
00493 
00494   end subroutine aero_mode_vol_conc
00495 
00496 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00497 
00498   !> Compute weighted sampled number concentrations.
00499   subroutine aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
00500        weighted_num_conc)
00501 
00502     !> Aerosol mode.
00503     type(aero_mode_t), intent(in) :: aero_mode
00504     !> Aerosol weight.
00505     type(aero_weight_t), intent(in) :: aero_weight
00506     !> Weighted number concentration
00507     real(kind=dp), intent(out) :: weighted_num_conc(:)
00508 
00509     integer :: i_sample
00510     real(kind=dp) :: x0, x1, xr
00511 
00512     call assert(256667423, aero_mode%type == AERO_MODE_TYPE_SAMPLED)
00513     call assert(878731017, &
00514          size(weighted_num_conc) == size(aero_mode%sample_num_conc))
00515 
00516     if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00517        weighted_num_conc = aero_mode%sample_num_conc
00518     elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00519          .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00520        do i_sample = 1,size(aero_mode%sample_num_conc)
00521           x0 = log(aero_mode%sample_radius(i_sample))
00522           x1 = log(aero_mode%sample_radius(i_sample + 1))
00523           xr = log(aero_weight%ref_radius)
00524           weighted_num_conc(i_sample) = aero_mode%sample_num_conc(i_sample) &
00525                / aero_weight%exponent &
00526                * (exp(- aero_weight%exponent * (x0 - xr)) &
00527                - exp(- aero_weight%exponent * (x1 - xr))) &
00528                / (x1 - x0)
00529        end do
00530     else
00531        call die_msg(576124393, "unknown aero_weight type: " &
00532             // trim(integer_to_string(aero_weight%type)))
00533     end if
00534 
00535   end subroutine aero_mode_weighted_sampled_num_conc
00536 
00537 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00538 
00539   !> Return the total number of computational particles for an \c aero_mode.
00540   real(kind=dp) function aero_mode_number(aero_mode, aero_weight)
00541 
00542     !> Aero_mode to sample radius from.
00543     type(aero_mode_t), intent(in) :: aero_mode
00544     !> Aero weight.
00545     type(aero_weight_t), intent(in) :: aero_weight
00546 
00547     real(kind=dp) :: x_mean_prime
00548     real(kind=dp), allocatable :: weighted_num_conc(:)
00549 
00550     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00551        if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00552           aero_mode_number = aero_mode%num_conc * aero_weight%comp_vol
00553        elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00554             .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00555             x_mean_prime = log10(aero_mode%char_radius) &
00556             - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
00557             * log(10d0)
00558           aero_mode_number = aero_mode%num_conc * aero_weight%comp_vol &
00559                * aero_weight%ref_radius**aero_weight%exponent &
00560                * exp((x_mean_prime**2 - log10(aero_mode%char_radius)**2) &
00561                / (2d0 * aero_mode%log10_std_dev_radius**2))
00562        else
00563           call die_msg(466668240, "unknown aero_weight type: " &
00564                // trim(integer_to_string(aero_weight%type)))
00565        end if
00566     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00567        if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00568           aero_mode_number = aero_mode%num_conc * aero_weight%comp_vol
00569        else
00570           call die_msg(822252601, &
00571                "cannot use exponential modes with weighting")
00572        end if
00573     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00574        aero_mode_number = aero_mode%num_conc &
00575             / aero_weight_num_conc_at_radius(aero_weight, &
00576             aero_mode%char_radius)
00577     elseif (aero_mode%type == AERO_MODE_TYPE_SAMPLED) then
00578        allocate(weighted_num_conc(size(aero_mode%sample_num_conc)))
00579        call aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
00580             weighted_num_conc)
00581        aero_mode_number = sum(weighted_num_conc) * aero_weight%comp_vol
00582        deallocate(weighted_num_conc)
00583     else
00584        call die_msg(901140225, "unknown aero_mode type: " &
00585             // trim(integer_to_string(aero_mode%type)))
00586     end if
00587 
00588   end function aero_mode_number
00589 
00590 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00591 
00592   !> Return a radius randomly sampled from the mode distribution.
00593   subroutine aero_mode_sample_radius(aero_mode, aero_weight, radius)
00594 
00595     !> Aero_mode to sample radius from.
00596     type(aero_mode_t), intent(in) :: aero_mode
00597     !> Aero weight.
00598     type(aero_weight_t), intent(in) :: aero_weight
00599     !> Sampled radius (m).
00600     real(kind=dp), intent(out) :: radius
00601 
00602     real(kind=dp) :: x_mean_prime, x0, x1, x, r, inv_nc0, inv_nc1, inv_nc
00603     integer :: i_sample
00604     real(kind=dp), allocatable :: weighted_num_conc(:)
00605 
00606     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00607        if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00608           x_mean_prime = log10(aero_mode%char_radius)
00609        elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00610             .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00611           x_mean_prime = log10(aero_mode%char_radius) &
00612                - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
00613                * log(10d0)
00614        else
00615           call die_msg(517376844, "unknown aero_weight type: " &
00616                // trim(integer_to_string(aero_weight%type)))
00617        end if
00618        radius = 10d0**rand_normal(x_mean_prime, &
00619             aero_mode%log10_std_dev_radius)
00620     elseif (aero_mode%type == AERO_MODE_TYPE_SAMPLED) then
00621        allocate(weighted_num_conc(size(aero_mode%sample_num_conc)))
00622        call aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
00623             weighted_num_conc)
00624        i_sample = sample_cts_pdf(weighted_num_conc)
00625        deallocate(weighted_num_conc)
00626        if ((aero_weight%type == AERO_WEIGHT_TYPE_NONE) &
00627             .or. (((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00628             .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) &
00629             .and. (aero_weight%exponent == 0d0))) then
00630           x0 = log(aero_mode%sample_radius(i_sample))
00631           x1 = log(aero_mode%sample_radius(i_sample + 1))
00632           r = pmc_random()
00633           x = (1d0 - r) * x0 + r * x1
00634           radius = exp(x)
00635        elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00636             .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00637           inv_nc0 = 1d0 / aero_weight_num_conc_at_radius(aero_weight, &
00638                aero_mode%sample_radius(i_sample))
00639           inv_nc1 = 1d0 / aero_weight_num_conc_at_radius(aero_weight, &
00640                aero_mode%sample_radius(i_sample + 1))
00641           r = pmc_random()
00642           inv_nc = (1d0 - r) * inv_nc0 + r * inv_nc1
00643           radius = aero_weight_radius_at_num_conc(aero_weight, 1d0 / inv_nc)
00644        else
00645           call die_msg(769131141, "unknown aero_weight type: " &
00646                // trim(integer_to_string(aero_weight%type)))
00647        end if
00648     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00649        if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00650           radius = vol2rad(- rad2vol(aero_mode%char_radius) &
00651                * log(pmc_random()))
00652        elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00653             .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00654           call die_msg(678481276, &
00655                "cannot use exponential modes with weighting")
00656        else
00657           call die_msg(301787712, "unknown aero_weight type: " &
00658                // trim(integer_to_string(aero_weight%type)))
00659        end if
00660     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00661        radius = aero_mode%char_radius
00662     else
00663        call die_msg(749122931, "Unknown aero_mode type: " &
00664             // trim(integer_to_string(aero_mode%type)))
00665     end if
00666 
00667   end subroutine aero_mode_sample_radius
00668 
00669 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00670 
00671   !> Read volume fractions from a data file.
00672   subroutine spec_file_read_vol_frac(file, aero_data, vol_frac)
00673 
00674     !> Spec file to read mass fractions from.
00675     type(spec_file_t), intent(inout) :: file
00676     !> Aero_data data.
00677     type(aero_data_t), intent(in) :: aero_data
00678     !> Aerosol species volume fractions.
00679     real(kind=dp), intent(inout) :: vol_frac(:)
00680 
00681     integer :: n_species, species, i
00682     character(len=SPEC_LINE_MAX_VAR_LEN), pointer :: species_name(:)
00683     real(kind=dp), pointer :: species_data(:,:)
00684     real(kind=dp) :: tot_vol_frac
00685 
00686     !> \page input_format_mass_frac Input File Format: Aerosol Mass Fractions
00687     !!
00688     !! An aerosol mass fractions file must consist of one line per
00689     !! aerosol species, with each line having the species name
00690     !! followed by the species mass fraction in each aerosol
00691     !! particle. The valid species names are those specfied by the
00692     !! \ref input_format_aero_data file, but not all species have to
00693     !! be listed. Any missing species will have proportions of
00694     !! zero. If the proportions do not sum to 1 then they will be
00695     !! normalized before use. For example, a mass fractions file file
00696     !! could contain:
00697     !! <pre>
00698     !! # species   proportion
00699     !! OC          0.3
00700     !! BC          0.7
00701     !! </pre>
00702     !! indicating that the diesel particles are 30% organic carbon and
00703     !! 70% black carbon.
00704     !!
00705     !! See also:
00706     !!   - \ref spec_file_format --- the input file text format
00707     !!   - \ref input_format_aero_dist --- the format for a complete
00708     !!     aerosol distribution with several modes
00709     !!   - \ref input_format_aero_mode --- the format for each mode
00710     !!     of an aerosol distribution
00711 
00712     ! read the aerosol data from the specified file
00713     allocate(species_name(0))
00714     allocate(species_data(0,0))
00715     call spec_file_read_real_named_array(file, 0, species_name, &
00716          species_data)
00717 
00718     ! check the data size
00719     n_species = size(species_data, 1)
00720     if (n_species < 1) then
00721        call die_msg(628123166, 'file ' // trim(file%name) &
00722             // ' must contain at least one line of data')
00723     end if
00724     if (size(species_data, 2) /= 1) then
00725        call die_msg(427666881, 'each line in file ' &
00726             // trim(file%name) // ' must contain exactly one data value')
00727     end if
00728 
00729     ! copy over the data
00730     vol_frac = 0d0
00731     do i = 1,n_species
00732        species = aero_data_spec_by_name(aero_data, species_name(i))
00733        if (species == 0) then
00734           call die_msg(775942501, 'unknown species ' &
00735                // trim(species_name(i)) // ' in file ' &
00736                // trim(file%name))
00737        end if
00738        vol_frac(species) = species_data(i,1)
00739     end do
00740     deallocate(species_name)
00741     deallocate(species_data)
00742 
00743     ! convert mass fractions to volume fractions
00744     vol_frac = vol_frac / aero_data%density
00745     
00746     ! normalize
00747     tot_vol_frac = sum(vol_frac)
00748     if ((minval(vol_frac) < 0d0) .or. (tot_vol_frac <= 0d0)) then
00749        call die_msg(356648030, 'vol_frac in ' // trim(file%name) &
00750             // ' is not positive')
00751     end if
00752     vol_frac = vol_frac / tot_vol_frac
00753 
00754   end subroutine spec_file_read_vol_frac
00755 
00756 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00757 
00758   !> Read a size distribution from a data file.
00759   subroutine spec_file_read_size_dist(file, sample_radius, sample_num_conc)
00760 
00761     !> Spec file to read size distribution from.
00762     type(spec_file_t), intent(inout) :: file
00763     !> Sample radius values (m).
00764     real(kind=dp), pointer :: sample_radius(:)
00765     !> Sample number concentrations (m^{-3}).
00766     real(kind=dp), pointer :: sample_num_conc(:)
00767 
00768     character(len=SPEC_LINE_MAX_VAR_LEN), pointer :: names(:)
00769     real(kind=dp), pointer :: data(:,:)
00770     integer :: n_sample, i_sample
00771 
00772     !> \page input_format_size_dist Input File Format: Size Distribution
00773     !!
00774     !! A size distribution file must consist of two lines:
00775     !! - the first line must begin with \c diam and be followed by
00776     !!   \f$N + 1\f$ space-separated real scalars, giving the radii
00777     !!   \f$r_1,\ldots,r_{N+1}\f$ of bin edges (m) --- these must be
00778     !!   in increasing order, so \f$r_i < r_{i+1}\f$
00779     !! - the second line must begin with \c num_conc and be followed
00780     !!   by \f$N\f$ space-separated real scalars, giving the number
00781     !!   concenrations \f$C_1,\ldots,C_N\f$ in each bin (#/m^3) ---
00782     !!   \f$C_i\f$ is the total number concentrations of particles
00783     !!   with diameters in \f$[r_i, r_{i+1}]\f$
00784     !!
00785     !! The resulting size distribution is taken to be piecewise
00786     !! constant in log-diameter coordinates.
00787     !!
00788     !! Example: a size distribution could be:
00789     !! <pre>
00790     !! diam 1e-7 1e-6 1e-5  # bin edge diameters (m)
00791     !! num_conc 1e9 1e8     # bin number concentrations (m^{-3})
00792     !! </pre>
00793     !! This distribution has 1e9 particles per cubic meter with
00794     !! diameters between 0.1 micron and 1 micron, and 1e8 particles
00795     !! per cubic meter with diameters between 1 micron and 10 micron.
00796     !!
00797     !! See also:
00798     !!   - \ref spec_file_format --- the input file text format
00799     !!   - \ref input_format_aero_dist --- the format for a complete
00800     !!     aerosol distribution with several modes
00801     !!   - \ref input_format_aero_mode --- the format for each mode
00802     !!     of an aerosol distribution
00803 
00804     ! read the data from the file
00805     allocate(names(0))
00806     allocate(data(0,0))
00807     call spec_file_read_real_named_array(file, 1, names, data)
00808     call spec_file_check_name(file, 'diam', names(1))
00809     n_sample = size(data,2) - 1
00810     call spec_file_assert_msg(669011124, file, n_sample >= 1, &
00811          'must have at least two diam values')
00812 
00813     deallocate(sample_radius)
00814     allocate(sample_radius(n_sample + 1))
00815     sample_radius = diam2rad(data(1,:))
00816     do i_sample = 1,n_sample
00817        call spec_file_assert_msg(528089871, file, &
00818             sample_radius(i_sample) < sample_radius(i_sample + 1), &
00819             'diam values must be strictly increasing')
00820     end do
00821     
00822     call spec_file_read_real_named_array(file, 1, names, data)
00823     call spec_file_check_name(file, 'num_conc', names(1))
00824 
00825     call spec_file_assert_msg(721029144, file, size(data, 2) == n_sample, &
00826          'must have one fewer num_conc than diam values')
00827 
00828     deallocate(sample_num_conc)
00829     allocate(sample_num_conc(n_sample))
00830     sample_num_conc = data(1,:)
00831     do i_sample = 1,n_sample
00832        call spec_file_assert_msg(356490397, file, &
00833             sample_num_conc(i_sample) >= 0d0, &
00834             'num_conc values must be non-negative')
00835     end do
00836 
00837   end subroutine spec_file_read_size_dist
00838 
00839 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00840 
00841   !> Read one mode of an aerosol distribution (number concentration,
00842   !> volume fractions, and mode shape).
00843   subroutine spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
00844 
00845     !> Spec file.
00846     type(spec_file_t), intent(inout) :: file
00847     !> Aero_data data.
00848     type(aero_data_t), intent(inout) :: aero_data
00849     !> Aerosol mode.
00850     type(aero_mode_t), intent(inout) :: aero_mode
00851     !> If eof instead of reading data.
00852     logical :: eof
00853 
00854     character(len=SPEC_LINE_MAX_VAR_LEN) :: tmp_str, mode_type
00855     character(len=SPEC_LINE_MAX_VAR_LEN) :: mass_frac_filename
00856     character(len=SPEC_LINE_MAX_VAR_LEN) :: size_dist_filename
00857     type(spec_line_t) :: line
00858     type(spec_file_t) :: mass_frac_file, size_dist_file
00859     real(kind=dp) :: diam
00860 
00861     ! note that doxygen's automatic list creation breaks on the list
00862     ! below for some reason, so we just use html format
00863 
00864     !> \page input_format_aero_mode Input File Format: Aerosol Distribution Mode
00865     !!
00866     !! An aerosol distribution mode has the parameters:
00867     !! <ul>
00868     !! <li> \b mode_name (string): the name of the mode (for
00869     !!      informational purposes only)
00870     !! <li> \b mass_frac (string): name of file from which to read the
00871     !!      species mass fractions --- the file format should
00872     !!      be \subpage input_format_mass_frac
00873     !! <li> \b mode_type (string): the functional form of the mode ---
00874     !!      must be one of: \c log_normal for a log-normal
00875     !!      distribution; \c exp for an exponential distribution; \c
00876     !!      mono for a mono-disperse distribution; or \c sampled for a
00877     !!      sampled distribution
00878     !! <li> if \c mode_type is \c log_normal then the mode distribution
00879     !!      shape is
00880     !!      \f[ n(\log D) {\rm d}\log D
00881     !!      = \frac{N_{\rm total}}{\sqrt{2\pi} \log \sigma_{\rm g}}
00882     !!      \exp\left(\frac{(\log D - \log D_{\rm gn})^2}{2 \log ^2
00883     !!      \sigma_{\rm g}}\right)
00884     !!      {\rm d}\log D \f]
00885     !!      and the following parameters are:
00886     !!      <ul>
00887     !!      <li> \b num_conc (real, unit 1/m^3): the total number
00888     !!           concentration \f$N_{\rm total}\f$ of the mode
00889     !!      <li> \b geom_mean_diam (real, unit m): the geometric mean
00890     !!           diameter \f$D_{\rm gn}\f$
00891     !!      <li> \b log10_geom_std_dev (real, dimensionless):
00892     !!           \f$\log_{10}\f$ of the geometric standard deviation
00893     !!           \f$\sigma_{\rm g}\f$ of the diameter
00894     !!      </ul>
00895     !! <li> if \c mode_type is \c exp then the mode distribution shape is
00896     !!      \f[ n(v) {\rm d}v = \frac{N_{\rm total}}{v_{\rm \mu}}
00897     !!      \exp\left(- \frac{v}{v_{\rm \mu}}\right)
00898     !!      {\rm d}v \f]
00899     !!      and the following parameters are:
00900     !!      <ul>
00901     !!      <li> \b num_conc (real, unit 1/m^3): the total number
00902     !!           concentration \f$N_{\rm total}\f$ of the mode
00903     !!      <li> \b diam_at_mean_vol (real, unit m): the diameter
00904     !!           \f$D_{\rm \mu}\f$ such that \f$v_{\rm \mu}
00905     !!           = \frac{\pi}{6} D^3_{\rm \mu}\f$
00906     !!      </ul>
00907     !! <li> if \c mode_type is \c mono then the mode distribution shape
00908     !!      is a delta distribution at diameter \f$D_0\f$ and the
00909     !!      following parameters are:
00910     !!      <ul>
00911     !!      <li> \b num_conc (real, unit 1/m^3): the total number
00912     !!           concentration \f$N_{\rm total}\f$ of the mode
00913     !!      <li> \b radius (real, unit m): the radius \f$R_0\f$ of the
00914     !!           particles, so that \f$D_0 = 2 R_0\f$
00915     !!      </ul>
00916     !! <li> if \c mode_type is \c sampled then the mode distribution
00917     !!      shape is piecewise constant (in log-diameter coordinates)
00918     !!      and the following parameters are:
00919     !!      <ul>
00920     !!      <li> \b size_dist (string): name of file from which to
00921     !!           read the size distribution --- the file format should
00922     !!           be \subpage input_format_size_dist
00923     !!      </ul>
00924     !! </ul>
00925     !!
00926     !! Example:
00927     !! <pre>
00928     !! mode_name diesel          # mode name (descriptive only)
00929     !! mass_frac comp_diesel.dat # mass fractions in each aerosol particle
00930     !! mode_type log_normal      # type of distribution
00931     !! num_conc 1.6e8            # particle number density (#/m^3)
00932     !! geom_mean_diam 2.5e-8     # geometric mean diameter (m)
00933     !! log10_geom_std_dev 0.24   # log_10 of geometric standard deviation
00934     !! </pre>
00935     !!
00936     !! See also:
00937     !!   - \ref spec_file_format --- the input file text format
00938     !!   - \ref input_format_aero_dist --- the format for a complete
00939     !!     aerosol distribution with several modes
00940     !!   - \ref input_format_mass_frac --- the format for the mass
00941     !!     fractions file
00942 
00943     call spec_line_allocate(line)
00944     call spec_file_read_line(file, line, eof)
00945     if (.not. eof) then
00946        call spec_file_check_line_name(file, line, "mode_name")
00947        call spec_file_check_line_length(file, line, 1)
00948        call aero_mode_deallocate(aero_mode)
00949        call aero_mode_allocate_size(aero_mode, aero_data%n_spec)
00950        tmp_str = line%data(1) ! hack to avoid gfortran warning
00951        aero_mode%name = tmp_str(1:AERO_MODE_NAME_LEN)
00952        aero_mode%source = aero_data_source_by_name(aero_data, aero_mode%name)
00953        call spec_file_read_string(file, 'mass_frac', mass_frac_filename)
00954        call spec_file_open(mass_frac_filename, mass_frac_file)
00955        call spec_file_read_vol_frac(mass_frac_file, aero_data, &
00956             aero_mode%vol_frac)
00957        call spec_file_close(mass_frac_file)
00958        call spec_file_read_string(file, 'mode_type', mode_type)
00959        if (trim(mode_type) == 'log_normal') then
00960           aero_mode%type = AERO_MODE_TYPE_LOG_NORMAL
00961           call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
00962           call spec_file_read_real(file, 'geom_mean_diam', diam)
00963           aero_mode%char_radius = diam2rad(diam)
00964           call spec_file_read_real(file, 'log10_geom_std_dev', &
00965                aero_mode%log10_std_dev_radius)
00966        elseif (trim(mode_type) == 'exp') then
00967           aero_mode%type = AERO_MODE_TYPE_EXP
00968           call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
00969           call spec_file_read_real(file, 'diam_at_mean_vol', diam)
00970           aero_mode%char_radius = diam2rad(diam)
00971        elseif (trim(mode_type) == 'mono') then
00972           aero_mode%type = AERO_MODE_TYPE_MONO
00973           call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
00974           call spec_file_read_real(file, 'diam', diam)
00975           aero_mode%char_radius = diam2rad(diam)
00976        elseif (trim(mode_type) == 'sampled') then
00977           aero_mode%type = AERO_MODE_TYPE_SAMPLED
00978           call spec_file_read_string(file, 'size_dist', size_dist_filename)
00979           call spec_file_open(size_dist_filename, size_dist_file)
00980           call spec_file_read_size_dist(size_dist_file, &
00981                aero_mode%sample_radius, aero_mode%sample_num_conc)
00982           call spec_file_close(size_dist_file)
00983        else
00984           call spec_file_die_msg(729472928, file, &
00985                "Unknown distribution mode type: " // trim(mode_type))
00986        end if
00987     end if
00988     call spec_line_deallocate(line)
00989 
00990   end subroutine spec_file_read_aero_mode
00991 
00992 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00993 
00994   !> Determines the number of bytes required to pack the given value.
00995   integer function pmc_mpi_pack_size_aero_mode(val)
00996 
00997     !> Value to pack.
00998     type(aero_mode_t), intent(in) :: val
00999 
01000     pmc_mpi_pack_size_aero_mode = &
01001          pmc_mpi_pack_size_string(val%name) &
01002          + pmc_mpi_pack_size_integer(val%type) &
01003          + pmc_mpi_pack_size_real(val%char_radius) &
01004          + pmc_mpi_pack_size_real(val%log10_std_dev_radius) &
01005          + pmc_mpi_pack_size_real_array(val%sample_radius) &
01006          + pmc_mpi_pack_size_real_array(val%sample_num_conc) &
01007          + pmc_mpi_pack_size_real(val%num_conc) &
01008          + pmc_mpi_pack_size_real_array(val%vol_frac) &
01009          + pmc_mpi_pack_size_integer(val%source)
01010 
01011   end function pmc_mpi_pack_size_aero_mode
01012 
01013 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01014 
01015   !> Packs the given value into the buffer, advancing position.
01016   subroutine pmc_mpi_pack_aero_mode(buffer, position, val)
01017 
01018     !> Memory buffer.
01019     character, intent(inout) :: buffer(:)
01020     !> Current buffer position.
01021     integer, intent(inout) :: position
01022     !> Value to pack.
01023     type(aero_mode_t), intent(in) :: val
01024 
01025 #ifdef PMC_USE_MPI
01026     integer :: prev_position
01027 
01028     prev_position = position
01029     call pmc_mpi_pack_string(buffer, position, val%name)
01030     call pmc_mpi_pack_integer(buffer, position, val%type)
01031     call pmc_mpi_pack_real(buffer, position, val%char_radius)
01032     call pmc_mpi_pack_real(buffer, position, val%log10_std_dev_radius)
01033     call pmc_mpi_pack_real_array(buffer, position, val%sample_radius)
01034     call pmc_mpi_pack_real_array(buffer, position, val%sample_num_conc)
01035     call pmc_mpi_pack_real(buffer, position, val%num_conc)
01036     call pmc_mpi_pack_real_array(buffer, position, val%vol_frac)
01037     call pmc_mpi_pack_integer(buffer, position, val%source)
01038     call assert(497092471, &
01039          position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
01040 #endif
01041 
01042   end subroutine pmc_mpi_pack_aero_mode
01043 
01044 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01045 
01046   !> Unpacks the given value from the buffer, advancing position.
01047   subroutine pmc_mpi_unpack_aero_mode(buffer, position, val)
01048 
01049     !> Memory buffer.
01050     character, intent(inout) :: buffer(:)
01051     !> Current buffer position.
01052     integer, intent(inout) :: position
01053     !> Value to pack.
01054     type(aero_mode_t), intent(inout) :: val
01055 
01056 #ifdef PMC_USE_MPI
01057     integer :: prev_position
01058 
01059     prev_position = position
01060     call pmc_mpi_unpack_string(buffer, position, val%name)
01061     call pmc_mpi_unpack_integer(buffer, position, val%type)
01062     call pmc_mpi_unpack_real(buffer, position, val%char_radius)
01063     call pmc_mpi_unpack_real(buffer, position, val%log10_std_dev_radius)
01064     call pmc_mpi_unpack_real_array(buffer, position, val%sample_radius)
01065     call pmc_mpi_unpack_real_array(buffer, position, val%sample_num_conc)
01066     call pmc_mpi_unpack_real(buffer, position, val%num_conc)
01067     call pmc_mpi_unpack_real_array(buffer, position, val%vol_frac)
01068     call pmc_mpi_unpack_integer(buffer, position, val%source)
01069     call assert(874467577, &
01070          position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
01071 #endif
01072 
01073   end subroutine pmc_mpi_unpack_aero_mode
01074 
01075 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01076 
01077 end module pmc_aero_mode