The pyspectral API

Blackbody radiation

Planck radiation equation.

pyspectral.blackbody.blackbody(wavel, temp)

Derive the Planck radiation as a function of wavelength.

SI units. blackbody(wavelength, temperature) wavel = Wavelength or a sequence of wavelengths (m) temp = Temperature (scalar) or a sequence of temperatures (K)

Output: The spectral radiance per meter (not micron!)
Unit = W/m^2 sr^-1 m^-1
pyspectral.blackbody.blackbody_rad2temp(wavelength, radiance)

Derive brightness temperatures from radiance using the Planck function.

Wavelength space. Assumes SI units as input and returns temperature in Kelvin

pyspectral.blackbody.blackbody_wn(wavenumber, temp)

Derive the Planck radiation as a function of wavenumber.

SI units. blackbody_wn(wavnum, temperature) wavenumber = A wavenumber (scalar) or a sequence of wave numbers (m-1) temp = A temperatfure (scalar) or a sequence of temperatures (K)

Output: The spectral radiance in Watts per square meter per steradian

per m-1: Unit = W/m^2 sr^-1 (m^-1)^-1 = W/m sr^-1

Converting from SI units to mW/m^2 sr^-1 (cm^-1)^-1: 1.0 W/m^2 sr^-1 (m^-1)^-1 = 0.1 mW/m^2 sr^-1 (cm^-1)^-1

pyspectral.blackbody.blackbody_wn_rad2temp(wavenumber, radiance)

Derive brightness temperatures from radiance using the Planck function.

Wavenumber space

pyspectral.blackbody.planck(wave, temperature, wavelength=True)

Derive the Planck radiation as a function of wavelength or wavenumber.

SI units. _planck(wave, temperature, wavelength=True) wave = Wavelength/wavenumber or a sequence of wavelengths/wavenumbers (m or m^-1) temp = Temperature (scalar) or a sequence of temperatures (K)

Output: Wavelength space: The spectral radiance per meter (not micron!)

Unit = W/m^2 sr^-1 m^-1

Wavenumber space: The spectral radiance in Watts per square meter per steradian per m-1: Unit = W/m^2 sr^-1 (m^-1)^-1 = W/m sr^-1

Converting from SI units to mW/m^2 sr^-1 (cm^-1)^-1: 1.0 W/m^2 sr^-1 (m^-1)^-1 = 1.0e5 mW/m^2 sr^-1 (cm^-1)^-1

Spectral responses

Reading the spectral responses in the internal pyspectral hdf5 format.

class pyspectral.rsr_reader.RSRDataBaseClass

Bases: object

Data container for the Relative Spectral Responses for all (supported) satellite sensors.

rsr_data_version_uptodate

Check whether RSR data are up to date.

class pyspectral.rsr_reader.RelativeSpectralResponse(platform_name=None, instrument=None, **kwargs)

Bases: pyspectral.rsr_reader.RSRDataBaseClass

Container for the relative spectral response functions for various satellite imagers.

convert()

Convert spectral response functions from wavelength to wavenumber.

get_number_of_detectors4bandname(h5f, bandname)

For a band name get the number of detectors, if any.

get_relative_spectral_responses(h5f)

Read the rsr data and add to the object.

integral(bandname)

Calculate the integral of the spectral response function for each detector.

load()

Read the internally formatet hdf5 relative spectral response data.

set_band_central_wavelength_per_detector(h5f, bandname, detector_name)

Set the central wavelength for the band and detector.

set_band_names(h5f)

Set the band names.

set_band_responses_per_detector(h5f, bandname, detector_name)

Set the RSR responses for the band and detector.

set_band_wavelengths_per_detector(h5f, bandname, detector_name)

Set the RSR wavelengths for the band and detector.

set_description(h5f)

Set the description.

set_instrument(h5f)

Set the instrument name.

set_platform_name(h5f)

Set the platform name.

pyspectral.rsr_reader.check_and_download(**kwargs)

Do a check for the version and attempt downloading only if needed.

pyspectral.rsr_reader.main()

Main.

Base class for reading raw instrument spectral responses.

class pyspectral.raw_reader.InstrumentRSR(bandname, platform_name, bandnames=None)

Bases: object

Base class for the raw (agency dependent) instrument response functions.

Solar irradiance

