PartMC 2.1.2
aero_dist.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_dist module.
00007 
00008 !> The aero_dist_t structure and associated subroutines.
00009 !!
00010 !! The initial size distributions are computed as number densities, so
00011 !! they can be used for both sectional and particle-resolved
00012 !! simulations. The routine dist_to_n() converts a number concentration
00013 !! distribution to an actual number of particles ready for a
00014 !! particle-resolved simulation.
00015 !!
00016 !! Initial distributions should be normalized so that <tt>sum(n_den) =
00017 !! 1/log_width</tt>.
00018 module pmc_aero_dist
00019 
00020   use pmc_bin_grid
00021   use pmc_util
00022   use pmc_constants
00023   use pmc_spec_file
00024   use pmc_aero_data
00025   use pmc_aero_mode
00026   use pmc_mpi
00027   use pmc_rand
00028 #ifdef PMC_USE_MPI
00029   use mpi
00030 #endif
00031 
00032   !> A complete aerosol distribution, consisting of several modes.
00033   type aero_dist_t
00034      !> Number of modes.
00035      integer :: n_mode
00036      !> Internally mixed modes [length \c n_mode].
00037      type(aero_mode_t), pointer :: mode(:)
00038   end type aero_dist_t
00039 
00040 contains
00041 
00042 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00043 
00044   !> Allocates an aero_dist.
00045   subroutine aero_dist_allocate(aero_dist)
00046 
00047     !> Aerosol distribution.
00048     type(aero_dist_t), intent(out) :: aero_dist
00049 
00050     integer :: i
00051 
00052     aero_dist%n_mode = 0
00053     allocate(aero_dist%mode(0))
00054 
00055   end subroutine aero_dist_allocate
00056 
00057 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00058 
00059   !> Allocates an aero_dist of the given size.
00060   subroutine aero_dist_allocate_size(aero_dist, n_mode, n_spec)
00061 
00062     !> Aerosol distribution.
00063     type(aero_dist_t), intent(out) :: aero_dist
00064     !> Number of modes.
00065     integer, intent(in) :: n_mode
00066     !> Number of species.
00067     integer, intent(in) :: n_spec
00068 
00069     integer :: i
00070 
00071     aero_dist%n_mode = n_mode
00072     allocate(aero_dist%mode(n_mode))
00073     do i = 1,n_mode
00074        call aero_mode_allocate_size(aero_dist%mode(i), n_spec)
00075     end do
00076 
00077   end subroutine aero_dist_allocate_size
00078 
00079 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00080 
00081   !> Free all storage.
00082   subroutine aero_dist_deallocate(aero_dist)
00083 
00084     !> Aerosol distribution.
00085     type(aero_dist_t), intent(inout) :: aero_dist
00086 
00087     integer :: i
00088 
00089     do i = 1,aero_dist%n_mode
00090        call aero_mode_deallocate(aero_dist%mode(i))
00091     end do
00092     deallocate(aero_dist%mode)
00093 
00094   end subroutine aero_dist_deallocate
00095 
00096 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00097 
00098   !> Copy an aero_dist.
00099   subroutine aero_dist_copy(aero_dist_from, aero_dist_to)
00100 
00101     !> Aero_dist original.
00102     type(aero_dist_t), intent(in) :: aero_dist_from
00103     !> Aero_dist copy.
00104     type(aero_dist_t), intent(inout) :: aero_dist_to
00105 
00106     integer :: n_spec, i
00107 
00108     if (aero_dist_from%n_mode > 0) then
00109        n_spec = size(aero_dist_from%mode(1)%vol_frac)
00110     else
00111        n_spec = 0
00112     end if
00113     call aero_dist_deallocate(aero_dist_to)
00114     call aero_dist_allocate_size(aero_dist_to, aero_dist_from%n_mode, n_spec)
00115     do i = 1,aero_dist_from%n_mode
00116        call aero_mode_copy(aero_dist_from%mode(i), aero_dist_to%mode(i))
00117     end do
00118 
00119   end subroutine aero_dist_copy
00120 
00121 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00122 
00123   !> Returns the total number concentration of a distribution. (#/m^3)
00124   real(kind=dp) function aero_dist_total_num_conc(aero_dist)
00125 
00126     !> Aerosol distribution.
00127     type(aero_dist_t), intent(in) :: aero_dist
00128 
00129     aero_dist_total_num_conc = sum(aero_dist%mode%num_conc)
00130 
00131   end function aero_dist_total_num_conc
00132 
00133 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00134 
00135   !> Returns the total weighted number concentration of a
00136   !> distribution. (#/m^3)
00137   real(kind=dp) function aero_dist_weighted_num_conc(aero_dist, &
00138        aero_weight)
00139 
00140     !> Aerosol distribution.
00141     type(aero_dist_t), intent(in) :: aero_dist
00142     !> Aerosol weight.
00143     type(aero_weight_t), intent(in) :: aero_weight
00144 
00145     integer :: i_mode
00146 
00147     aero_dist_weighted_num_conc = 0d0
00148     do i_mode = 1,aero_dist%n_mode
00149        aero_dist_weighted_num_conc = aero_dist_weighted_num_conc &
00150             + aero_mode_weighted_num_conc(aero_dist%mode(i_mode), &
00151             aero_weight)
00152     end do
00153 
00154   end function aero_dist_weighted_num_conc
00155 
00156 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00157 
00158   !> Return the binned number concentration for an aero_dist.
00159   subroutine aero_dist_num_conc(aero_dist, bin_grid, aero_data, &
00160        num_conc)
00161 
00162     !> Aero dist for which to compute number concentration.
00163     type(aero_dist_t), intent(in) :: aero_dist
00164     !> Bin grid.
00165     type(bin_grid_t), intent(in) :: bin_grid
00166     !> Aerosol data.
00167     type(aero_data_t), intent(in) :: aero_data
00168     !> Number concentration (#(ln(r))d(ln(r))).
00169     real(kind=dp), intent(out) :: num_conc(bin_grid%n_bin)
00170 
00171     integer :: i_mode
00172     real(kind=dp) :: mode_num_conc(size(num_conc, 1))
00173 
00174     num_conc = 0d0
00175     do i_mode = 1,aero_dist%n_mode
00176        call aero_mode_num_conc(aero_dist%mode(i_mode), bin_grid, &
00177             aero_data, mode_num_conc)
00178        num_conc = num_conc + mode_num_conc
00179     end do
00180 
00181   end subroutine aero_dist_num_conc
00182 
00183 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00184 
00185   !> Return the binned per-species volume concentration for an
00186   !> aero_dist.
00187   subroutine aero_dist_vol_conc(aero_dist, bin_grid, aero_data, &
00188        vol_conc)
00189 
00190     !> Aero dist for which to compute volume concentration.
00191     type(aero_dist_t), intent(in) :: aero_dist
00192     !> Bin grid.
00193     type(bin_grid_t), intent(in) :: bin_grid
00194     !> Aerosol data.
00195     type(aero_data_t), intent(in) :: aero_data
00196     !> Volume concentration (V(ln(r))d(ln(r))).
00197     real(kind=dp), intent(out) :: vol_conc(bin_grid%n_bin, aero_data%n_spec)
00198 
00199     integer :: i_mode
00200     real(kind=dp) :: mode_vol_conc(size(vol_conc, 1), size(vol_conc, 2))
00201 
00202     vol_conc = 0d0
00203     do i_mode = 1,aero_dist%n_mode
00204        call aero_mode_vol_conc(aero_dist%mode(i_mode), bin_grid, &
00205             aero_data, mode_vol_conc)
00206        vol_conc = vol_conc + mode_vol_conc
00207     end do
00208 
00209   end subroutine aero_dist_vol_conc
00210 
00211 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00212 
00213   !> Determine the current aero_dist and rate by interpolating at the
00214   !> current time with the lists of aero_dists and rates.
00215   subroutine aero_dist_interp_1d(aero_dist_list, time_list, &
00216          rate_list, time, aero_dist, rate)
00217 
00218     !> Gas states.
00219     type(aero_dist_t), intent(in) :: aero_dist_list(:)
00220     !> Times (s).
00221     real(kind=dp), intent(in) :: time_list(size(aero_dist_list))
00222     !> Rates (s^{-1}).
00223     real(kind=dp), intent(in) :: rate_list(size(aero_dist_list))
00224     !> Current time (s).
00225     real(kind=dp), intent(in) :: time
00226     !> Current gas state.
00227     type(aero_dist_t), intent(inout) :: aero_dist
00228     !> Current rate (s^{-1}).
00229     real(kind=dp), intent(out) :: rate
00230 
00231     integer :: n, p, n_bin, n_spec, i, i_new
00232     real(kind=dp) :: y, alpha
00233 
00234     n = size(aero_dist_list)
00235     p = find_1d(n, time_list, time)
00236     if (p == 0) then
00237        ! before the start, just use the first state and rate
00238        call aero_dist_copy(aero_dist_list(1), aero_dist)
00239        rate = rate_list(1)
00240     elseif (p == n) then
00241        ! after the end, just use the last state and rate
00242        call aero_dist_copy(aero_dist_list(n), aero_dist)
00243        rate = rate_list(n)
00244     else
00245        ! in the middle, use the previous dist
00246        call aero_dist_copy(aero_dist_list(p), aero_dist)
00247        rate = rate_list(p)
00248     end if
00249 
00250   end subroutine aero_dist_interp_1d
00251 
00252 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00253 
00254   !> Read continuous aerosol distribution composed of several modes.
00255   subroutine spec_file_read_aero_dist(file, aero_data, aero_dist)
00256 
00257     !> Spec file to read data from.
00258     type(spec_file_t), intent(inout) :: file
00259     !> Aero_data data.
00260     type(aero_data_t), intent(inout) :: aero_data
00261     !> Aerosol dist.
00262     type(aero_dist_t), intent(inout) :: aero_dist
00263 
00264     type(aero_mode_t), pointer :: new_aero_mode_list(:)
00265     type(aero_mode_t) :: aero_mode
00266     integer :: i, j
00267     logical :: eof
00268 
00269     ! note that the <p> is needed below to force correct paragraph
00270     ! breaking by doxygen
00271     
00272     !> \page input_format_aero_dist Input File Format: Aerosol Distribution
00273     !!
00274     !! <p>An aerosol distribution file consists of zero or more modes,
00275     !! each in the format described by \subpage input_format_aero_mode
00276     !!
00277     !! See also:
00278     !!   - \ref spec_file_format --- the input file text format
00279     !!   - \ref input_format_aero_mode --- the format for each mode
00280     !!     of an aerosol distribution
00281 
00282     call aero_dist_deallocate(aero_dist)
00283     call aero_dist_allocate(aero_dist)
00284     call aero_mode_allocate(aero_mode)
00285     call spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
00286     do while (.not. eof)
00287        aero_dist%n_mode = aero_dist%n_mode + 1
00288        allocate(new_aero_mode_list(aero_dist%n_mode))
00289        do i = 1,aero_dist%n_mode
00290           call aero_mode_allocate(new_aero_mode_list(i))
00291        end do
00292        call aero_mode_copy(aero_mode, &
00293             new_aero_mode_list(aero_dist%n_mode))
00294        do i = 1,(aero_dist%n_mode - 1)
00295           call aero_mode_copy(aero_dist%mode(i), &
00296                new_aero_mode_list(i))
00297           call aero_mode_deallocate(aero_dist%mode(i))
00298        end do
00299        deallocate(aero_dist%mode)
00300        aero_dist%mode => new_aero_mode_list
00301        nullify(new_aero_mode_list)
00302        call spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
00303     end do
00304     call aero_mode_deallocate(aero_mode)
00305 
00306   end subroutine spec_file_read_aero_dist
00307 
00308 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00309 
00310   !> Read an array of aero_dists with associated times and rates from
00311   !> the given file.
00312   subroutine spec_file_read_aero_dists_times_rates(file, aero_data, &
00313        bin_grid, times, rates, aero_dists)
00314 
00315     !> Spec file to read data from.
00316     type(spec_file_t), intent(inout) :: file
00317     !> Aero data.
00318     type(aero_data_t), intent(inout) :: aero_data
00319     !> Bin grid.
00320     type(bin_grid_t), intent(in) :: bin_grid
00321     !> Times (s).
00322     real(kind=dp), pointer :: times(:)
00323     !> Rates (s^{-1}).
00324     real(kind=dp), pointer :: rates(:)
00325     !> Aero dists.
00326     type(aero_dist_t), pointer :: aero_dists(:)
00327 
00328     type(spec_line_t) :: aero_dist_line
00329     type(spec_file_t) :: aero_dist_file
00330     integer :: n_time, i_time
00331     character(len=SPEC_LINE_MAX_VAR_LEN), pointer :: names(:)
00332     real(kind=dp), pointer :: data(:,:)
00333 
00334     !> \page input_format_aero_dist_profile Input File Format: Aerosol Distribution Profile
00335     !!
00336     !! An aerosol distribution profile input file must consist of
00337     !! three lines:
00338     !! - the first line must begin with \c time and should be followed
00339     !!   by \f$N\f$ space-separated real scalars, giving the times (in
00340     !!   s after the start of the simulation) of the aerosol
00341     !!   distrbution set points --- the times must be in increasing
00342     !!   order
00343     !! - the second line must begin with \c rate and should be
00344     !!   followed by \f$N\f$ space-separated real scalars, giving the
00345     !!   scalings at the corresponding times
00346     !! - the third line must begin with \c dist and should be followed
00347     !!   by \f$N\f$ space-separated filenames, each specifying an
00348     !!   aerosol distribution in the format \ref input_format_aero_dist
00349     !!   at the corresponding time
00350     !!
00351     !! The units of the \c rate line depend on the type of aerosol
00352     !! distribution profile:
00353     !! - Emissions aerosol profiles have rates with units m/s ---
00354     !!   the aerosol distribution number concentrations are multiplied
00355     !!   by the rate to give an emission rate with unit #/(m^2 s)
00356     !!   which is then divided by the current mixing layer height
00357     !!   to give a per-volume emission rate.
00358     !! - Background aerosol profiles have rates with units \f$\rm
00359     !!   s^{-1}\f$, which is the dilution rate between the background
00360     !!   and the simulated air parcel. That is, if the simulated
00361     !!   number concentration is \f$N\f$ and the background number
00362     !!   concentration is \f$N_{\rm back}\f$, then dilution is modeled
00363     !!   as \f$\dot{N} = r N_{\rm back} - r N\f$, where \f$r\f$ is the
00364     !!   rate.
00365     !!
00366     !! Between the specified times the aerosol profile is interpolated
00367     !! step-wise and kept constant at its last value. That is, if the
00368     !! times are \f$t_i\f$, the rates are \f$r_i\f$, and the aerosol
00369     !! distributions are \f$a_i\f$ (all with \f$i = 1,\ldots,n\f$),
00370     !! then between times \f$t_i\f$ and \f$t_{i+1}\f$ the aerosol
00371     !! state is constant at \f$r_i a_i\f$. Before time \f$t_1\f$ the
00372     !! aerosol state is \f$r_1 a_1\f$, while after time \f$t_n\f$ it
00373     !! is \f$r_n a_n\f$.
00374     !!
00375     !! Example: an emissions aerosol profile could be:
00376     !! <pre>
00377     !! time  0          600        1800       # time (in s) after simulation start
00378     !! rate  1          0.5        1          # scaling factor in m/s
00379     !! dist  dist1.dat  dist2.dat  dist3.dat  # aerosol distribution files
00380     !! </pre>
00381     !! Here the emissions between 0&nbsp;min and 10&nbsp;min are given
00382     !! by <tt>dist1.dat</tt> (with the number concentration
00383     !! interpreted as having units 1/(m^2 s)), the emissions between
00384     !! 10&nbsp;min and 30&nbsp;min are given by <tt>dist2.dat</tt>
00385     !! (scaled by 0.5), while the emissions after 30&nbsp;min are
00386     !! given by <tt>dist3.dat</tt>.
00387     !!
00388     !! See also:
00389     !!   - \ref spec_file_format --- the input file text format
00390     !!   - \ref input_format_aero_data --- the aerosol species list
00391     !!     and material data
00392     !!   - \ref input_format_aero_dist --- the format of the
00393     !!     instantaneous aerosol distribution files
00394 
00395     ! read the data from the file
00396     allocate(names(0))
00397     allocate(data(0,0))
00398     call spec_file_read_real_named_array(file, 2, names, data)
00399     call spec_line_allocate(aero_dist_line)
00400     call spec_file_read_line_no_eof(file, aero_dist_line)
00401     call spec_file_check_line_name(file, aero_dist_line, "dist")
00402     call spec_file_check_line_length(file, aero_dist_line, size(data, 2))
00403 
00404     ! check the data size
00405     if (trim(names(1)) /= 'time') then
00406        call die_msg(570205795, 'row 1 in ' // trim(file%name) &
00407             // ' must start with: time not: ' // trim(names(1)))
00408     end if
00409     if (trim(names(2)) /= 'rate') then
00410        call die_msg(221270915, 'row 2 in ' // trim(file%name) &
00411             // ' must start with: rate not: ' // trim(names(1)))
00412     end if
00413     n_time = size(data, 2)
00414     if (n_time < 1) then
00415        call die_msg(457229710, 'each line in ' // trim(file%name) &
00416             // ' must contain at least one data value')
00417     end if
00418 
00419     ! copy over the data
00420     do i_time = 1,size(aero_dists)
00421        call aero_dist_deallocate(aero_dists(i_time))
00422     end do
00423     deallocate(aero_dists)
00424     deallocate(times)
00425     deallocate(rates)
00426     allocate(aero_dists(n_time))
00427     allocate(times(n_time))
00428     allocate(rates(n_time))
00429     do i_time = 1,n_time
00430        call spec_file_open(aero_dist_line%data(i_time), aero_dist_file)
00431        call aero_dist_allocate(aero_dists(i_time))
00432        call spec_file_read_aero_dist(aero_dist_file, &
00433             aero_data, aero_dists(i_time))
00434        call spec_file_close(aero_dist_file)
00435        times(i_time) = data(1,i_time)
00436        rates(i_time) = data(2,i_time)
00437     end do
00438     deallocate(names)
00439     deallocate(data)
00440     call spec_line_deallocate(aero_dist_line)
00441 
00442   end subroutine spec_file_read_aero_dists_times_rates
00443 
00444 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00445 
00446   !> Determines the number of bytes required to pack the given value.
00447   integer function pmc_mpi_pack_size_aero_dist(val)
00448 
00449     !> Value to pack.
00450     type(aero_dist_t), intent(in) :: val
00451 
00452     integer :: i, total_size
00453 
00454     total_size = pmc_mpi_pack_size_integer(val%n_mode)
00455     do i = 1,size(val%mode)
00456        total_size = total_size + pmc_mpi_pack_size_aero_mode(val%mode(i))
00457     end do
00458     pmc_mpi_pack_size_aero_dist = total_size
00459 
00460   end function pmc_mpi_pack_size_aero_dist
00461 
00462 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00463 
00464   !> Packs the given value into the buffer, advancing position.
00465   subroutine pmc_mpi_pack_aero_dist(buffer, position, val)
00466 
00467     !> Memory buffer.
00468     character, intent(inout) :: buffer(:)
00469     !> Current buffer position.
00470     integer, intent(inout) :: position
00471     !> Value to pack.
00472     type(aero_dist_t), intent(in) :: val
00473 
00474 #ifdef PMC_USE_MPI
00475     integer :: prev_position, i
00476 
00477     prev_position = position
00478     call pmc_mpi_pack_integer(buffer, position, val%n_mode)
00479     do i = 1,size(val%mode)
00480        call pmc_mpi_pack_aero_mode(buffer, position, val%mode(i))
00481     end do
00482     call assert(440557910, &
00483          position - prev_position <= pmc_mpi_pack_size_aero_dist(val))
00484 #endif
00485 
00486   end subroutine pmc_mpi_pack_aero_dist
00487 
00488 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00489 
00490   !> Unpacks the given value from the buffer, advancing position.
00491   subroutine pmc_mpi_unpack_aero_dist(buffer, position, val)
00492 
00493     !> Memory buffer.
00494     character, intent(inout) :: buffer(:)
00495     !> Current buffer position.
00496     integer, intent(inout) :: position
00497     !> Value to pack.
00498     type(aero_dist_t), intent(inout) :: val
00499 
00500 #ifdef PMC_USE_MPI
00501     integer :: prev_position, i
00502 
00503     call aero_dist_deallocate(val)
00504     prev_position = position
00505     call pmc_mpi_unpack_integer(buffer, position, val%n_mode)
00506     allocate(val%mode(val%n_mode))
00507     do i = 1,size(val%mode)
00508        call aero_mode_allocate(val%mode(i))
00509        call pmc_mpi_unpack_aero_mode(buffer, position, val%mode(i))
00510     end do
00511     call assert(742535268, &
00512          position - prev_position <= pmc_mpi_pack_size_aero_dist(val))
00513 #endif
00514 
00515   end subroutine pmc_mpi_unpack_aero_dist
00516 
00517 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00518 
00519 end module pmc_aero_dist