aerosol package¶
aerosol.functions module¶
Aerosol number-size distribution¶
In the below functions of this package aerosol number-size distribution is assumed to be a pandas DataFrame where
- indexpandas DatetimeIndex
timestamps
- columnsfloat
size bin geomean diameters, m
- valuesfloat
normalized concentration, dN/dlogDp, cm-3
- aerosol.functions.air_density(temp, pres)[source]¶
Calculate air density
- Parameters:
- tempfloat or series of lenght n
absolute temperature (K)
- presfloat or series of length n
absolute pressure (Pa)
- Returns:
- float or series of length n
air density (kg/m3)
- aerosol.functions.air_viscosity(temp)[source]¶
Calculate air viscosity using Enskog-Chapman theory
- Parameters:
- tempfloat or series of length n
air temperature, unit: K
- Returns:
- float or series of length n
viscosity of air, unit: m2 s-1
- aerosol.functions.beta(dp, temp, pres, diffusivity, molar_mass)[source]¶
Calculate Fuchs Sutugin correction factor
Sutugin et al. (1971): https://doi.org/10.1016/0021-8502(71)90061-9
- Parameters:
- dpfloat or series of lenght m
aerosol particle diameter(s), unit: m
- tempfloat or series of lenght n
temperature, unit: K
- presfloat or series of lenght n
pressure, unit: Pa
- diffusivityfloat or series of length n
diffusivity of the gas that is condensing, unit: m2/s
- molar_massfloat
molar mass of the condensing gas, unit: g/mol
- Returns:
- float or dataframe of shape (n,m)
Fuchs Sutugin correction factor for each particle diameter and temperature/pressure unit: m2/s
- aerosol.functions.binary_diffusivity(temp, pres, Ma, Mb, Va, Vb)[source]¶
Binary diffusivity in a mixture of gases a and b
Fuller et al. (1966): https://doi.org/10.1021/ie50677a007
- Parameters:
- tempfloat or series of length n
temperature, unit: K
- presfloat or series of length n
pressure, unit: Pa
- Mafloat
relative molecular mass of gas a, unit: dimensionless
- Mbfloat
relative molecular mass of gas b, unit: dimensionless
- Vafloat
diffusion volume of gas a, unit: dimensionless
- Vbfloat
diffusion volume of gas b, unit: dimensionless
- Returns:
- float or series of length n
binary diffusivity, unit: m2 s-1
- aerosol.functions.calc_bin_edges(dp)[source]¶
Calculate bin edges given bin centers
- Parameters:
- dppandas series of lenght n
bin center diameters
- Returns:
- pandas series of lenght n+1
log bin edges
- aerosol.functions.calc_coags(df, dp, temp=293.15, pres=101325.0, dp_start=None)[source]¶
Calculate coagulation sink
Kulmala et al (2012): doi:10.1038/nprot.2012.091
- Parameters:
- dfdataframe
Aerosol number size distribution
- dpfloat or series of length m
Particle diameter(s) for which you want to calculate the CoagS, unit: m
- tempfloat or series indexed by DatetimeIndex
Ambient temperature corresponding to the data, unit: K If single value given it is used for all data
- presfloat or series indexed by DatetimeIndex
Ambient pressure corresponding to the data, unit: Pa If single value given it is used for all data
- dp_startfloat or None
The smallest size that you consider as part of the coagulation sink If None (default) then the smallest size is from dp
- Returns:
- float or dataframe
Coagulation sink for the given diamater(s), unit: s-1
- aerosol.functions.calc_conc(df, dmin, dmax, frac=0.5)[source]¶
Calculate particle number concentration from aerosol number-size distribution by adding whole bins
- Parameters:
- dfdataframe
Aerosol number-size distribution
- dminfloat or series of length n
Size range lower diameter(s), unit: m
- dmaxfloat or series of length n
Size range upper diameter(s), unit: m
- fracfloat
Minimum fraction of available data when calculating a concentration point
- Returns:
- dataframe
Number concentration in the given size range(s), unit: cm-3
- aerosol.functions.calc_conc_interp(df, dmin, dmax, threshold=0.0)[source]¶
Calculate particle number concentration from aerosol number-size distribution using integration and linear approximation
- Parameters:
- dfdataframe
Aerosol number-size distribution
- dminfloat or series of length n
Size range lower diameter(s), unit: m
- dmaxfloat or series of length n
Size range upper diameter(s), unit: m
- thresholdfloat
fraction of nans accepted per row
- Returns:
- dataframe
Number concentration in the given size range(s), unit: cm-3 in the index of the original dataframe
- aerosol.functions.calc_cs(df, temp=293.15, pres=101325.0)[source]¶
Calculate condensation sink, assuming that the condensing gas is sulfuric acid in air with aerosol particles.
Kulmala et al (2012): doi:10.1038/nprot.2012.091
- Parameters:
- dfpandas.DataFrame
aerosol number size distribution (dN/dlogDp)
- temppandas.Series or float
Ambient temperature corresponding to the data, unit: K If single value given it is used for all data
- prespandas.Series or float
Ambient pressure corresponding to the data, unit: Pa If single value given it is used for all data
- Returns:
- pandas.Series
condensation sink, unit: s-1
- aerosol.functions.calc_formation_rate(df, dp1, dp2, gr, coags=None, temp=293.15, pres=101325.0)[source]¶
Calculate particle formation rate
Kulmala et al (2012): doi:10.1038/nprot.2012.091
- Parameters:
- dfdataframe
Aerosol particle number size distribution
- dp1float or series of length n
Lower diameter of the size range(s) Unit m
- dp2float or series of length n
Upper diameter of the size range(s) Unit m
- grfloat or series of length n
Growth rate for particles out of the size range(s), unit nm h-1
- coagsdataframe or None
Coagulation sink of representative diameter, If None it is calculated from the number size distribution Unit: s-1
- tempfloat or series
Ambient temperature corresponding to the data, unit K
- presfloat or series
Ambient pressure corresponding to the data unit Pa
- Returns:
- dataframe
Particle formation rate(s) for the diameter range(s) Unit cm3 s-1
- aerosol.functions.calc_ion_formation_rate(df_particles, df_negions, df_posions, dp1, dp2, gr_negions, gr_posions, temp=293.15, pres=101325.0)[source]¶
Calculate ion formation rate
Kulmala et al (2012): doi:10.1038/nprot.2012.091
- Parameters:
- df_particlesdataframe
Aerosol particle number size distribution
- df_negionsdataframe
Negative ion number size distribution
- df_posionsdataframe
Positive ion number size distribution
- dp1float or series of length n
Lower diameter of the size range(s), unit: m
- dp2float or series of length n
Upper diameter of the size range(s), unit: m
- gr_negionsfloat or series of length n
Growth rate for negative ions out of the size range(s), unit: nm h-1
- gr_posionsfloat or series of length n
Growth rate for positive ions out of the size range(s), unit: nm h-1
- tempfloat or series
Ambient temperature corresponding to the data, unit: K
- presor series
Ambient pressure corresponding to the data, unit: Pa
- Returns:
- dataframe
Negative ion formation rate(s), unit : cm3 s-1
- dataframe
Positive ion formation rate(s), unit: cm3 s-1
- aerosol.functions.calc_ion_production_rate(df_ions, df_particles, temp=293.15, pres=101325.0)[source]¶
Calculate the ion production rate from measurements
- Parameters:
- df_ionsdataframe of shape (n,m)
negative or positive ion number size distribution unit cm-3
- df_particlesdataframe of shape (n,m)
particle number size distribution unit cm-3
- tempfloat or series of length n
ambient temperature unit K
- presfloat or series of length n
ambient pressure unit Pa
- Returns:
- series of lenght n
ion production rate in cm-3 s-1
- aerosol.functions.calc_ldsa(df)[source]¶
Calculate total LDSA from number size distribution data
ICRP, 1994. Human respiratory tract model for radiological protection. A report of a task group of the international commission on radiological protection. Ann. ICRP 24 (1-3), 1-482
- Parameters:
- dfpandas.DataFrame
Aerosol number-size distribution
- Returns:
- pandas.DataFrame
Total LDSA for alveoli (“al”), trachea/bronchi (“tb”) head-airways (“ha”) and all combiend (“tot”) unit: um2 cm-3
- aerosol.functions.calc_lung_df(dp)[source]¶
Calculate lung deposition fractions for particle diameters
ICRP, 1994. Human respiratory tract model for radiological protection. A report of a task group of the international commission on radiological protection. Ann. ICRP 24 (1-3), 1-482
- Parameters:
- dppandas.Series
aerosol particle diameters unit: m
- Returns:
- pandas.DataFrame
Lung deposition fractions for alveoli (“DF_al”), trachea/bronchi (“DF_tb”) head-airways (“DF_ha”) and all combiend (“DF_tot”)
- aerosol.functions.calc_tube_residence_time(tube_diam, tube_length, flowrate)[source]¶
Calculate residence time in a circular tube
- Parameters:
- tube_diamfloat or series of length m
Inner diameter of the tube (m)
- tube_lengthfloat or series of length m
Length of the tube (m)
- flowratefloat or series of length n
Volumetric flow rate (lpm)
- Returns:
- float or dataframe of shape (n,m)
Average residence time in seconds
- aerosol.functions.coagulation_coef(dp1, dp2, temp=293.15, pres=101325.0)[source]¶
Calculate Brownian coagulation coefficient (Fuchs)
- Parameters:
- dp1float
first particle diameter, unit: m
- dp2float or series of lenght m
second particle diameter, unit: m
- tempfloat or series of lenght n
air temperature, unit: K
- presfloat or series of lenght n
air pressure, unit: Pa
- Returns:
- float or dataframe of shape (n,m)
Brownian coagulation coefficient (Fuchs),
If dataframe is returned the columns correspond to diameter pairs (dp1,dp2) and are labeled by elements in dp2.
unit m3 s-1
- aerosol.functions.conical_dma_mob2volts(Q, R1_max, R2, L, alpha, Z)[source]¶
Conical DMA voltage corresponding to mobility
- Parameters:
- Qfloat
sheath flow rate, unit lpm
- R1_maxfloat
inner electrode radius at outlet at distance L, unit m
- R2float
outer electrode radius, unit m
- Lfloat
effective electrode length, unit m
- alphafloat
tapering angle, unit degrees
- Zfloat
mobility, unit cm2 s-1 V-1
- Returns:
- float
DMA voltage, unit V
- aerosol.functions.cross_corr_gr(df, dmin, dmax, window_width_hours=3.0, row_threshold=0.0, col_threshold=0.0)[source]¶
Calculate GR using cross-correlation method
- Parameters:
- dfpandas dataframe
Aerosol number size distribution
- dminfloat
Lower size limit for GR
- dmaxfloat
Upper size limit for GR
- smoothing_windowfloat
Window length used in smoothing the data in hours
- row_thresholdfloat
Maximum fraction of NaNs present in the rows
- col_thersholdfloat
Maximum fraction of NaNs present in the columns
- Returns:
- dictionary
Results in a dictionary:
gr: growth rate in nm/h
lags: time shifts applied in seconds
corrs: correlation coefficicients for the lags
max_corr_lag: time shift with maximum correlation
max_corr: the value of the maximum correlation
channels: concentration in the two size channels normalized between 0 and 1
- aerosol.functions.cs2coags(cs, dp, m=-1.6)[source]¶
Estimate coagulation sink from condensation sink
- Parameters:
- cspandas.Series
The condensation sink time series: unit s-1
- dpfloat
Particle diameter for which CoagS is calculated, unit: nm
- mfloat
Exponent in the equation
- Returns:
- coagspandas.Series
Coagulation sink time series for size dp
References
Kulmala et al (2012), doi:10.1038/nprot.2012.091
- aerosol.functions.datenum2datetime(datenum)[source]¶
Convert from matlab datenum to python datetime
- Parameters:
- datenumfloat or int
A serial date number representing the whole and fractional number of days from 1-Jan-0000 to a specific date (MATLAB datenum)
- Returns:
- pandas.Timestamp
- aerosol.functions.datetime2datenum(dt)[source]¶
Convert from python datetime to matlab datenum
- Parameters:
- dtdatetime object
- Returns:
- float
A serial date number representing the whole and fractional number of days from 1-Jan-0000 to a specific date (MATLAB datenum)
- aerosol.functions.denan(df)[source]¶
Interpolate away any nans and drop nan tails
- Parameters:
- dfpandas datafarme
- Returns:
- pandas dataframe
- aerosol.functions.diam2mob(dp, temp=293.15, pres=101325.0, ne=1)[source]¶
Convert electrical mobility diameter to electrical mobility in air
- Parameters:
- dpfloat
particle diameter(s), unit : nm
- tempfloat
ambient temperature default 20 C unit: K
- presfloat
ambient pressure, default 1 atm unit: Pa
- neint
number and polarity of charges on the aerosol particle default 1
- Returns:
- float
particle electrical mobility, unit: cm2 s-1 V-1
- aerosol.functions.dma_mob2volts(Q, R1, R2, L, Z)[source]¶
Cylindrical DMA voltage corresponding to mobility
- Parameters:
- Qfloat
sheath flow rate, unit lpm
- R1float
inner electrode radius, unit m
- R2float
outer electrode radius, unit m
- Lfloat
effective electrode length, unit m
- Zfloat
mobility, unit cm2 s-1 V-1
- Returns:
- float
DMA voltage, unit V
- aerosol.functions.dma_volts2mob(Q, R1, R2, L, V)[source]¶
Theoretical selected mobility from cylindrical DMA
- Parameters:
- Qfloat
sheath flow rate, unit lpm
- R1float
inner electrode radius, unit m
- R2float
outer electrode radius, unit m
- Lfloat
effective electrode length, unit m
- Vfloat or series
applied voltage, unit V
- Returns:
- float or series
selected mobility, unit cm2 s-1 V-1
- aerosol.functions.dndlogdp2dn(df)[source]¶
Convert from normalized number concentrations to unnormalized number concentrations.
- Parameters:
- dfdataframe
Aerosol number-size distribution (dN/dlogDp)
- Returns:
- dataframe
Aerosol number size distribution (dN)
- aerosol.functions.eq_charge_frac(dp, N)[source]¶
Calculate equilibrium charge fraction using Wiedensohler (1988) approximation
- Parameters:
- dpfloat
Particle diameter (m)
- Nint
Amount of elementary charge in range [-2,2]
- Returns:
- float
Fraction of particles of diameter dp having N elementary charges
- aerosol.functions.flow_velocity_in_pipe(tube_diam, flowrate)[source]¶
Calculate fluid speed from the flow rate in circular tube
- Parameters:
- tube_diamfloat or series of lenght m
Diameter of circular tube (m)
- flowratefloat or series of lenght n
Volumetric flow rate (lpm)
- Returns:
- float or dataframe of shape (n,m)
Speed of fluid (m/s)
- aerosol.functions.ions2particles(neg_ions, pos_ions, temp=293.15, mob_ratio=1.0)[source]¶
Estimate particle number size distribution from ions using Li et al. (2022)
- Parameters:
- neg_ionspandas dataframe of shape (n,m)
negative ion number size distribution
- pos_ionspandas dataframe of shape (n,m)
positive ion number size distribution
- tempfloat or series of length n
ambient temperature in K
- mob_ratiofloat
mobility ratio to be used default 1.0
- Returns:
- pandas dataframe of shape (n,m)
estimated particle number size distribution
References
Li et al. (2022), https://doi.org/10.1080/02786826.2022.2060795
- aerosol.functions.mean_free_path(temp, pres)[source]¶
Calculate mean free path in air
- Parameters:
- tempfloat or series of length n
air temperature, unit: K
- presfloat or series of length n
air pressure, unit: Pa
- Returns:
- float or series of length n
mean free path in air, unit: m
- aerosol.functions.mob2diam(Zp, temp=293.15, pres=101325.0, ne=1, tol=0.001, maxiter=20)[source]¶
Convert electrical mobility to electrical mobility diameter in air
- Parameters:
- Zpfloat
particle electrical mobility or mobilities, unit: cm2 s-1 V-1
- tempfloat
ambient temperature, unit: K
- presfloat
ambient pressure, unit: Pa
- neinteger
number and polarity of elementary charges on the aerosol particle
- Returns:
- float
particle diameter, unit: m
- aerosol.functions.nanoranking(df, dmin, dmax, row_threshold=0, col_threshold=0)[source]¶
Simplified method of calculating the nanorank
- Parameters:
- dfpandas dataframe
number size distribution
- dminfloat
lower diameter limit
- dmaxfloat
upper diameter limit
- row_thresholdfloat
maximum fraction of nans when calculating the number concentartion
- col_thresholdfloat
maximum fraction of nans in the conentration timeseries
- Returns:
- dictionary
the result dictionary has the following keys:
norm_conc: Concentration with removed background
rank: Value of the peak normalized concentration
rank_time: Time where the peak occurs
Notes
The nanorank is calculated for one day day. See Aliaga et al 2023
- aerosol.functions.particle_diffusivity(dp, temp=293.15, pres=101325.0)[source]¶
Particle brownian diffusivity in air
- Parameters:
- dpfloat or series of lenght m
particle diameter, unit: m
- tempfloat or series of lenght n
air temperature, unit: K
- presfloat or series of lenght n
air pressure, unit: Pa
- Returns:
- float or dataframe of shape (n,m)
Particle Brownian diffusivity in air unit m2 s-1
- aerosol.functions.particle_mean_free_path(dp, temp=293.15, pres=101325.0)[source]¶
Particle mean free path in air
- Parameters:
- dpfloat or series of length m
particle diameter, unit: m
- tempfloat or series of length n
air temperature, unit: K
- presfloat or series of length n
air pressure, unit: Pa
- Returns:
- float or dataframe of shape (n,m)
Particle mean free path, unit: m
- aerosol.functions.particle_thermal_speed(dp, temp)[source]¶
Particle thermal speed
- Parameters:
- dpfloat or series
particle diameter, unit: m
- tempfloat or series
air temperature, unit: K
- Returns:
- float or dataframe
Particle thermal speed point, unit: m s-1
- aerosol.functions.pipe_reynolds(tube_diam, flowrate, temp=293.15, pres=101325.0)[source]¶
Calculate Reynolds number in a tube
- Parameters:
- tube_diamfloat or series of length m
Inner diameter of the tube (m)
- flowratefloat or series of lenght n
Volumetric flow rate (lpm)
- tempfloat or series of length n
Temperature in K
- presfloat or series of length n
Pressure in Pa
- Returns:
- float or dataframe of shape (n,m)
Reynolds number
- aerosol.functions.sample_from_dist(x, y, n)[source]¶
Draw n samples from empirical distribution defined by points (x,y)
- Parameters:
- xnumpy array
x-data points
- ynumpy array
y-data points
- nint
number of samples to draw
- Returns:
- numpy array
samples drawn
- aerosol.functions.slipcorr(dp, temp=293.15, pres=101325.0)[source]¶
Slip correction factor in air
- Parameters:
- dpfloat or series of lenght m
particle diameter, unit m
- tempfloat or series of length n
air temperature, unit K
- presfloat or series of lenght n
air pressure, unit Pa
- Returns:
- float or dataframe fo shape (n,m)
For dataframe the index is taken from temperature or pressure series. Columns are particle diameters. unit dimensionless
Notes
Correction is done according to Mäkelä et al. (1996)
- aerosol.functions.surf_dist(df)[source]¶
Calculate the aerosol surface area size distribution
- Parameters:
- dfpandas.DataFrame
Aerosol number-size distribution
- Returns:
- pandas.DataFrame
Aerosol surface area-size distribution unit: m2 cm-3
- aerosol.functions.thab_dp2volts(thab_voltage, dp)[source]¶
Convert particle diameters to DMA voltages
- Parameters:
- thab_voltagefloat
Voltage at THA+ peak (V)
- dpfloat or series
Particle diameters (nm)
- Returns:
- float or series:
DMA voltage (V) corresponding to dp
Notes
See https://doi.org/10.1016/j.jaerosci.2005.02.009
Assumptions:
Sheath flow is air
Mobility standard used is THA+ monomer
T = 293.15 K and p = 101325 Pa
- aerosol.functions.thab_volts2dp(thab_voltage, dma_voltage)[source]¶
Convert DMA voltages to particle diameters
- Parameters:
- thab_voltagefloat
Voltage at THA+ peak (V)
- dma_voltagefloat
DMA voltage (V)
- Returns:
- float:
particle diameter corresponding to DMA voltage (nm)
Notes
See https://doi.org/10.1016/j.jaerosci.2005.02.009
Assumptions:
Sheath flow is air
Mobility standard used is THA+ monomer
T = 293.15 K and p = 101325 Pa
- aerosol.functions.tubeloss(diam, flowrate, tubelength, temp=293.15, pres=101325.0)[source]¶
Calculate diffusional particle losses to walls of straight cylindrical tube assuming a laminar flow regime
- Parameters:
- diamfloat or series of length m
Particle diameters for which to calculate the losses, unit: m
- flowratefloat or series of length n
unit: L/min
- tubelengthfloat
Length of the cylindrical tube unit: m
- tempfloat or series of length n
temperature unit: K
- presfloat or series of lenght n
air pressure unit: Pa
- Returns:
- float or dataframe of shape (n,m)
Fraction of particles passing through. Each column represents diameter and each each row represents different temperature pressure and flowrate value
- aerosol.functions.tubeloss_turbulent(diam, flowrate, tube_length, tube_diam, temp=293.15, pres=101325.0)[source]¶
Calculate particle losses to walls of a straight cylindrical tube assuming a turbulent flow regime and air as the carrier gas.
- Parameters:
- diamfloat or series of length m
Particle diameters for which to calculate the losses, unit: m
- flowratefloat or series of length n
unit: L/min
- tube_lengthfloat
Length of the cylindrical tube unit: m
- tube_diamfloat
Diameter of the cylindrical tube unit: m
- tempfloat or series of length n
temperature unit: K
- presfloat or series of lenght n
air pressure unit: Pa
- Returns:
- float or dataframe of shape (n,m)
Fraction of particles passing through. Each column represents diameter and each each row represents different temperature pressure and flowrate value
aerosol.plotting module¶
- aerosol.plotting.generate_log_ticks(min_exp, max_exp)[source]¶
Generate ticks and ticklabels for log axis
- Parameters:
- min_expint
The exponent in the smallest power of ten
- max_expint
The exponent in the largest power of ten
- Returns:
- numpy.array
minor tick values
- numpy.array
major tick values
- list of strings
major tick labels (powers of ten)
- aerosol.plotting.generate_timeticks(t_min, t_max, minortick_interval, majortick_interval, ticklabel_format)[source]¶
- Parameters:
- t_minpandas timestamp
- t_maxpandas timestamp
- majortick_intervalpandas date frequency string
See for all options here: https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases
- minortick_intervalpandas date frequency string
- ticklabel_formatpython date format string
See for all options here: https://docs.python.org/3/library/datetime.html#strftime-and-strptime-format-code
- Returns:
- pandas DatetimeIndex
minor tick values
- pandas DatetimeIndex
major tick values
- pandas Index containing strings
major tick labels
- aerosol.plotting.plot_aerosol_dist(v, ax, norm='log', clim=None, cmap='turbo', xmajortick_interval=None, xminortick_interval=None, xticklabel_format='%Y-%m-%d %H:%M')[source]¶
Plot aerosol particle number-size distribution surface plot
- Parameters:
- vpandas.DataFrame or list of pandas.DataFrames
Aerosol number size distribution (continuous index)
- axaxes object
axis on which to plot the data
- normstring
Define how to normalize the colors. “linear” or “log”
- climlist with two elements or None
Minimum and maximum value in colorbar, None calculates automatic limits
- xminortick_intervalpandas date frequency string
See for all options here: https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases
- xmajortick_intervalpandas date frequency string
- xticklabel_formatstr
See for all options here: https://docs.python.org/3/library/datetime.html#strftime-and-strptime-format-code
- Returns:
- image handle
- colorbar handle
- aerosol.plotting.plot_one_to_one_line(ax=None, color='black', linestyle='--', linewidth=1)[source]¶
Draws a one-to-one line (y = x) on the visible part of the plot.
- Parameters:
- axmatplotlib.axes.Axes (optional)
The axes to draw the line on. If None, uses current axes.
- colorstr
Color of the line.
- linestylestr
Style of the line (e.g., ‘–’ for dashed)
- linewidthfloat
Width of the line
- aerosol.plotting.rotate_xticks(ax, degrees)[source]¶
- Parameters:
- axmatplotlib axes
- degreesint or float
number of degrees to rotate the xticklabels
- aerosol.plotting.set_legend_outside(ax, handles=None, labels=None, coords=(1, 1), **kwargs)[source]¶
Put legend outside axes (upper right corner)
- Parameters:
- axAxes
axes to add the legend to
- handleslist of handles
list of lines or points in the legend
- labelslist of strings
labels for the legend entries
- coordstuple of floats or ints
anchor point for the legend
- kwargs
other parameters passed to legend
- Returns:
- Legend
- aerosol.plotting.show_matrix_values(ax, matrix, text_format='%d', text_color='white')[source]¶
Plot numerical values on top of the cells when visualizing a matrix with imshow()
- Parameters:
- axmatplotlib.axes
- matrixnumpy 2d-array
- text_formatstr
- text_colorstr
- aerosol.plotting.subplot_aerosol_dist(vlist, grid, norm='log', vmin=10, vmax=10000, xminortick_interval='1H', xmajortick_interval='2H', xticklabel_format='%H:%M', keep_inner_ticklabels=False, hspace_padding=None, vspace_padding=None, subplot_labels=None, label_color='black', label_size=10, column_titles=None, fill_order='row', **kwargs)[source]¶
Plot aerosol size distributions (subplots)
- Parameters:
- vlistlist of pandas.DataFrames
Aerosol size distributions (continuous index)
- gridtuple (rows,columns)
define number of rows and columns
- normstring
Define how to normalize the colors. “linear” or “log”
- vminfloat or int
Minimum value in colorbar
- vmaxfloat or int
Maximum value in colorbar
- xminortick_intervalstr
A pandas date frequency string. See for all options here: https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases
- xmajortick_intervalstr
A pandas date frequency string
- xticklabel_formatstr
Date format string. See for all options here: https://docs.python.org/3/library/datetime.html#strftime-and-strptime-format-code
- keep_inner_ticklabelsbool
If True, use ticklabels in all subplots. If False, use ticklabels only on outer subplots.
- subplot_paddingnumber or None
Adjust space between subplots
- subplot_labelslist of str or None
The labels to put to labels the subplots with
- label_colorstr
- label_sizefloat
- column_titleslist of strings or None
- fill_orderstr
“rows” fills the subplots row by row “columns” fills the subplots column by column
- **kwargsoptional parameters passed to matplotlib imshow()
- Returns:
- figure object
- array of axes objects
- list of image handles
- colorbar handle
aerosol.fitting module¶
- aerosol.fitting.fit_multimode(x, y, timestamp, n_modes=None, n_samples=10000)[source]¶
Fit multimodal Gaussian to aerosol number-size distribution
- Parameters:
- x1d numpy array
log10 of bin diameters in nm.
- y1d numpy array
Number size distribution
- timestamppandas Timestamp
timestamp associated with the number size distributions
- n_modesint or None
number of modes to fit, if None the number is determined using automatic method
- n_samplesint
Number of samples to draw from the distribution during the fitting process.
- Returns:
- dict:
Fit results
- aerosol.fitting.fit_multimodes(df, n_modes=None, n_samples=10000)[source]¶
Fit multimodal Gaussian to a aerosol number size distribution (dataframe)
- Parameters:
- dfpandas DataFrame
Aerosol number size distribution
- n_modesint or None
number of modes to fit, if None the number is determined using automatic method for each timestamp
- n_samplesint
Number of samples to draw from the distribution during the fitting process.
- Returns:
- list:
List of fit results