Read solar irradiances and calculate solar flux.

Module to read solar irradiance spectra and calculate the solar flux over various instrument bands given their relative spectral response functions

class pyspectral.solar.SolarIrradianceSpectrum(filename, **options)

Bases: object

Total Top of Atmosphere (TOA) Solar Irradiance Spectrum.

Wavelength is in units of microns (10^-6 m). The spectral Irradiance in the file TOTAL_IRRADIANCE_SPECTRUM_2000ASTM is in units of W/m^2/micron

convert2wavenumber()

Convert from wavelengths to wavenumber.

Units:
Wavelength: micro meters (1e-6 m) Wavenumber: cm-1
inband_solarflux(rsr, scale=1.0, **options)

Get the in band solar flux.

Derive the inband solar flux for a given instrument relative spectral response valid for an earth-sun distance of one AU.

inband_solarirradiance(rsr, scale=1.0, **options)

Get the in band solar irradiance.

Derive the inband solar irradiance for a given instrument relative spectral response valid for an earth-sun distance of one AU.

(Same as the in band solar flux).

interpolate(**options)

Interpolate Irradiance to a specified evenly spaced resolution/grid.

This is necessary to make integration and folding (with a channel relative spectral response) straightforward.

dlambda = wavelength interval in microns start = Start of the wavelength interval (left/lower) end = End of the wavelength interval (right/upper end) options: dlambda: Delta wavelength used when interpolating/resampling ival_wavelength: Tuple. The start and end interval in wavelength space, defining where to integrate/convolute the spectral response curve on the spectral irradiance data.

plot(plotname=None, **options)

Plot the data.

solar_constant()

Calculate the solar constant.

Near-Infrared reflectance

Derive NIR reflectances of a a band in the 3-4 micron window region.

Derive the Near-Infrared reflectance of a given band in the solar and thermal range (usually the 3.7-3.9 micron band) using a thermal atmospheric window channel (usually around 11-12 microns).

class pyspectral.near_infrared_reflectance.Calculator(platform_name, instrument, band, detector='det-1', wavespace='wavelength', solar_flux=None, sunz_threshold=85.0, masking_limit=85.0)

Bases: pyspectral.radiance_tb_conversion.RadTbConverter

A thermal near-infrared (~3.7 micron) band reflectance calculator.

Given the relative spectral response of the NIR band, the solar zenith angle, and the brightness temperatures of the NIR and the Thermal bands, derive the solar reflectance for the NIR band removing the thermal (terrestrial) part. The in-band solar flux over the NIR band is optional. If not provided, it will be calculated here!

The relfectance calculated is without units and should be between 0 and 1.

derive_rad39_corr(bt11, bt13, method='rosenfeld')

Derive the CO2 correction to be applied to the 3.9 channel.

Derive the 3.9 radiance correction factor to account for the attenuation of the emitted 3.9 radiance by CO2 absorption. Requires the 11 micron window band and the 13.4 CO2 absorption band, as e.g. available on SEVIRI. Currently only supports the Rosenfeld method

emissive_part_3x(tb=True)

Get the emissive part of the 3.x band.

reflectance_from_tbs(sun_zenith, tb_near_ir, tb_thermal, **kwargs)

Derive reflectances from Tb’s in the 3.x band.

The relfectance calculated is without units and should be between 0 and 1.

Inputs:

sun_zenith: Sun zenith angle for every pixel - in degrees

tb_near_ir: The 3.7 (or 3.9 or equivalent) IR Tb’s at every pixel
(Kelvin)
tb_thermal: The 10.8 (or 11 or 12 or equivalent) IR Tb’s at every
pixel (Kelvin)
tb_ir_co2: The 13.4 micron channel (or similar - co2 absorption band)
brightness temperatures at every pixel. If None, no CO2 absorption correction will be applied.
pyspectral.near_infrared_reflectance.get_as_array(variable)

Return variable as a Dask or Numpy array.

Variable may be a scalar, a list or a Numpy/Dask array.

Rayleigh scattering

Atmospheric correction of shortwave imager bands in the wavelength range 400 to 800 nm.

exception pyspectral.rayleigh.BandFrequencyOutOfRange

Bases: ValueError

Exception when the band frequency is out of the visible range.

