Scattering Truncation Corrections#

The study of aerosol optical properties, such as Single Scattering Albedo (SSA), plays a crucial role in understanding the effects of aerosols on climate and air quality. SSA, a measure of the fraction of light scattered by particles relative to the total light extinction (scattering plus absorption), is essential for assessing aerosol radiative forcing. Instruments like the Cavity Attenuated Phase Shift (CAPS) SSA instrument are pivotal in making these measurements with high accuracy and reliability. However, one of the critical challenges in accurately determining SSA, especially with the CAPS instrument, is the phenomenon of scattering truncation.

The Importance of Correcting for Scattering Truncation#

Scattering truncation occurs due to the limited angular range over which scattering measurements can be made, leading to an underestimation of the total scattering by aerosol particles. This limitation is particularly pronounced in instruments like the CAPS, where the design inherently restricts the detection of scattered light to a finite angular range (missing some back scatter and forward scattering light). The consequence of this truncation is a potential bias (low) in the measured SSA values, which can significantly affect the interpretation of aerosol optical properties and, by extension, their climatic and environmental impacts.

Size-Dependent Nature of Scattering Truncation#

The extent of scattering truncation is not uniform across all particle sizes; rather, it exhibits a pronounced size dependency. This variation arises because the scattering efficiency and the angular distribution of scattered light are functions of particle size relative to the wavelength of the incident light. Smaller particles tend to scatter light more isotropically, while larger particles preferentially scatter light in the forward direction. As a result, instruments with limited angular detection ranges may miss a substantial portion of the forward-scattered light from larger particles, leading to more significant truncation effects for these particle sizes.

Addressing the Challenge#

Correcting for scattering truncation is thus vital for ensuring the accuracy of SSA measurements and, by extension, our understanding of aerosol optical properties. This notebook focuses on methodologies for correcting scattering truncation effects in SSA measurements made using the CAPS instrument. By implementing these corrections, we aim to achieve more accurate and representative SSA values, enhancing our ability to model and predict aerosol impacts on atmospheric processes.

The following sections will delve into scattering truncation, explore its size-dependent characteristics, and introduce correction techniques tailored to the CAPS SSA instrument.

import numpy as np
import matplotlib.pyplot as plt

# particula imports
from particula.util import convert, distribution_discretization
from particula.data.process import scattering_truncation

Addressing Scattering Truncation for Single Particle Measurements#

In the quest to accurately assess aerosol optical properties using the Cavity Attenuated Phase Shift (CAPS) SSA instrument, accounting for scattering truncation is paramount. Scattering truncation arises due to the instrument’s limited angular range in capturing scattered light, necessitating corrections to obtain true scattering coefficients. This section introduces the methodology to correct for scattering truncation for a single aerosol particle using the scattering_truncation.trunc_mono function.

Implementing Truncation Correction with Python#

The process begins by defining the optical properties of the particle (refractive index and wavelength) and calculating both the ideal (untruncated) and the actual (truncated) scattering efficiencies. We then determine the correction factor needed to adjust for the truncation effect observed with the CAPS instrument.

Insights from Truncation Correction#

  • Ideal vs. Truncated Scattering Efficiency: The qsca_ideal output represents the scattering efficiency as if measured by a perfect instrument with no angular limitations. In contrast, qsca_trunc reflects the efficiency as captured by the CAPS instrument, which is inherently limited by scattering truncation.

  • Correction Factor: The trunc_corr value indicates the factor by which the measured (truncated) scattering efficiency needs to be adjusted to align with the ideal, untruncated scenario. This correction is crucial for accurately interpreting SSA measurements and understanding aerosol scattering behaviors.

This correction approach for a single particle lays the groundwork for more comprehensive analyses.

Note: The calibrated_trunc boolean applies an additional correction factor of 1.0224, so at a 150 nm diameter the truncation correction is 1. This factor is numerically determined and it is not clear where this error is introduced in the angular truncation calculation.

# Define the refractive index of the aerosol particle and the light wavelength
m_sphere = 1.5  # Refractive index of the particle
wavelength = 450.0  # Wavelength of incident light in nanometers (nm)

# Perform the truncation correction for a single particle of a given diameter
trunc_corr, z_axis, qsca_trunc, qsca_ideal, theta1, theta2 = scattering_truncation.trunc_mono(
    m_sphere=m_sphere,  # Refractive index of the particle
    wavelength=wavelength,  # Wavelength of incident light
    diameter=100,  # Diameter of the particle in nanometers
    full_output=True,  # Request full output for detailed analysis
    calibrated_trunc=True
)

