PartMC 2.1.3
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 
00037   !> An aerosol size distribution mode.
00038   !!
00039   !! Each mode is assumed to be fully internally mixed so that every
00040   !! particle has the same composition. The \c num_conc array then
00041   !! stores the number concentration distribution.
00042   type aero_mode_t
00043      !> Mode name, used to track particle sources.
00044      character(len=AERO_MODE_NAME_LEN) :: name
00045      !> Mode type (given by module constants).
00046      integer :: type
00047      !> Characteristic radius, with meaning dependent on mode type (m).
00048      real(kind=dp) :: char_radius
00049      !> Log base 10 of geometric standard deviation of radius, if
00050      !> necessary (m).
00051      real(kind=dp) :: log10_std_dev_radius
00052      !> Total number concentration of mode (#/m^3).
00053      real(kind=dp) :: num_conc
00054      !> Species fractions by volume [length \c aero_data%%n_spec] (1).
00055      real(kind=dp), pointer :: vol_frac(:)
00056      !> Source number.
00057      integer :: source
00058   end type aero_mode_t
00059 
00060 contains
00061 
00062 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00063 
00064   !> Return a string representation of a kernel type.
00065   character(len=AERO_MODE_TYPE_LEN) function aero_mode_type_to_string(type)
00066 
00067     !> Aero mode type.
00068     integer, intent(in) :: type
00069    
00070     if (type == AERO_MODE_TYPE_INVALID) then
00071        aero_mode_type_to_string = "invalid"
00072     elseif (type == AERO_MODE_TYPE_LOG_NORMAL) then
00073        aero_mode_type_to_string = "log_normal"
00074     elseif (type == AERO_MODE_TYPE_EXP) then
00075        aero_mode_type_to_string = "exp"
00076     elseif (type == AERO_MODE_TYPE_MONO) then
00077        aero_mode_type_to_string = "mono"
00078     else
00079        aero_mode_type_to_string = "unknown"
00080     end if
00081 
00082   end function aero_mode_type_to_string
00083 
00084 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00085 
00086   !> Allocates an aero_mode.
00087   subroutine aero_mode_allocate(aero_mode)
00088 
00089     !> Aerosol mode.
00090     type(aero_mode_t), intent(out) :: aero_mode
00091 
00092     allocate(aero_mode%vol_frac(0))
00093 
00094   end subroutine aero_mode_allocate
00095 
00096 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00097 
00098   !> Allocates an aero_mode of the given size.
00099   subroutine aero_mode_allocate_size(aero_mode, n_spec)
00100 
00101     !> Aerosol mode.
00102     type(aero_mode_t), intent(out) :: aero_mode
00103     !> Number of species.
00104     integer, intent(in) :: n_spec
00105 
00106     allocate(aero_mode%vol_frac(n_spec))
00107 
00108   end subroutine aero_mode_allocate_size
00109 
00110 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00111 
00112   !> Free all storage.
00113   subroutine aero_mode_deallocate(aero_mode)
00114 
00115     !> Aerosol mode.
00116     type(aero_mode_t), intent(inout) :: aero_mode
00117 
00118     deallocate(aero_mode%vol_frac)
00119 
00120   end subroutine aero_mode_deallocate
00121 
00122 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00123 
00124   !> Copy an aero_mode.
00125   subroutine aero_mode_copy(aero_mode_from, aero_mode_to)
00126 
00127     !> Aerosol mode original.
00128     type(aero_mode_t), intent(in) :: aero_mode_from
00129     !> Aerosol mode copy.
00130     type(aero_mode_t), intent(inout) :: aero_mode_to
00131 
00132     call aero_mode_deallocate(aero_mode_to)
00133     call aero_mode_allocate_size(aero_mode_to, size(aero_mode_from%vol_frac))
00134     aero_mode_to%name = aero_mode_from%name
00135     aero_mode_to%type = aero_mode_from%type
00136     aero_mode_to%char_radius = aero_mode_from%char_radius
00137     aero_mode_to%log10_std_dev_radius = aero_mode_from%log10_std_dev_radius
00138     aero_mode_to%num_conc = aero_mode_from%num_conc
00139     aero_mode_to%vol_frac = aero_mode_from%vol_frac
00140     aero_mode_to%source = aero_mode_from%source
00141 
00142   end subroutine aero_mode_copy
00143 
00144 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00145 
00146   !> Compute a log-normal distribution, normalized so that
00147   !> sum(num_conc(k) * log_width) = 1
00148   subroutine num_conc_log_normal(geom_mean_radius, log10_sigma_g, &
00149        bin_grid, num_conc)
00150     
00151     !> Geometric mean radius (m).
00152     real(kind=dp), intent(in) :: geom_mean_radius
00153     !> log_10(geom. std. dev.) (1).
00154     real(kind=dp), intent(in) :: log10_sigma_g
00155     !> Bin grid.
00156     type(bin_grid_t), intent(in) :: bin_grid
00157     !> Number concentration (#(ln(r))d(ln(r))).
00158     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00159     
00160     integer :: k
00161     
00162     do k = 1,bin_grid%n_bin
00163        num_conc(k) = 1d0 / (sqrt(2d0 * const%pi) * log10_sigma_g) * &
00164             dexp(-(dlog10(bin_grid%center_radius(k)) &
00165             - dlog10(geom_mean_radius))**2d0 &
00166             / (2d0 * log10_sigma_g**2d0)) / dlog(10d0)
00167     end do
00168     
00169     ! The formula above was originally for a distribution in
00170     ! log_10(r), while we are using log_e(r) for our bin grid. The
00171     ! division by dlog(10) at the end corrects for this.
00172 
00173     ! Remember that log_e(r) = log_10(r) * log_e(10).
00174     
00175   end subroutine num_conc_log_normal
00176   
00177 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00178 
00179   !> Compute a log-normal distribution in volume.
00180   subroutine vol_conc_log_normal(geom_mean_radius, log10_sigma_g, &
00181        bin_grid, vol_conc)
00182     
00183     !> Geometric mean radius (m).
00184     real(kind=dp), intent(in) :: geom_mean_radius
00185     !> log_10(geom. std. dev.) (1).
00186     real(kind=dp), intent(in) :: log10_sigma_g
00187     !> Bin grid.
00188     type(bin_grid_t), intent(in) :: bin_grid
00189     !> Volume concentration (V(ln(r))d(ln(r))).
00190     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00191     
00192     real(kind=dp) :: num_conc(bin_grid%n_bin)
00193 
00194     call num_conc_log_normal(geom_mean_radius, log10_sigma_g, &
00195          bin_grid, num_conc)
00196     vol_conc = num_conc * rad2vol(bin_grid%center_radius)
00197     
00198   end subroutine vol_conc_log_normal
00199   
00200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00201 
00202   !> Exponential distribution in volume
00203   !> \f[ n(v) = \frac{1}{\rm mean-vol} \exp(- v / {\rm mean-vol}) \f]
00204   !> Normalized so that sum(num_conc(k) * log_width) = 1
00205   subroutine num_conc_exp(radius_at_mean_vol, bin_grid, num_conc)
00206     
00207     !> Radius at mean volume (m).
00208     real(kind=dp), intent(in) :: radius_at_mean_vol
00209     !> Bin grid.
00210     type(bin_grid_t), intent(in) :: bin_grid
00211     !> Number concentration (#(ln(r))d(ln(r))).
00212     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00213     
00214     integer :: k
00215     real(kind=dp) :: mean_vol, num_conc_vol
00216     
00217     mean_vol = rad2vol(radius_at_mean_vol)
00218     do k = 1,bin_grid%n_bin
00219        num_conc_vol = 1d0 / mean_vol &
00220             * exp(-(rad2vol(bin_grid%center_radius(k)) / mean_vol))
00221        call vol_to_lnr(bin_grid%center_radius(k), num_conc_vol, num_conc(k))
00222     end do
00223     
00224   end subroutine num_conc_exp
00225   
00226 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00227 
00228   !> Exponential distribution in volume.
00229   subroutine vol_conc_exp(radius_at_mean_vol, bin_grid, vol_conc)
00230     
00231     !> Radius at mean volume (m).
00232     real(kind=dp), intent(in) :: radius_at_mean_vol
00233     !> Bin grid.
00234     type(bin_grid_t), intent(in) :: bin_grid
00235     !> Volume concentration (V(ln(r))d(ln(r))).
00236     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00237     
00238     real(kind=dp) :: num_conc(bin_grid%n_bin)
00239 
00240     call num_conc_exp(radius_at_mean_vol, bin_grid, num_conc)
00241     vol_conc = num_conc * rad2vol(bin_grid%center_radius)
00242     
00243   end subroutine vol_conc_exp
00244   
00245 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00246 
00247   !> Mono-disperse distribution.
00248   !> Normalized so that sum(num_conc(k) * log_width) = 1
00249   subroutine num_conc_mono(radius, bin_grid, num_conc)
00250     
00251     !> Radius of each particle (m^3).
00252     real(kind=dp), intent(in) :: radius
00253     !> Bin grid.
00254     type(bin_grid_t), intent(in) :: bin_grid
00255     !> Number concentration (#(ln(r))d(ln(r))).
00256     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00257     
00258     integer :: k
00259 
00260     num_conc = 0d0
00261     k = bin_grid_particle_in_bin(bin_grid, radius)
00262     num_conc(k) = 1d0 / bin_grid%log_width
00263     
00264   end subroutine num_conc_mono
00265   
00266 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00267 
00268   !> Mono-disperse distribution in volume.
00269   subroutine vol_conc_mono(radius, bin_grid, vol_conc)
00270     
00271     !> Radius of each particle (m^3).
00272     real(kind=dp), intent(in) :: radius
00273     !> Bin grid.
00274     type(bin_grid_t), intent(in) :: bin_grid
00275     !> Volume concentration (V(ln(r))d(ln(r))).
00276     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin)
00277     
00278     integer :: k
00279 
00280     vol_conc = 0d0
00281     k = bin_grid_particle_in_bin(bin_grid, radius)
00282     vol_conc(k) = 1d0 / bin_grid%log_width * rad2vol(radius)
00283     
00284   end subroutine vol_conc_mono
00285   
00286 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00287 
00288   !> Return the binned number concentration for an aero_mode.
00289   subroutine aero_mode_num_conc(aero_mode, bin_grid, aero_data, &
00290        num_conc)
00291 
00292     !> Aero mode for which to compute number concentration.
00293     type(aero_mode_t), intent(in) :: aero_mode
00294     !> Bin grid.
00295     type(bin_grid_t), intent(in) :: bin_grid
00296     !> Aerosol data.
00297     type(aero_data_t), intent(in) :: aero_data
00298     !> Number concentration (#(ln(r))d(ln(r))).
00299     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00300 
00301     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00302        call num_conc_log_normal(aero_mode%char_radius, &
00303             aero_mode%log10_std_dev_radius, bin_grid, num_conc)
00304     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00305        call num_conc_exp(aero_mode%char_radius, bin_grid, num_conc)
00306     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00307        call num_conc_mono(aero_mode%char_radius, bin_grid, num_conc)
00308     else
00309        call die_msg(719625922, "unknown aero_mode type: " &
00310             // trim(integer_to_string(aero_mode%type)))
00311     end if
00312     num_conc = num_conc * aero_mode%num_conc
00313 
00314   end subroutine aero_mode_num_conc
00315 
00316 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00317 
00318   !> Return the binned per-species volume concentration for an
00319   !> aero_mode.
00320   subroutine aero_mode_vol_conc(aero_mode, bin_grid, aero_data, &
00321        vol_conc)
00322 
00323     !> Aero mode for which to compute volume concentration.
00324     type(aero_mode_t), intent(in) :: aero_mode
00325     !> Bin grid.
00326     type(bin_grid_t), intent(in) :: bin_grid
00327     !> Aerosol data.
00328     type(aero_data_t), intent(in) :: aero_data
00329     !> Volume concentration (V(ln(r))d(ln(r))).
00330     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin, aero_data%n_spec)
00331 
00332     integer :: i_spec
00333     real(kind=dp) :: vol_conc_total(bin_grid%n_bin)
00334 
00335     if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00336        call vol_conc_log_normal(aero_mode%char_radius, &
00337             aero_mode%log10_std_dev_radius, bin_grid, vol_conc_total)
00338     elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00339        call vol_conc_exp(aero_mode%char_radius, bin_grid, &
00340             vol_conc_total)
00341     elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00342        call vol_conc_mono(aero_mode%char_radius, bin_grid, &
00343             vol_conc_total)
00344     else
00345        call die_msg(314169653, "Unknown aero_mode type: " &
00346             // trim(integer_to_string(aero_mode%type)))
00347     end if
00348     vol_conc_total = vol_conc_total * aero_mode%num_conc
00349     do i_spec = 1,aero_data%n_spec
00350        vol_conc(:,i_spec) = vol_conc_total &
00351             * aero_mode%vol_frac(i_spec)
00352     end do
00353 
00354   end subroutine aero_mode_vol_conc
00355 
00356 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00357 
00358   !> Return the weighted number concentration for an \c aero_mode.
00359   real(kind=dp) function aero_mode_weighted_num_conc(aero_mode, aero_weight)
00360 
00361     !> Aero_mode to sample radius from.
00362     type(aero_mode_t), intent(in) :: aero_mode
00363     !> Aero weight.
00364     type(aero_weight_t), intent(in) :: aero_weight
00365 
00366     real(kind=dp) :: x_mean_prime
00367 
00368     if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00369        aero_mode_weighted_num_conc = aero_mode%num_conc
00370     elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00371          .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00372        if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00373           x_mean_prime = log10(aero_mode%char_radius) &
00374                - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
00375                * log(10d0)
00376           aero_mode_weighted_num_conc = aero_mode%num_conc &
00377                * aero_weight%ref_radius**aero_weight%exponent &
00378                * exp((x_mean_prime**2 - log10(aero_mode%char_radius)**2) &
00379                / (2d0 * aero_mode%log10_std_dev_radius**2))
00380        elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00381           call die_msg(822252601, &
00382                "cannot use exponential modes with weighting")
00383        elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00384           aero_mode_weighted_num_conc = aero_mode%num_conc &
00385                / aero_weight_value(aero_weight, aero_mode%char_radius)
00386        else
00387           call die_msg(901140225, "unknown aero_mode type: " &
00388                // trim(integer_to_string(aero_mode%type)))
00389        end if
00390     else
00391        call die_msg(742383510, "unknown aero_weight type: " &
00392             // trim(integer_to_string(aero_weight%type)))
00393     end if
00394 
00395   end function aero_mode_weighted_num_conc
00396 
00397 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00398 
00399   !> Return a radius randomly sampled from the mode distribution.
00400   subroutine aero_mode_sample_radius(aero_mode, aero_weight, radius)
00401 
00402     !> Aero_mode to sample radius from.
00403     type(aero_mode_t), intent(in) :: aero_mode
00404     !> Aero weight.
00405     type(aero_weight_t), intent(in) :: aero_weight
00406     !> Sampled radius (m).
00407     real(kind=dp), intent(out) :: radius
00408 
00409     real(kind=dp) :: x_mean_prime
00410 
00411     if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00412        if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00413           radius = 10d0**rand_normal(log10(aero_mode%char_radius), &
00414                aero_mode%log10_std_dev_radius)
00415        elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00416           radius = vol2rad(- rad2vol(aero_mode%char_radius) &
00417                * log(pmc_random()))
00418        elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00419           radius = aero_mode%char_radius
00420        else
00421           call die_msg(749122931, "Unknown aero_mode type: " &
00422                // trim(integer_to_string(aero_mode%type)))
00423        end if
00424     elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00425          .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00426        if (aero_mode%type == AERO_MODE_TYPE_LOG_NORMAL) then
00427           x_mean_prime = log10(aero_mode%char_radius) &
00428                - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
00429                * log(10d0)
00430           radius = 10d0**rand_normal(x_mean_prime, &
00431                aero_mode%log10_std_dev_radius)
00432        elseif (aero_mode%type == AERO_MODE_TYPE_EXP) then
00433           call die_msg(111024862, &
00434                "cannot use exponential modes with weighting")
00435        elseif (aero_mode%type == AERO_MODE_TYPE_MONO) then
00436           radius = aero_mode%char_radius
00437        else
00438           call die_msg(886417976, "unknown aero_mode type: " &
00439                // trim(integer_to_string(aero_mode%type)))
00440        end if
00441     else
00442        call die_msg(863127819, "unknown aero_weight type: " &
00443             // trim(integer_to_string(aero_weight%type)))
00444     end if
00445 
00446   end subroutine aero_mode_sample_radius
00447 
00448 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00449 
00450   !> Read volume fractions from a data file.
00451   subroutine spec_file_read_vol_frac(file, aero_data, vol_frac)
00452 
00453     !> Spec file to read mass fractions from.
00454     type(spec_file_t), intent(inout) :: file
00455     !> Aero_data data.
00456     type(aero_data_t), intent(in) :: aero_data
00457     !> Aerosol species volume fractions.
00458     real(kind=dp), intent(inout) :: vol_frac(:)
00459 
00460     integer :: n_species, species, i
00461     character(len=SPEC_LINE_MAX_VAR_LEN), pointer :: species_name(:)
00462     real(kind=dp), pointer :: species_data(:,:)
00463     real(kind=dp) :: tot_vol_frac
00464 
00465     !> \page input_format_mass_frac Input File Format: Aerosol Mass Fractions
00466     !!
00467     !! An aerosol mass fractions file must consist of one line per
00468     !! aerosol species, with each line having the species name
00469     !! followed by the species mass fraction in each aerosol
00470     !! particle. The valid species names are those specfied by the
00471     !! \ref input_format_aero_data file, but not all species have to
00472     !! be listed. Any missing species will have proportions of
00473     !! zero. If the proportions do not sum to 1 then they will be
00474     !! normalized before use. For example, a mass fractions file file
00475     !! could contain:
00476     !! <pre>
00477     !! # species   proportion
00478     !! OC          0.3
00479     !! BC          0.7
00480     !! </pre>
00481     !! indicating that the diesel particles are 30% organic carbon and
00482     !! 70% black carbon.
00483     !!
00484     !! See also:
00485     !!   - \ref spec_file_format --- the input file text format
00486     !!   - \ref input_format_aero_dist --- the format for a complete
00487     !!     aerosol distribution with several modes
00488     !!   - \ref input_format_aero_mode --- the format for each mode
00489     !!     of an aerosol distribution
00490 
00491     ! read the aerosol data from the specified file
00492     allocate(species_name(0))
00493     allocate(species_data(0,0))
00494     call spec_file_read_real_named_array(file, 0, species_name, &
00495          species_data)
00496 
00497     ! check the data size
00498     n_species = size(species_data, 1)
00499     if (n_species < 1) then
00500        call die_msg(427666881, 'file ' // trim(file%name) &
00501             // ' must contain at least one line of data')
00502     end if
00503     if (size(species_data, 2) /= 1) then
00504        call die_msg(427666881, 'each line in file ' &
00505             // trim(file%name) // ' must contain exactly one data value')
00506     end if
00507 
00508     ! copy over the data
00509     vol_frac = 0d0
00510     do i = 1,n_species
00511        species = aero_data_spec_by_name(aero_data, species_name(i))
00512        if (species == 0) then
00513           call die_msg(775942501, 'unknown species ' &
00514                // trim(species_name(i)) // ' in file ' &
00515                // trim(file%name))
00516        end if
00517        vol_frac(species) = species_data(i,1)
00518     end do
00519     deallocate(species_name)
00520     deallocate(species_data)
00521 
00522     ! convert mass fractions to volume fractions
00523     vol_frac = vol_frac / aero_data%density
00524     
00525     ! normalize
00526     tot_vol_frac = sum(vol_frac)
00527     if ((minval(vol_frac) < 0d0) .or. (tot_vol_frac <= 0d0)) then
00528        call die_msg(356648030, 'vol_frac in ' // trim(file%name) &
00529             // ' is not positive')
00530     end if
00531     vol_frac = vol_frac / tot_vol_frac
00532 
00533   end subroutine spec_file_read_vol_frac
00534 
00535 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00536 
00537   !> Read one mode of an aerosol distribution (number concentration,
00538   !> volume fractions, and mode shape).
00539   subroutine spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
00540 
00541     !> Spec file.
00542     type(spec_file_t), intent(inout) :: file
00543     !> Aero_data data.
00544     type(aero_data_t), intent(inout) :: aero_data
00545     !> Aerosol mode.
00546     type(aero_mode_t), intent(inout) :: aero_mode
00547     !> If eof instead of reading data.
00548     logical :: eof
00549 
00550     character(len=SPEC_LINE_MAX_VAR_LEN) :: tmp_str, mode_type
00551     character(len=SPEC_LINE_MAX_VAR_LEN) :: mass_frac_filename
00552     type(spec_line_t) :: line
00553     type(spec_file_t) :: mass_frac_file
00554     real(kind=dp) :: diam
00555 
00556     ! note that doxygen's automatic list creation breaks on the list
00557     ! below for some reason, so we just use html format
00558 
00559     !> \page input_format_aero_mode Input File Format: Aerosol Distribution Mode
00560     !!
00561     !! An aerosol distribution mode has the parameters:
00562     !! <ul>
00563     !! <li> \b mode_name (string): the name of the mode (for
00564     !!      informational purposes only)
00565     !! <li> \b mass_frac (string): name of file from which to read the
00566     !!      species mass fractions --- the file format should
00567     !!      be \subpage input_format_mass_frac
00568     !! <li> \b num_conc (real, unit 1/m^3): the total number
00569     !!      concentration \f$N_{\rm total}\f$ of the mode
00570     !! <li> \b mode_type (string): the functional form of the mode ---
00571     !!      must be one of: \c log_normal for a log-normal distribution;
00572     !!      \c exp for an exponential distribution; or \c mono for a
00573     !!      mono-disperse distribution
00574     !! <li> if \c mode_type is \c log_normal then the mode distribution
00575     !!      shape is
00576     !!      \f[ n(\log D) {\rm d}\log D
00577     !!      = \frac{N_{\rm total}}{\sqrt{2\pi} \log \sigma_{\rm g}}
00578     !!      \exp\left(\frac{(\log D - \log D_{\rm gn})^2}{2 \log ^2
00579     !!      \sigma_{\rm g}}\right)
00580     !!      {\rm d}\log D \f]
00581     !!      and the following parameters are:
00582     !!      <ul>
00583     !!      <li> \b geom_mean_diam (real, unit m): the geometric mean
00584     !!           diameter \f$D_{\rm gn}\f$
00585     !!      <li> \b log10_geom_std_dev (real, dimensionless):
00586     !!           \f$\log_{10}\f$ of the geometric standard deviation
00587     !!           \f$\sigma_{\rm g}\f$ of the diameter
00588     !!      </ul>
00589     !! <li> if \c mode_type is \c exp then the mode distribution shape is
00590     !!      \f[ n(v) {\rm d}v = \frac{N_{\rm total}}{v_{\rm \mu}}
00591     !!      \exp\left(- \frac{v}{v_{\rm \mu}}\right)
00592     !!      {\rm d}v \f]
00593     !!      and the following parameters are:
00594     !!      <ul>
00595     !!      <li> \b diam_at_mean_vol (real, unit m): the diameter
00596     !!           \f$D_{\rm \mu}\f$ such that \f$v_{\rm \mu}
00597     !!           = \frac{\pi}{6} D^3_{\rm \mu}\f$
00598     !!      </ul>
00599     !! <li> if \c mode_type is \c mono then the mode distribution shape
00600     !!      is a delta distribution at diameter \f$D_0\f$ and the
00601     !!      following parameters are:
00602     !!      <ul>
00603     !!      <li> \b radius (real, unit m): the radius \f$R_0\f$ of the
00604     !!           particles, so that \f$D_0 = 2 R_0\f$
00605     !!      </ul>
00606     !! </ul>
00607     !!
00608     !! Example:
00609     !! <pre>
00610     !! mode_name diesel          # mode name (descriptive only)
00611     !! mass_frac comp_diesel.dat # mass fractions in each aerosol particle
00612     !! num_conc 1.6e8            # particle number density (#/m^3)
00613     !! mode_type log_normal      # type of distribution
00614     !! geom_mean_diam 2.5e-8     # geometric mean diameter (m)
00615     !! log10_geom_std_dev 0.24   # log_10 of geometric standard deviation
00616     !! </pre>
00617     !!
00618     !! See also:
00619     !!   - \ref spec_file_format --- the input file text format
00620     !!   - \ref input_format_aero_dist --- the format for a complete
00621     !!     aerosol distribution with several modes
00622     !!   - \ref input_format_mass_frac --- the format for the mass
00623     !!     fractions file
00624 
00625     call spec_line_allocate(line)
00626     call spec_file_read_line(file, line, eof)
00627     if (.not. eof) then
00628        call spec_file_check_line_name(file, line, "mode_name")
00629        call spec_file_check_line_length(file, line, 1)
00630        call aero_mode_deallocate(aero_mode)
00631        call aero_mode_allocate_size(aero_mode, aero_data%n_spec)
00632        tmp_str = line%data(1) ! hack to avoid gfortran warning
00633        aero_mode%name = tmp_str(1:AERO_MODE_NAME_LEN)
00634        aero_mode%source = aero_data_source_by_name(aero_data, aero_mode%name)
00635        call spec_file_read_string(file, 'mass_frac', mass_frac_filename)
00636        call spec_file_open(mass_frac_filename, mass_frac_file)
00637        call spec_file_read_vol_frac(mass_frac_file, aero_data, &
00638             aero_mode%vol_frac)
00639        call spec_file_close(mass_frac_file)
00640        call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
00641        call spec_file_read_string(file, 'mode_type', mode_type)
00642        if (trim(mode_type) == 'log_normal') then
00643           aero_mode%type = AERO_MODE_TYPE_LOG_NORMAL
00644           call spec_file_read_real(file, 'geom_mean_diam', diam)
00645           aero_mode%char_radius = diam2rad(diam)
00646           call spec_file_read_real(file, 'log10_geom_std_dev', &
00647                aero_mode%log10_std_dev_radius)
00648        elseif (trim(mode_type) == 'exp') then
00649           aero_mode%type = AERO_MODE_TYPE_EXP
00650           call spec_file_read_real(file, 'diam_at_mean_vol', diam)
00651           aero_mode%char_radius = diam2rad(diam)
00652        elseif (trim(mode_type) == 'mono') then
00653           aero_mode%type = AERO_MODE_TYPE_MONO
00654           call spec_file_read_real(file, 'diam', diam)
00655           aero_mode%char_radius = diam2rad(diam)
00656        else
00657           call spec_file_die_msg(729472928, file, &
00658                "Unknown distribution mode type: " // trim(mode_type))
00659        end if
00660     end if
00661     call spec_line_deallocate(line)
00662 
00663   end subroutine spec_file_read_aero_mode
00664 
00665 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00666 
00667   !> Determines the number of bytes required to pack the given value.
00668   integer function pmc_mpi_pack_size_aero_mode(val)
00669 
00670     !> Value to pack.
00671     type(aero_mode_t), intent(in) :: val
00672 
00673     pmc_mpi_pack_size_aero_mode = &
00674          pmc_mpi_pack_size_string(val%name) &
00675          + pmc_mpi_pack_size_integer(val%type) &
00676          + pmc_mpi_pack_size_real(val%char_radius) &
00677          + pmc_mpi_pack_size_real(val%log10_std_dev_radius) &
00678          + pmc_mpi_pack_size_real(val%num_conc) &
00679          + pmc_mpi_pack_size_real_array(val%vol_frac) &
00680          + pmc_mpi_pack_size_integer(val%source)
00681 
00682   end function pmc_mpi_pack_size_aero_mode
00683 
00684 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00685 
00686   !> Packs the given value into the buffer, advancing position.
00687   subroutine pmc_mpi_pack_aero_mode(buffer, position, val)
00688 
00689     !> Memory buffer.
00690     character, intent(inout) :: buffer(:)
00691     !> Current buffer position.
00692     integer, intent(inout) :: position
00693     !> Value to pack.
00694     type(aero_mode_t), intent(in) :: val
00695 
00696 #ifdef PMC_USE_MPI
00697     integer :: prev_position
00698 
00699     prev_position = position
00700     call pmc_mpi_pack_string(buffer, position, val%name)
00701     call pmc_mpi_pack_integer(buffer, position, val%type)
00702     call pmc_mpi_pack_real(buffer, position, val%char_radius)
00703     call pmc_mpi_pack_real(buffer, position, val%log10_std_dev_radius)
00704     call pmc_mpi_pack_real(buffer, position, val%num_conc)
00705     call pmc_mpi_pack_real_array(buffer, position, val%vol_frac)
00706     call pmc_mpi_pack_integer(buffer, position, val%source)
00707     call assert(579699255, &
00708          position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
00709 #endif
00710 
00711   end subroutine pmc_mpi_pack_aero_mode
00712 
00713 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00714 
00715   !> Unpacks the given value from the buffer, advancing position.
00716   subroutine pmc_mpi_unpack_aero_mode(buffer, position, val)
00717 
00718     !> Memory buffer.
00719     character, intent(inout) :: buffer(:)
00720     !> Current buffer position.
00721     integer, intent(inout) :: position
00722     !> Value to pack.
00723     type(aero_mode_t), intent(inout) :: val
00724 
00725 #ifdef PMC_USE_MPI
00726     integer :: prev_position
00727 
00728     prev_position = position
00729     call pmc_mpi_unpack_string(buffer, position, val%name)
00730     call pmc_mpi_unpack_integer(buffer, position, val%type)
00731     call pmc_mpi_unpack_real(buffer, position, val%char_radius)
00732     call pmc_mpi_unpack_real(buffer, position, val%log10_std_dev_radius)
00733     call pmc_mpi_unpack_real(buffer, position, val%num_conc)
00734     call pmc_mpi_unpack_real_array(buffer, position, val%vol_frac)
00735     call pmc_mpi_unpack_integer(buffer, position, val%source)
00736     call assert(874467577, &
00737          position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
00738 #endif
00739 
00740   end subroutine pmc_mpi_unpack_aero_mode
00741 
00742 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00743 
00744 end module pmc_aero_mode