class pyspectral.rayleigh.Rayleigh(platform_name, sensor, **kwargs)

Bases: pyspectral.rayleigh.RayleighConfigBaseClass

Container for the atmospheric correction of satellite imager bands.

This class removes background contributions of Rayleigh scattering of molecules and Mie scattering and absorption by aerosols.

get_effective_wavelength(bandname)

Get the effective wavelength with Rayleigh scattering in mind.

get_reflectance(sun_zenith, sat_zenith, azidiff, bandname, redband=None)

Get the reflectance from the three sun-sat angles.

get_reflectance_lut()

Get reflectance LUT.

Read the LUT with reflectances as a function of wavelength, satellite zenith secant, azimuth difference angle, and sun zenith secant.

class pyspectral.rayleigh.RayleighConfigBaseClass(aerosol_type, atm_type='us-standard')

Bases: object

A base class for the Atmospheric correction, handling the configuration and LUT download.

lutfiles_version_uptodate

Tell whether LUT file is up to date or not.

pyspectral.rayleigh.check_and_download(**kwargs)

Download atm correction LUT tables if they are not up to date already.

Do a check for the version of the atmospheric correction LUTs and attempt downloading only if needed.

pyspectral.rayleigh.get_reflectance_lut(filename)

Get reflectance LUT.

Read the LUT with reflectances as a function of wavelength, satellite zenith secant, azimuth difference angle, and sun zenith secant

Utils

Utility functions.

class pyspectral.utils.NullHandler(level=0)

Bases: logging.Handler

Empty handler.

emit(record)

Record a message.

pyspectral.utils.bytes2string(var)

Decode a bytes variable and return a string.

pyspectral.utils.convert2hdf5(ClassIn, platform_name, bandnames, scale=1e-06)

Retrieve original RSR data and convert to internal hdf5 format.

scale is the number which has to be multiplied to the wavelength data in order to get it in the SI unit meter

pyspectral.utils.convert2str(value)

Convert a value to string. :param value: Either a str, bytes or 1-element numpy array

pyspectral.utils.convert2wavenumber(rsr)

Convert Spectral Responses from wavelength to wavenumber space.

Take rsr data set with all channels and detectors for an instrument each with a set of wavelengths and normalised responses and convert to wavenumbers and responses

Rsr:Relative Spectral Response function (all bands)
Returns:retv: Relative Spectral Responses in wave number space :info: Dictionary with scale (to go convert to SI units) and unit
pyspectral.utils.debug_on()

Turn debugging logging on.

pyspectral.utils.download_luts(**kwargs)

Download the luts from internet.

pyspectral.utils.download_rsr(**kwargs)

Download the relative spectral response functions.

Download the pre-compiled hdf5 formatet relative spectral response functions from the internet

pyspectral.utils.get_bandname_from_wavelength(sensor, wavelength, rsr, epsilon=0.1, multiple_bands=False)

Get the bandname from h5 rsr provided the approximate wavelength.

pyspectral.utils.get_central_wave(wav, resp, weight=1.0)

Calculate the central wavelength or the central wavenumber.

Calculate the central wavelength or the central wavenumber, depending on which parameters is input. On default the weighting funcion is f(lambda)=1.0, but it is possible to add a custom weight, e.g. f(lambda) = 1./lambda**4 for Rayleigh scattering calculations

pyspectral.utils.get_logger(name)

Return logger with null handle.

pyspectral.utils.get_wave_range(in_chan, threshold=0.15)

Return central, min and max wavelength in an RSR greater than threshold.

An RSR function will generally start near zero, increase to a maximum and then drop back to near zero. This function takes advantage of this to find the first and last points where the RSR is greater than a threshold. These points are then defined as the minimum and maximum wavelengths for a given channel, and can be used, for example, in Satpy reader YAML files.

pyspectral.utils.logging_off()

Turn logging off.

pyspectral.utils.logging_on(level=30)

Turn logging on.

pyspectral.utils.np2str(value)

Convert an numpy.string_ to str. :param value: scalar or 1-element numpy array to convert :type value: ndarray

Raises:ValueError – if value is array larger than 1-element or it is not of type numpy.string_ or it is not a numpy array
pyspectral.utils.sort_data(x_vals, y_vals)

Sort the data so that x is monotonically increasing and contains no duplicates.