PartMC 2.1.2
aero_weight.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2010-2011 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_weight module.
00007 
00008 !> The aero_weight_t structure and associated subroutines.
00009 module pmc_aero_weight
00010 
00011   use pmc_util
00012   use pmc_constants
00013   use pmc_rand
00014   use pmc_spec_file
00015   use pmc_netcdf
00016   use pmc_mpi
00017 #ifdef PMC_USE_MPI
00018   use mpi
00019 #endif
00020 
00021   !> Type code for an undefined or invalid weighting.
00022   integer, parameter :: AERO_WEIGHT_TYPE_INVALID  = 0
00023   !> Type code for no (or flat) weighting.
00024   integer, parameter :: AERO_WEIGHT_TYPE_NONE     = 1
00025   !> Type code for power function weighting.
00026   integer, parameter :: AERO_WEIGHT_TYPE_POWER    = 2
00027   !> Type code for MFA weighting.
00028   integer, parameter :: AERO_WEIGHT_TYPE_MFA      = 3
00029 
00030   !> An aerosol size distribution weighting function.
00031   type aero_weight_t
00032      !> Weight type (given by module constants).
00033      integer :: type
00034      !> Reference radius at which the weight is 1 (hard-coded at present).
00035      real(kind=dp) :: ref_radius
00036      !> Exponent for "power" weight.
00037      real(kind=dp) :: exponent
00038   end type aero_weight_t
00039 
00040 contains
00041 
00042 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00043 
00044   !> Allocates an aero_weight.
00045   subroutine aero_weight_allocate(aero_weight)
00046 
00047     !> Aerosol weight.
00048     type(aero_weight_t), intent(out) :: aero_weight
00049 
00050     call aero_weight_zero(aero_weight)
00051 
00052   end subroutine aero_weight_allocate
00053 
00054 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00055 
00056   !> Allocates an aero_weight of the given size.
00057   subroutine aero_weight_allocate_size(aero_weight)
00058 
00059     !> Aerosol weight.
00060     type(aero_weight_t), intent(out) :: aero_weight
00061 
00062     call aero_weight_zero(aero_weight)
00063 
00064   end subroutine aero_weight_allocate_size
00065 
00066 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00067 
00068   !> Free all storage.
00069   subroutine aero_weight_deallocate(aero_weight)
00070 
00071     !> Aerosol weight.
00072     type(aero_weight_t), intent(in) :: aero_weight
00073 
00074   end subroutine aero_weight_deallocate
00075 
00076 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00077 
00078   !> Zeros the contents of the \c aero_weight.
00079   subroutine aero_weight_zero(aero_weight)
00080 
00081     !> Aerosol weight.
00082     type(aero_weight_t), intent(inout) :: aero_weight
00083 
00084     aero_weight%type = AERO_WEIGHT_TYPE_INVALID
00085     aero_weight%ref_radius = 0d0
00086     aero_weight%exponent = 0d0
00087 
00088   end subroutine aero_weight_zero
00089 
00090 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00091 
00092   !> Copy an aero_weight.
00093   subroutine aero_weight_copy(aero_weight_from, aero_weight_to)
00094 
00095     !> Aerosol weight original.
00096     type(aero_weight_t), intent(in) :: aero_weight_from
00097     !> Aerosol weight copy.
00098     type(aero_weight_t), intent(inout) :: aero_weight_to
00099 
00100     aero_weight_to%type = aero_weight_from%type
00101     aero_weight_to%ref_radius = aero_weight_from%ref_radius
00102     aero_weight_to%exponent = aero_weight_from%exponent
00103 
00104   end subroutine aero_weight_copy
00105 
00106 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00107 
00108   !> Give the value of an aero_weight at a specific radius.
00109   real(kind=dp) function aero_weight_value(aero_weight, radius)
00110 
00111     !> Aerosol weight.
00112     type(aero_weight_t), intent(in) :: aero_weight
00113     !> Radius to compute weight at (m).
00114     real(kind=dp), intent(in) :: radius
00115 
00116     if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then
00117        aero_weight_value = 1d0
00118     elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) &
00119          .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then
00120        aero_weight_value &
00121             = (radius / aero_weight%ref_radius)**aero_weight%exponent
00122     else
00123        call die_msg(700421478, "unknown aero_weight type: " &
00124             // trim(integer_to_string(aero_weight%type)))
00125     end if
00126 
00127   end function aero_weight_value
00128 
00129 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00130 
00131   !> Read an aero_weight from a spec file.
00132   subroutine spec_file_read_aero_weight(file, aero_weight)
00133 
00134     !> Spec file.
00135     type(spec_file_t), intent(inout) :: file
00136     !> Aerosol weight.
00137     type(aero_weight_t), intent(inout) :: aero_weight
00138 
00139     character(len=SPEC_LINE_MAX_VAR_LEN) :: weight_type
00140 
00141     !> \page input_format_aero_weight Input File Format: Aerosol Weighting Function
00142     !!
00143     !! For efficiency the aerosol population can be weighted so that
00144     !! the true number distribution \f$n(D)\f$ is given by
00145     !! \f[ n(D) = w(D) c(D) \f]
00146     !! where \f$w(D)\f$ is a fixed weighting function, \f$c(D)\f$ is
00147     !! the computational (simulated) number distribution, and \f$D\f$
00148     !! is the diameter. Thus a large value of \f$w(D)\f$ means that
00149     !! relatively few computational particles are used at diameter
00150     !! \f$D\f$, while a small value of \f$w(D)\f$ means that
00151     !! relatively many computational particles will be used at that
00152     !! diameter.
00153     !!
00154     !! The aerosol weighting function is specified by the parameters:
00155     !!   - \b weight (string): the type of weighting function --- must
00156     !!     be one of: "none" for no weighting (\f$w(D) = 1\f$);
00157     !!     "power" for a power-law weighting (\f$w(D) = D^\alpha\f$),
00158     !!     or "mfa" for the mass flow algorithm weighting (\f$w(D) =
00159     !!     D^{-3}\f$ with dependent coagulation particle removal)
00160     !!   - if the \c weight is \c power then the next parameter is:
00161     !!     - \b exponent (real, dimensionless): the exponent
00162     !!       \f$\alpha\f$ in the power law relationship --- setting
00163     !!       the \c exponent to 0 is equivalent to no weighting, while
00164     !!       setting the exponent negative uses more computational
00165     !!       particles at larger diameters and setting the exponent
00166     !!       positive uses more computational particles at smaller
00167     !!       diameters; in practice exponents between 0 and -3 are
00168     !!       most useful
00169     !!
00170     !! See also:
00171     !!   - \ref spec_file_format --- the input file text format
00172 
00173     call spec_file_read_string(file, 'weight', weight_type)
00174     if (trim(weight_type) == 'none') then
00175        aero_weight%type = AERO_WEIGHT_TYPE_NONE
00176     elseif (trim(weight_type) == 'power') then
00177        aero_weight%type = AERO_WEIGHT_TYPE_POWER
00178        aero_weight%ref_radius = 1d0
00179        call spec_file_read_real(file, 'exponent', aero_weight%exponent)
00180     elseif (trim(weight_type) == 'mfa') then
00181        aero_weight%type = AERO_WEIGHT_TYPE_MFA
00182        aero_weight%ref_radius = 1d0
00183        aero_weight%exponent = -3d0
00184     else
00185        call spec_file_die_msg(456342050, file, "unknown weight_type: " &
00186             // trim(weight_type))
00187     end if
00188 
00189   end subroutine spec_file_read_aero_weight
00190 
00191 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00192 
00193   !> Determines the number of bytes required to pack the given value.
00194   integer function pmc_mpi_pack_size_aero_weight(val)
00195 
00196     !> Value to pack.
00197     type(aero_weight_t), intent(in) :: val
00198 
00199     pmc_mpi_pack_size_aero_weight = &
00200          pmc_mpi_pack_size_integer(val%type) &
00201          + pmc_mpi_pack_size_real(val%ref_radius) &
00202          + pmc_mpi_pack_size_real(val%exponent)
00203 
00204   end function pmc_mpi_pack_size_aero_weight
00205 
00206 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00207 
00208   !> Packs the given value into the buffer, advancing position.
00209   subroutine pmc_mpi_pack_aero_weight(buffer, position, val)
00210 
00211     !> Memory buffer.
00212     character, intent(inout) :: buffer(:)
00213     !> Current buffer position.
00214     integer, intent(inout) :: position
00215     !> Value to pack.
00216     type(aero_weight_t), intent(in) :: val
00217 
00218 #ifdef PMC_USE_MPI
00219     integer :: prev_position
00220 
00221     prev_position = position
00222     call pmc_mpi_pack_integer(buffer, position, val%type)
00223     call pmc_mpi_pack_real(buffer, position, val%ref_radius)
00224     call pmc_mpi_pack_real(buffer, position, val%exponent)
00225     call assert(579699255, &
00226          position - prev_position <= pmc_mpi_pack_size_aero_weight(val))
00227 #endif
00228 
00229   end subroutine pmc_mpi_pack_aero_weight
00230 
00231 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00232 
00233   !> Unpacks the given value from the buffer, advancing position.
00234   subroutine pmc_mpi_unpack_aero_weight(buffer, position, val)
00235 
00236     !> Memory buffer.
00237     character, intent(inout) :: buffer(:)
00238     !> Current buffer position.
00239     integer, intent(inout) :: position
00240     !> Value to pack.
00241     type(aero_weight_t), intent(inout) :: val
00242 
00243 #ifdef PMC_USE_MPI
00244     integer :: prev_position
00245 
00246     prev_position = position
00247     call pmc_mpi_unpack_integer(buffer, position, val%type)
00248     call pmc_mpi_unpack_real(buffer, position, val%ref_radius)
00249     call pmc_mpi_unpack_real(buffer, position, val%exponent)
00250     call assert(874467577, &
00251          position - prev_position <= pmc_mpi_pack_size_aero_weight(val))
00252 #endif
00253 
00254   end subroutine pmc_mpi_unpack_aero_weight
00255 
00256 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00257 
00258   !> Write full state.
00259   subroutine aero_weight_output_netcdf(aero_weight, ncid)
00260 
00261     !> Aero weight to write.
00262     type(aero_weight_t), intent(in) :: aero_weight
00263     !> NetCDF file ID, in data mode.
00264     integer, intent(in) :: ncid
00265 
00266     !> \page output_format_aero_weight Output File Format: Aerosol Weighting Function
00267     !!
00268     !! The aerosol weighting function NetCDF variables are:
00269     !!   - \b weight_type (dimensionless integer): the type of the
00270     !!     weighting function, with 0 = invalid weight, 1 = no weight
00271     !!     (\f$w(D) = 1\f$), 2 = power weight (\f$w(D) =
00272     !!     (D/D_0)^\alpha\f$), 3 = MFA weight (\f$w(D) =
00273     !!     (D/D_0)^{-3}\f$)
00274     !!   - \b weight_exponent (dimensionless): the exponent
00275     !!     \f$\alpha\f$ for the power \c weight_type, set to -3
00276     !!     for the MFA \c weight_type, and zero otherwise
00277     !!
00278     !! See also:
00279     !!   - \ref input_format_aero_weight --- the corresponding input format
00280 
00281     call pmc_nc_write_integer(ncid, aero_weight%type, "weight_type", &
00282          description="type of aerosol weighting function: 0 = invalid, " &
00283          // "1 = none (w(D) = 1), 2 = power (w(D) = (D/D_0)^alpha), " &
00284          // "3 = MFA (mass flow) (w(D) = (D/D_0)^(-3))")
00285     call pmc_nc_write_real(ncid, aero_weight%exponent, "weight_exponent", &
00286          unit="1", &
00287          description="exponent alpha for the power weight_type, " &
00288          // "set to -3 for MFA, and zero otherwise")
00289 
00290   end subroutine aero_weight_output_netcdf
00291 
00292 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00293 
00294   !> Read full state.
00295   subroutine aero_weight_input_netcdf(aero_weight, ncid)
00296 
00297     !> Environment state to read.
00298     type(aero_weight_t), intent(inout) :: aero_weight
00299     !> NetCDF file ID, in data mode.
00300     integer, intent(in) :: ncid
00301 
00302     call pmc_nc_read_integer(ncid, aero_weight%type, "weight_type")
00303     aero_weight%ref_radius = 1d0
00304     call pmc_nc_read_real(ncid, aero_weight%exponent, "weight_exponent")
00305 
00306   end subroutine aero_weight_input_netcdf
00307 
00308 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00309 
00310 end module pmc_aero_weight