# Output the calculated scattering efficiencies and truncation correction
# factor
print(f"Ideal Q_sca: {qsca_ideal} for a perfect instrument")
print(f"Truncated Q_sca: {qsca_trunc} for a truncated CAPS instrument")
print(f"Truncation correction factor: {trunc_corr}")
Ideal Q_sca: 0.5478174528225518 for a perfect instrument
Truncated Q_sca: 0.5363833503036401 for a truncated CAPS instrument
Truncation correction factor: 0.9988859309373502

Full Size Dependence of Truncation#

The size-dependent nature of scattering truncation necessitates a thorough understanding of how the phenomenon varies across different particle sizes. This section explores the full size dependence of scattering truncation and introduces the scattering_truncation.truncation_for_diameters function to correct for truncation effects across a range of particle sizes.

# Generating diameters in lin space from 50 nm to 1000 nm with 200 points
diameters = np.linspace(50, 2000, 100)

truncation_array = scattering_truncation.truncation_for_diameters(
    m_sphere=m_sphere,
    wavelength=wavelength,
    diameter_sizes=diameters,
    discretize=False,
    calibrated_trunc=True
)

# Plot the truncation correction factor as a function of particle diameter
plt.figure()
plt.plot(diameters, truncation_array, '.')
plt.xlabel('Particle diameter (nm)')
plt.ylabel('Truncation correction factor')
plt.title('Truncation correction factor as a function of particle diameter')
plt.show()
../../_images/be37a5ba8eb411b78988ce73dd19b27f500e33702cbbcc31dab78a3c89fb8c4e.png

Correcting Scattering Truncation 250 nm Distribution#

Aerosol populations in the atmosphere exhibit a broad range of particle sizes, each contributing distinctively to the overall optical scattering properties. This section demonstrates the application of the scattering_truncation.correction_for_distribution function, which facilitates the adjustment for scattering truncation across an aerosol size distributions measurement.

Generating a Representative Aerosol Distribution#

For this analysis, we focus on a distribution centered around 250 nm, a size range significant for many atmospheric aerosol studies. We simulate this aerosol population using a log-normal distribution, a common representation for atmospheric aerosol size distributions.

# Define the refractive index and wavelength for the aerosols
m_sphere = 1.5
wavelength = 450.0  # in nanometers

# Generate diameters using a log-spaced array for better representation
diameters = np.logspace(np.log10(50), np.log10(800), 100)  # From 50 nm to 800 nm

# Parameters for the log-normal size distribution
sigma = 1.4  # Geometric standard deviation
modes = 250  # Peak diameter for monomodal distribution
number_total = 1e3  # Total number of particles

# Create a log-normal size distribution
pdf_dist = distribution_discretization.discretize(
    interval=diameters,
    disttype="lognormal",
    gsigma=sigma,
    mode=modes,
    nparticles=number_total
).m

# Convert the size distribution to a probability mass function (PMF)
number_conc_pms = convert.distribution_convert_pdf_pms(
    x_array=diameters,
    distribution=pdf_dist,
    to_pdf=False
) * number_total

# Calculate the truncation correction for the entire size distribution
trunc_correction = scattering_truncation.correction_for_distribution(
    m_sphere=m_sphere,
    wavelength=wavelength,
    diameter_sizes=diameters,
    number_per_cm3=number_conc_pms,
    discretize=True
)

# Visualize the size distribution and truncation correction
plt.figure()
plt.plot(diameters, number_conc_pms, label="Size Distribution")
plt.xlabel('Particle Diameter (nm)')
plt.ylabel('Number Concentration (#/cm³)')
plt.title(f"Truncation Correction: {trunc_correction:.5f} for the entire distribution")
plt.legend()
plt.show()

print(f"Overall truncation correction for the distribution: {trunc_correction}")
../../_images/c7274dd6d208cb96b67bd81bfd40ccaeeffde66f63d454d6f9f7ed1f566281d7.png
Overall truncation correction for the distribution: 1.041606234228776

Correcting Scattering Truncation 1000 nm Distribution#

In addition to the 250 nm distribution, we also explore the correction of scattering truncation for a larger aerosol population centered around 1000 nm. This size range is particularly relevant for understanding the scattering properties of coarse-mode aerosols, which play a crucial role in atmospheric processes.

The truncation factor is much larger for the 1000 nm distribution, indicating the more significant impact of truncation for larger particles. This underscores the importance of accounting for scattering truncation across the full range of aerosol sizes to obtain accurate SSA measurements.

# Define the refractive index and wavelength for the aerosols
m_sphere = 1.5
wavelength = 450.0  # in nanometers

# Generate diameters using a log-spaced array for better representation
diameters = np.logspace(
    np.log10(200),
    np.log10(2500),
    100)  # From 50 nm to 800 nm

