Sun

class sbpy.calib.Sun(source, description=None, bibcode=None)[source]

Bases: SpectralStandard

Solar spectrum.

Most functionality requires synphot. The exception is retrieval of flux densities by filter name via the solar_fluxd context manager.

Parameters:
waveQuantity

Spectral wavelengths.

fluxdQuantity

Solar spectral flux density at 1 au.

descriptionstring, optional

Brief description of the source spectrum.

bibcodestring, optional

Bibliography code for sbpy.bib.register.

metadict, optional

Any additional meta data, passed on to SourceSpectrum.

Examples

Get the default solar spectrum:

>>> sun = Sun.from_default()

Create solar standard from synphot.SourceSpectrum:

>>> import astropy.constants as const
>>> import synphot
>>> source = (synphot.SourceSpectrum(synphot.BlackBody1D, temperature=5770)
...           * (3.14159 * const.R_sun**2 / const.au**2).decompose())
>>> sun = Sun(source, description='5770 K blackbody Sun')

Create solar standard from a file:

>>> sun = Sun.from_file('filename')        

Create solar standard from an array:

>>> import numpy as np
>>> import astropy.units as u
>>> import astropy.constants as const
>>> from astropy.modeling.models import BlackBody
>>> from sbpy.calib import Sun
>>> B = BlackBody(
...     temperature=5770 * u.K,
...     scale=((const.R_sun / const.au)**2
...            * u.W / u.m**2 / u.um / u.sr)
... )
>>> wave = np.logspace(-1, 2, 300) * u.um
>>> fluxd = B(wave) * np.pi * u.sr
>>> sun = Sun.from_array(wave, fluxd, description='5770 K blackbody Sun')

Interpolate to / evaluate at 0.62 μm:

>>> print(sun(0.62 * u.um))                
1611.7080046473307 W / (um m2)

Observe as through a spectrometer:

>>> import numpy as np
>>> import astropy.units as u
>>> sun = Sun.from_default()
>>> wave = np.linspace(1, 2.5) * u.um
>>> fluxd = sun.observe(wave)              

Observe as through a filter, integrating with the bandpass transmission:

>>> from sbpy.photometry import bandpass
>>> sun = Sun.from_default()
>>> v = bandpass('Johnson V')
>>> print(sun.observe(v))                  
1839.93273227  W / (um m2)

Observe through a filter, using sbpy’s filter calibration system:

>>> from sbpy.calib import solar_fluxd
>>> import sbpy.units as sbu
>>> solar_fluxd.set({
...     'V': -26.76 * sbu.VEGAmag,
...     'V(lambda eff)': 548 * u.nm,
...     'V(lambda pivot)': 551 * u.nm
... })                                     
>>> sun = Sun.from_default()
>>> print(sun.observe('V'))
-26.76 mag(VEGA)
Attributes:
wave - Wavelengths of the source spectrum.
fluxd - Source spectrum.
description - Brief description of the source spectrum.
meta - Meta data.

Attributes Summary

description

Description of the source spectrum.

fluxd

The source spectrum.

meta

source

wave

Wavelengths of the source spectrum.

Methods Summary

__call__(wave_or_freq[, unit])

Evaluate/interpolate the source spectrum.

color_index(wfb, unit)

Color index (magnitudes) and effective wavelengths.

from_array(wave, fluxd[, meta])

Create standard from arrays.

from_builtin(name)

Spectrum from a built-in sbpy source.

from_default()

Initialize new spectral standard from current default.

from_file(filename[, wave_unit, flux_unit, ...])

Load the source spectrum from a file.

observe(wfb[, unit, interpolate])

Observe as through filters or spectrometer.

observe_bandpass(bp[, unit])

Observe through a bandpass.

observe_filter_name(filt[, unit])

Flux density through this filter.

observe_spectrum(wave_or_freq[, unit])

Observe source as through a spectrometer.

redden(S)

Redden the spectrum.

show_builtin([print])

List built-in spectra.

Attributes Documentation

description

Description of the source spectrum.

fluxd

The source spectrum.

meta
source
wave

Wavelengths of the source spectrum.

Methods Documentation

__call__(wave_or_freq, unit=None)

Evaluate/interpolate the source spectrum.

Parameters:
wave_or_freqQuantity

Requested wavelengths or frequencies of the resulting spectrum.

unitstring, Unit, optional

Spectral units of the output (flux density). If None, the default depends on wave_or_freq: W/(m2 μm) for wavelengths, Jy for frequencies.

Returns:
fluxdQuantity

The spectrum evaluated/interpolated to the requested wavelengths or frequencies.

color_index(wfb, unit)

Color index (magnitudes) and effective wavelengths.

Parameters:
wfbQuantity, or tuple/list of SectralElement, string

Two wavelengths, frequencies, or bandpasses.

unitstring or MagUnit

Units for the calculation, e.g., astropy.units.ABmag or sbpy.units.VEGAmag.

Returns:
eff_waveQuantity

Effective wavelengths for each wfb.

ciQuantity

Color index, m_0 - m_1, where 0 and 1 are element indexes for wfb.

classmethod from_array(wave, fluxd, meta=None, **kwargs)

Create standard from arrays.

Parameters:
waveQuantity

The spectral wavelengths.

fluxdQuantity

The spectral flux densities.

metadict, optional

Meta data.

**kwargs

Passed to object initialization.

classmethod from_builtin(name)

Spectrum from a built-in sbpy source.

Parameters:
namestring

The name of the spectrum. See show_builtin() for available sources.

classmethod from_default()

Initialize new spectral standard from current default.

The spectrum will be None if synphot is not available.

classmethod from_file(filename, wave_unit=None, flux_unit=None, cache=True, **kwargs)

Load the source spectrum from a file.

NaNs are dropped.

Parameters:
filenamestring

The name of the file. See from_file for details.

wave_unit, flux_unitstr or Unit, optional

Wavelength and flux units in the file.

cachebool, optional

If True, cache the contents of URLs.

**kwargs

Passed to object initialization.

observe(wfb, unit=None, interpolate=False, **kwargs)

Observe as through filters or spectrometer.

Calls observe_bandpass, observe_spectrum, or self(), as appropriate.

Parameters:
wfbQuantity, SpectralElement, string

Wavelengths, frequencies, or bandpasses. Bandpasses may be a filter name (string). May also be a list of SpectralElement or strings.

unitstring, Unit, optional

Units of the output (spectral flux density).

interpolatebool, optional

For wavelengths/frequencies, set to True for interpolation instead of rebinning. Use this when the spectral resolution of the source is close to that of the requested wavelengths.

**kwargs

Additional keyword arguments for Observation, e.g., force.

Returns:
fluxdQuantity
observe_bandpass(bp, unit=None, **kwargs)

Observe through a bandpass.

Parameters:
bpSpectralElement, list, or tuple

Bandpass.

unitstring, Unit, optional

Spectral flux density units for the output. The default is W/(m2 μm).

**kwargs

Additional keyword arguments for Observation, e.g., force.

Returns:
lambda_effQuantity

Effective wavelength(s) of the observation(s).

fluxdQuantity

The spectrum rebinned.

observe_filter_name(filt, unit=None)

Flux density through this filter.

Does not use the spectrum, but instead the flux density calibration manager. If the name of the filter is BP, then the expected keys are:

BP : flux density BP(lambda eff) : effective wavelength, optional BP(lambda pivot) : pivot wavelength, optional

Parameters:
filtstring

Name of the filter.

unitstring, Unit, optional

Spectral flux density units for the output.

Returns:
lambda_eff: Quantity

Effective wavelength. None if it is not provided.

lambda_pivot: Quantity

Pivot wavelength. None if it is not provided.

fluxdQuantity

Spectral flux density.

Raises:
FilterLookupError if the filter is not defined.
observe_spectrum(wave_or_freq, unit=None, **kwargs)

Observe source as through a spectrometer.

Important

This method works best when the requested spectral resolution is lower than the spectral resolution of the internal data. If the requested wavelengths/frequencies are exactly the same as the internal spectrum, then the internal spectrum will be returned without binning. This special case does not work for subsets of the wavelengths.

Parameters:
wave_or_freqQuantity

Wavelengths or frequencies of the spectrum. Spectral bins will be centered at these values. The length must be larger than 1.

unitstring, Unit, optional

Spectral flux density units for the output. If None, the default is W/(m2 μm) for wavelengths, Jy for frequencies.

**kwargs

Additional keyword arguments for Observation, e.g., force.

Returns:
fluxdQuantity

The spectrum rebinned.

Raises:
SinglePointSpectrumError - If requested wavelengths or

frequencies has only one value.

Notes

Method for spectra adapted from AstroBetter post by Jessica Lu: https://www.astrobetter.com/blog/2013/08/12/python-tip-re-sampling-spectra-with-pysynphot/

redden(S)

Redden the spectrum.

Parameters:
SSpectralGradient

The spectral gradient to redden.

Returns:
specSpectralSource

Reddened spectrum

classmethod show_builtin(print=True)

List built-in spectra.