# Parameters for the log-normal size distribution
sigma = 1.4  # Geometric standard deviation
modes = 1000  # Peak diameter for monomodal distribution
number_total = 1e3  # Total number of particles

# Create a log-normal size distribution
pdf_dist = distribution_discretization.discretize(
    interval=diameters,
    disttype="lognormal",
    gsigma=sigma,
    mode=modes,
    nparticles=number_total
).m

# Convert the size distribution to a probability mass function (PMF)
number_conc_pms = convert.distribution_convert_pdf_pms(
    x_array=diameters,
    distribution=pdf_dist,
    to_pdf=False
) * number_total

# Calculate the truncation correction for the entire size distribution
trunc_correction = scattering_truncation.correction_for_distribution(
    m_sphere=m_sphere,
    wavelength=wavelength,
    diameter_sizes=diameters,
    number_per_cm3=number_conc_pms,
    discretize=True
)

# Visualize the size distribution and truncation correction
plt.figure()
plt.plot(diameters, number_conc_pms, label="Size Distribution")
plt.xlabel('Particle Diameter (nm)')
plt.ylabel('Number Concentration (#/cm³)')
plt.title(f"Truncation Correction: {trunc_correction:.5f} for the entire distribution")
plt.legend()
plt.show()

print(f"Overall truncation correction for the distribution: {trunc_correction}")
../../_images/ec2a2e9171c72324027d63ebd8e84cd99c19057e91e3114f9d0d793f03fda3e0.png
Overall truncation correction for the distribution: 1.1821845913839404

Summary#

Understanding the climatic and environmental impacts of aerosols necessitates accurate determinations of their optical properties, especially Single Scattering Albedo (SSA). The challenge of scattering truncation, however, complicates the precise measurement of SSA, a dilemma particularly pronounced with the use of the Cavity Attenuated Phase Shift (CAPS) SSA instrument. Throughout this notebook, we have explored and applied methodologies to correct for the effects of scattering truncation, delving into its size-dependent characteristics and the broader implications for accurately assessing SSA.

Addressing Humidified Aerosol Measurements#

Our exploration does extend to humidified aerosol measurements, where changes in particle refractive index and diameter under varying humidity levels introduce additional complexity. The scattering_truncation module, through its correction_for_humidified and correction_for_humidified_looped functions, offers robust solutions for incorporating these dynamic factors into SSA calculations, ensuring that measurements remain reflective of actual atmospheric conditions.

Leveraging Truncation Corrections for Broader Applications#

The adaptability of these correction techniques to a wide range of real-world aerosol distributions highlights their significance. Aerosols in the atmosphere exhibit a vast diversity in size, composition, and behavior, necessitating flexible and accurate measurement and analysis methods. By enhancing the precision of SSA and other optical property measurements, we lay a stronger foundation for atmospheric modeling, contributing to more reliable predictions of aerosol impacts on climate and environmental health.

help(scattering_truncation.trunc_mono)
Help on _lru_cache_wrapper in module particula.data.process.scattering_truncation:

trunc_mono(m_sphere: Union[complex, float], wavelength: float, diameter: float, full_output: bool = False, calibrated_trunc: bool = True, discretize: bool = True) -> Union[float, Tuple[float, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray]]
    Calculates the single scattering albedo (SSA) correction due to truncation
    for monodisperse aerosol measurements using the CAPS-PM-SSA instrument. The
    correction accounts for the incomplete angular range of scattering
    measurements due to the instrument's geometry.
    
    Parameters
    ----------
    m_sphere : Union[complex, float]
        Complex or real refractive index of the aerosol.
    wavelength : float
        Wavelength of light in nanometers used in the CAPS instrument.
    diameter : Union[float, np.float64]
        Diameter of the monodisperse aerosol in nanometers.
    full_output : bool, optional
        If True, additional details about the calculation are returned,
        including z-axis values, angles of integration, and both truncated and
        ideal scattering efficiencies.
    calibrated_trunc : bool, optional
        If True, applies a numberical calibration factor to the truncation
        correction, so 150 is 1. Default is True.
    discretize : bool, optional
        If True, discretizes the input parameters for potentially improved
        stability/performance in scattering function calculations. Can not be
        done for full_output=True
    
    Returns
    -------
    Union[float, Tuple[float, np.ndarray, np.ndarray, np.ndarray, np.ndarray,
    np.ndarray]]
    If fullOutput is False, returns only the truncation correction factor.
    If fullOutput is True, returns a tuple containing the truncation correction
    factor, z-axis positions, truncated scattering efficiency, ideal scattering
    efficiency, forward scattering angle, and backward scattering angle.
help(scattering_truncation.truncation_for_diameters)
Help on function truncation_for_diameters in module particula.data.process.scattering_truncation:

truncation_for_diameters(m_sphere: Union[complex, float], wavelength: float, diameter_sizes: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], discretize: bool = True, calibrated_trunc: bool = True) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]
    Calculates the truncation correction for an array of particle diameters
    given a specific refractive index and wavelength. This function is
    particularly useful for aerosol optical property measurements where
    truncation effects due to instrument geometry need to be accounted for.
    
    Parameters
    ----------
    m_sphere : Union[complex, float]
        The complex or real refractive index of the particles.
    wavelength : float
        The wavelength of the incident light in nanometers (nm).
    diameter_sizes : NDArray[np.float64]
        An array of particle diameters in nanometers (nm) for which the
        truncation correction will be calculated.
    discretize : bool, optional
        A flag indicating whether to discretize the input parameters for
        potentially improved calculation performance. Default is True.
    calibrated_trunc : bool, optional
        If True, applies a numberical calibration factor to the truncation
        correction, so 150 is 1. Default is True.
    
    Returns
    -------
    NDArray[np.float64]
        An array of truncation corrections corresponding to the input array of
        particle diameters.
help(scattering_truncation.correction_for_distribution)
Help on function correction_for_distribution in module particula.data.process.scattering_truncation:

correction_for_distribution(m_sphere: Union[complex, float], wavelength: float, diameter_sizes: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], number_per_cm3: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], discretize: bool = True) -> Union[float, numpy.float64]
    Calculates the correction factor for scattering measurements due to
    truncation effects in aerosol size distribution measurements. This
    correction factor is used to adjust the measured scattering
    coefficient, accounting for the limited angular range of the instrument.
    
    Parameters
    ----------
    m_sphere : Union[complex, float]
        The complex or real refractive index of the particles.
    wavelength : float
        The wavelength of the incident light in nanometers (nm).
    diameter_sizes : NDArray[np.float64]
        An array of particle diameters in nanometers (nm) corresponding to the
        size distribution.
    number_per_cm3 : NDArray[np.float64]
        An array of particle number concentrations (#/cm^3) for each diameter
        in the size distribution.
    discretize : bool, optional
        If True, the calculation will use discretized values for the
        refractive index, wavelength, and diameters to potentially improve
        computation performance. Default is True.
    
    Returns
    -------
    float
        The correction factor for scattering measurements. This factor is
        dimensionless and is used to correct the measured scattering
        coefficient for truncation effects, calculated as the ratio of
        the ideal (full angular range) to truncated scattering coefficient.
    
    Example
    -------
    b_sca_corrected = b_sca_measured * bsca_correction
help(scattering_truncation.correction_for_humidified)
Help on function correction_for_humidified in module particula.data.process.scattering_truncation:

correction_for_humidified(kappa: Union[float, numpy.float64], number_per_cm3: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], diameter: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], water_activity_sizer: numpy.float64, water_activity_sample: numpy.float64, refractive_index_dry: Union[complex, float] = 1.45, water_refractive_index: Union[complex, float] = 1.33, wavelength: float = 450, discretize: bool = True) -> numpy.float64
    Calculates the scattering correction for humidified aerosol measurements,
    accounting for water uptake by adjusting the aerosol's refractive index.
    This function requires the kappa values for the particles, which describe
    their hygroscopic growth.
    
    Parameters
    ----------
    kappa : Union[float, NDArray[np.float64]]
        Hygroscopicity parameter kappa, indicating the water uptake capability
        of the particles.
    number_per_cm3 : NDArray[np.float64]
        Number concentration of particles per cubic centimeter (#/cm³) for
        each size bin.
    diameter : NDArray[np.float64]
        Array of particle diameters in nanometers (nm).
    water_activity_sizer : NDArray[np.float64]
        Water activity (relative humidity/100) of the air sample used for
        sizing.
    water_activity_sample : NDArray[np.float64]
        Water activity (relative humidity/100) of the air sample in optical
        measurements.
    refractive_index_dry : Union[complex, float], optional
        Refractive index of the dry particles. Default is 1.45.
    water_refractive_index : Union[complex, float], optional
        Refractive index of water. Default is 1.33.
    wavelength : float, optional
        Wavelength of the incident light in nanometers (nm). Default is 450.
    discretize : bool, optional
        If True, calculation uses discretized values for refractive index,
        wavelength, and diameters to improve performance. Default is True.
    
    Returns
    -------
    np.float64
        A numpy array of scattering correction factors for each particle size
        in the distribution, to be applied to measured backscatter
        coefficients to account for truncation effects due to humidity.