firefly instrument

This class handles instrumental properties that influence the FIREFLY fitting process.

Provides a set of functions to handle instrumental effects.

log_rebin() has been pulled from mangadap.contrib.ppxf_util.py and modified.

Source location:
$MANGADAP_DIR/python/mangadap/util/instrument.py

Python2/3 compliance:

from __future__ import division
from __future__ import print_function
from __future__ import absolute_import

import sys
if sys.version > '3':
    long = int

Imports:

import warnings
import numpy
from scipy.interpolate import InterpolatedUnivariateSpline
import astropy.constants
Revision history:
27 May 2015: Original implementation by K. Westfall (KBW) based on downgrader_MANGA.f provided by D. Thomas, O. Steele, D. Wilkinson, D. Goddard. 13 June 2015: D.Wilkinson edit to not calculate unimportant convolution terms -> runs 5x faster.
class firefly_instrument.convolution_integral_element(y, sigma, ye=None)[source]

Support class for variable sigma convolution. See convolution_variable_sigma().

Args:

y (np.ndarray): Vector to convolve sigma (np.ndarray): Coordinate-dependent standard deviation of the

Gaussian kernel

ye (np.ndarray): (Optional) Error in the vector to convolve

Raises:
Exception: Raised if y is not a 1D vector, or if the shape of
y and sigma (and ye if provided) are different.
Attributes:

x (np.ndarray): Pixel coordinate vector y (np.ndarray): Vector to convolve ye (np.ndarray): Error in the vector to convolve sigma (np.ndarray): Coordinate-dependent standard deviation of the

Gaussian kernel
norm (np.ndarray): Gaussian normalization; calculated once for
efficiency
error(xc)[source]

Calculates the error in the weighted mean of y using nominal error propagation. The weights are defined by a Gaussian with standard deviation sigma and centered at xc.

Args:
xc (float): Center for the Gaussian weighting function
Returns:
float: The error in the weighted mean of y
firefly_instrument.convolution_variable_sigma(y, sigma, ye=None)[source]

Convolve a discretely sampled function \(y(x)\) with a Gaussian kernel, \(g\), where the standard deviation of the kernel is a function of \(x\), \(\sigma(x)\). Nominal calculations can be performed to propagate the error in the result; these calculations do not include the covariance between the pixels, which will mean that the calculations likely have significant error!

The convolution is defined as:

\[\begin{split}(y\ast g)(x) &= \int_{-\infty}^{\infty} y(X)\ g(\sigma,x-X)\ dX \\ &= \int_{-\infty}^{\infty} \frac{y(X)}{\sqrt{2\pi}\ \sigma(X)}\ \exp\left(-\frac{(x-X)^2}{2\ \sigma(X)^2}\right) dX .\end{split}\]

To minimize edge effects and account for the censoring of the data (finite range in \(x\)), the convolution is actually calculated as a definite integral and normalized as follows:

\[(y\ast g)(x) \sim\frac{ \int_{x-n_\sigma*\sigma(x)}^{x+n_\sigma*\sigma(x)} y(X)\ g(\sigma,x-X)\ dX}{ \int_{x-n_\sigma*\sigma(x)}^{x+n_\sigma*\sigma(x)} g(\sigma,x-X)\ dX} .\]

The above is identical to getting the weighted mean of \(y\) at each position \(x\), where the weights are defined by a Gaussian kernel centered at \(x\) with a variable dispersion.

Use of this function requires:
  • y and sigma must be 1D vectors
  • y and sigma must be uniformly sampled on the same grid
  • sigma must be in pixel units.
Args:

y (np.ndarray): A uniformly sampled function to convolve. sigma (np.ndarray): The standard deviation of the Gaussian

kernel sampled at the same positions as y. The units of sigma must be in pixels.
ye (np.ndarray): (Optional) Errors in the function
\(y(x)\).
Returns:
np.ndarray: Arrays with the convolved function :math:`(yast
g)(x)` sampled at the same positions as the input \(x\) vector and its error. The second array will be returned as None if the error vector is not provided.
Raises:
Exception: Raised if trying to calculate the errors because they haven’t been implemented yet.
firefly_instrument.downgrade(wave, flux, deltal_in, sigma_galaxy, wave_instrument, r_instrument)[source]

Adapted from the manga DAP downgrader from Kyle Westfall.

Downgrades an input spectrum to a given galaxy velocity dispersion using the input SEDs resolution and the resolution of the observation.

Returns flux of downgraded SED.

firefly_instrument.log_rebin(lamRange, spec, oversample=None, velscale=None, flux=False, log10=False, newRange=None, wave_in_ang=False, unobs=0.0)[source]

Note

Copyright (C) 2001-2014, Michele Cappellari E-mail: cappellari_at_astro.ox.ac.uk

This software is provided as is without any warranty whatsoever. Permission to use, for non-commercial purposes is granted. Permission to modify for personal or internal use is granted, provided this copyright and disclaimer are included unchanged at the beginning of the file. All other rights are reserved.

Logarithmically rebin a spectrum, while rigorously conserving the flux. Basically the photons in the spectrum are simply ridistributed according to a new grid of pixels, with non-uniform size in the spectral direction.

This routine makes the `standard’ zero-order assumption that the spectrum is constant within each pixels. It is possible to perform log-rebinning by assuming the spectrum is represented by a piece-wise polynomial of higer degree, while still obtaining a uniquely defined linear problem, but this reduces to a deconvolution and amplifies noise.

Warning

This assumption can be poor for sharp features in the spectrum. Beware if resampling spectra with strong, marginally sampled features!

This same routine can be used to compute approximate errors of the log-rebinned spectrum. To do this type the command

>>> err2New, logLam, velscale = log_rebin(lamRange, np.square(err))

and the desired errors will be given by np.sqrt(err2New).

Warning

This rebinning of the error-spectrum is very approximate as it does not consider the correlation introduced by the rebinning!

Args:

lamRange (np.ndarray): two elements vector containing the
central wavelength of the first and last pixels in the spectrum, which is assumed to have constant wavelength scale! E.g. from the values in the standard FITS keywords: LAMRANGE = CRVAL1 + [0,CDELT1*(NAXIS1-1)]. It must be LAMRANGE[0] < LAMRANGE[1].

spec (np.ndarray): Input spectrum. oversample (int): (Optional) Oversampling can be done, not to

loose spectral resolution, especally for extended wavelength ranges and to avoid aliasing. Default is to provide the same number of output pixels as input.
velscale (float): Velocity scale in km/s per pixels. If this
variable is not defined, then it will contain in output the velocity scale. If this variable is defined by the user it will be used to set the output number of pixels and wavelength scale.
flux (bool): (Optional) Set this keyword to preserve total flux.

In this case the log rebinning changes the pixels flux in proportion to their dLam so the following command will show large differences beween the spectral shape before and after log_rebin():

# Plot log-rebinned spectrum
pyplot.plot(exp(logLam), specNew)
pyplot.plot(np.arange(lamRange[0],lamRange[1],spec.size), spec, 'g')
pyplot.show()

By default, when this keyword is not set, the above two lines produce two spectra that almost perfectly overlap each other.

log10 (bool): (Optional) Flag that the spectrum should be binned
in units of base-10 log wavelength, instead of natural log
newRange (np.ndarray): (Optional) Force the spectrum to be
sampled to a new spectral range (lamRange is the existing spectral range).
wave_in_ang (bool): (Optional) Return the wavelength coordinates
in angstroms, not log(angstroms)
unobs (float): (Optional) Default value for unobserved spectral
regions.
Returns:
np.ndarray, float
: Returns three variables: logarithmically
rebinned spectrum, the log of the wavelength at the geometric center of each pixel, and the velocity scale of each pixel in km/s.
Raises:
ValueError: Raised if the input spectrum is not a
one-dimensional np.ndarray.
Modification History:
V1.0.0: Using interpolation. Michele Cappellari, Leiden, 22 October 2001
V2.0.0: Analytic flux conservation. MC, Potsdam, 15 June 2003
V2.1.0: Allow a velocity scale to be specified by the user. MC, Leiden, 2 August 2003
V2.2.0: Output the optional logarithmically spaced wavelength at the geometric mean of the wavelength at the border of each pixel. Thanks to Jesus Falcon-Barroso. MC, Leiden, 5 November 2003
V2.2.1: Verify that lamRange[0] < lamRange[1]. MC, Vicenza, 29 December 2004
V2.2.2: Modified the documentation after feedback from James Price. MC, Oxford, 21 October 2010
V2.3.0: By default now preserve the shape of the spectrum, not the total flux. This seems what most users expect from the procedure. Set the keyword /FLUX to preserve flux like in previous version. MC, Oxford, 30 November 2011
V3.0.0: Translated from IDL into Python. MC, Santiago, 23 November 2013
V3.1.0: Fully vectorized log_rebin. Typical speed up by two orders of magnitude. MC, Oxford, 4 March 2014
05 Jun 2015: (K. Westfall, KBW) Pulled from ppxf_util.py. Conform to mangadap documentation standard. Transcribe edits made to IDL version that provides for the log10 and newRange arguments. Add option to return wavelength in angstroms, not log(angstroms). Break out determination of input and output spectrum pixel coordinates to a new function, log_rebin_pix(). Added default value for unobserved pixels. Default behavior unchanged.
firefly_instrument.log_rebin_pix(lamRange, n, oversample=None, velscale=None, log10=False, newRange=None)[source]

Determine the number of new pixels and their coordinate step when rebinning a spectrum in geometrically stepped bins. The input spectrum must be sampled linearly in wavelength. This is primarily a support routine for log_rebin().

Although breaking this out from the main log_rebin() function leads to a few repeat calculations in that function, the use of this function is in determine a common wavelength range for a large number of spectra before resampling the spectra themselves. See mangadap.TplLibrary .

Args:
lamRange (np.ndarray): two elements vector containing the
central wavelength of the first and last pixels in the spectrum, which is assumed to have constant wavelength scale! E.g. from the values in the standard FITS keywords: LAMRANGE = CRVAL1 + [0,CDELT1*(NAXIS1-1)]. It must be LAMRANGE[0] < LAMRANGE[1].

n (int): Number of pixels in the original spectrum. oversample (int): (Optional) Oversampling can be done, not to

loose spectral resolution, especally for extended wavelength ranges and to avoid aliasing. Default is to provide the same number of output pixels as input.
velscale (float): Velocity scale in km/s per pixels. If this
variable is not defined, then it will contain in output the velocity scale. If this variable is defined by the user it will be used to set the output number of pixels and wavelength scale.
log10 (bool): (Optional) Flag that the spectrum should be binned
in units of base-10 log wavelength, instead of natural log
newRange (np.ndarray): (Optional) Force the spectrum to be
sampled to a new spectral range (lamRange is the existing spectral range).
Returns:
float, int
: Returns (1) the linear wavelength step of each
pixel in the input spectrum, (2) the number of pixels for the rebinned spectrum, (3) the log-linear wavelength step for each pixel in the new spectrum, (4) the velocity step for each pixel in the new spectrum.
Raises:
ValueError: Raised if the input wavelength range (lamRange or
newRange) does not have two elements or is not sorted.
firefly_instrument.match_spectral_resolution(wave, flux, sres, new_sres_wave, new_sres, ivar=None, mask=None, min_sig_pix=0.0, no_offset=True, variable_offset=False, log10=False, new_log10=False)[source]

Adjust the existing spectral resolution of to a lower resolution as best as possible. The primary functionality is in spectral_resolution, which determines the Gaussian kernel parameters needed to match the resolution, and convolve_variable_sigma(), which actually performs the convolution to match the resolution.

In particular, see spectral_resolution.GaussianKernelDifference() for a description of how the kernel parameters are determined and how regions where the target resolution is higher than the existing resolution. In this case, one of the options is to adopt an offset of the resolution (in km/s) that could be corrected for in subsequent analysis. In this case, setting variable_offset to True allows the offset to be different for all input spectra. If one expects to combine the spectra, the default behavior should be used, forcing all the spectra to have a constant offset.

Args:
wave (np.ndarray): A 1D or 2D (:math:`N_{rm spec}times
N_{rm pix}`) array with the wavelength in angstroms for a set of spectra. The sampling may be either in linear steps of wavelength or \(\log_{10}\) steps, as set using log10.
flux (np.ndarray): A 1D or 2D (:math:`N_{rm spec}times
N_{rm pix}`) array with the flux sampled at the provided wavelengths.
sres (np.ndarray): A 1D or 2D (:math:`N_{rm spec}times
N_{rm pix}`) array with the spectral resolution, \(R\), at the provided wavelengths.
new_sres_wave (np.ndarray): A 1D vector with the wavelength
in angstroms at which the new resolution of the input spectra has been sampled. The sampling may be either in linear steps of wavelength or \(\log_{10}\) steps, as set using new_log10.
new_sres (np.ndarray): A 1D vector with the new resolution
for the input spectra.
ivar (np.ndarray): (Optional) A 1D or 2D (:math:`N_{rm
spec}times N_{rm pix}`) array with the inverse variance of the flux sampled at the provided wavelengths. Currently never used and, if provided, the output inverse variances are simply a carbon copy of this array.
mask (np.ndarray): (Optional) A 1D or 2D (:math:`N_{rm
spec}times N_{rm pix}`) array with a uint mask for the flux sampled at the provided wavelengths.
no_offset (bool): (Optional) Force \(\sigma^2_{v,o} = 0\) by
masking regions with \(\sigma_{p,d} < -\epsilon_\sigma\); i.e., the value of this arguments selects Option 1 (True) or Option 2 (False). See spectral_resolution.GaussianKernelDifference().
min_sig_pix (float): (Optional) Minimum value of the standard
deviation in pixels allowed before assuming the kernel is a Delta function.
variable_offset (bool): (Optional) Flag to allow the offset
applied to each spectrum (when the input contains more than one spectraum) to be tailored to each spectrum. Otherwise (variable_offset=False) the offset is forced to be the same for all spectra.
log10 (bool): (Optional) Flag that the spectrum has been binned
logarithmically (base 10) in wavelength
new_log10 (bool): (Optional) Flag that the coordinates of the
new spectral resolution are spectrum as been binned logarithmically (base 10) in wavelength.

Returns:

np.ndarray : Four or Five arrays are returned:

  • A 1D or 2D (\(N_{\rm spec}\times N_{\rm pix}\)) array with the resolution-matched flux sampled at the input wavelengths.
  • A 1D or 2D (\(N_{\rm spec}\times N_{\rm pix}\)) array with the spectral resolution, \(R\), of the resolution-matched spectra at the provided wavelengths.
  • A 1D vector with any constant offset in resolution in km/s between the targetted value and the result. See spectral_resolution.GaussianKernelDifference().
  • A 1D or 2D (\(N_{\rm spec}\times N_{\rm pix}\)) array with a uint mask for the resolution-matched flux sampled at the input wavelengths. This is returned regardless of whether an input mask was provided. Any pixel that had a resolution that was lower than the target resolution (up to some tolerance defined by min_sig_pix) is returned as masked.
  • (Optional) A 1D or 2D (\(N_{\rm spec}\times N_{\rm pix}\)) array with the inverse variance of the resolution-matched flux sampled at the input wavelengths. Only output if *ivar* is provided, and in that case this is just a carbon copy of the input array.
Raises:

Exception: Raised if:

  • the input wave array is 2D and the sres array is not; a 1D wavelength array is allowed for a 2D sres array but not vice versa
  • the number of spectral pixels in wave, flux, and sres is not the same
  • the shape of the flux, mask (if provided), and ivar (if provided) are not the same
  • the shape of the new_sres_wave and new_sres arrays are not the same and/or not 1D
class firefly_instrument.spectral_resolution(wave, sres, log10=False, interp_ext='extrapolate')[source]

Container class for the resolution, \(R = \lambda/\Delta\lambda\), of a spectrum. The primary functionality is to determine the parameters necessary to match the resolution of one spectrum to another. It can also be used as a function to interpolate the spectral resolution at a given wavelenth.

Args:
wave (np.ndarray): A 1D vector with the wavelength in
angstroms. The sampling may be either in linear steps of wavelength or \(\log_{10}\) steps.
sres (np.ndarray): A 1D vector with the spectral resolution,
\(R\), sampled at the positions of the provided wavelength vector.
log10 (bool): (Optional) Flag that the spectrum has been binned
logarithmically (base 10) in wavelength
interp_ext (int or str): (Optional) The value to pass as ext
to the interpolator, which defines its behavior when attempting to sample the spectral resolution beyond where it is defined. See scipy.interpolate.InterpolatedUnivariateSpline. Default is to extrapolate.
Raises:
Exception: Raised if wave is not a 1D vector or if wave and
sres do not have the same shape.
Attributes:
interpolator
(scipy.interpolate.InterpolatedUnivariateSpline): An object used to interpolate the spectral resolution at any given wavelength. The interpolation is hard-wired to be linear and its extrapolation behavior is defined by interp_ext. The wavelength and resolution vectors are held by this object for later reference if needed.
log10 (bool): Flag that the spectrum has been binned
logarithmically (base 10) in wavelength
cnst (mangadap.util.constants): Used to define the
conversion factor between a Gaussian sigma and FWHM

c (float): The speed of light; provided by astropy.constants. dv (float): The velocity step per pixel in km/s. Defined using

spectrum_velocity_scale() if log10 is True; otherwise set to None.
dw (float): The wavelength step per pixel in angstroms. Defined
as the wavelength step between the first two pixels if log10 is False; otherwise set to None.
min_sig (float): Minimum standard deviation allowed for the
kernel used to match two spectral resolutions. See GaussianKernelDifference().
sig_pd (np.ndarray): The standard deviation in pixels
required to match the spectral resolution of this object to the resolution defined by a different spectral_resolution object. See GaussianKernelDifference().
sig_mask (np.ndarray): A uint vector used to identify
measurements of sig_pd that should not be used to match the spectral resolution of this object to the resolution defined by a different spectral_resolution object. See GaussianKernelDifference().
sig_vo (float): A constant offset of the kernal standard
deviation in km/s that has been applied to sig_pd. See GaussianKernelDifference().

Warning

The default behavior of the interpolator is to extrapolate the input spectral resolution vector when trying to sample from regions beyond where it is sampled. Use interp_ext change this; see scipy.interpolate.InterpolatedUnivariateSpline.

GaussianKernelDifference(new_sres, no_offset=True, min_sig_pix=0.0)[source]

Determine the parameters of a Gaussian kernel required to convert the resolution of this object to the resolution of a different the spectral_resolution object, new_sres.

The spectral resolution is defined as \(R = \lambda/\Delta\lambda\), where \(\Delta\lambda\) is the FWHM of the spectral resolution element. The standard deviation of the resolution element in angstroms is then

\[\sigma_\lambda = \frac{\lambda}{f R}, \ \ {\rm where} \ \ f = \frac{{\rm FWHM_\lambda}}{\sigma_\lambda}.\]

Assuming a Gaussian (in angstroms) line-spread function:

\[\sigma^2_{\lambda,2} = \sigma^2_{\lambda,1} + \sigma^2_{\lambda,d}\]

such that

\[\sigma^2_{\lambda,d} = \left(\frac{\lambda}{f}\right)^2 (R^{-2}_2 - R^{-2}_1)\]

is the defining parameter of the Gaussian kernel needed to take a spectrum of resolution \(R_1\) to one with a resolution of \(R_2\).

For input to convolution_variable_sigma(), the wavelength-dependent parameter of the Gaussian kernel is converted to pixels. This function allows for the input spectra to be linearly sampled in angstroms or log10(angstroms). For the former (log10=False),

\[\sigma^2_{p,d} = \left(\frac{\lambda}{f\ \delta\lambda}\right)^2 (R^{-2}_2 - R^{-2}_1)\]

where \(\delta\lambda\) is the size of the pixel in angstroms. If the units are log10(angstrom) (log10=True), we approximate the velocity width of each pixel to be \(\delta v = c \ln(10.0) (\log\lambda[1]-\log\lambda[0])\), such that

\[\begin{split}\sigma^2_{p,d} &= \left(\frac{c}{ \delta v \lambda}\right)^2 \sigma^2_{\lambda,d} \\ &= \left(\frac{c}{ f\ \delta v}\right)^2 (R^{-2}_2 - R^{-2}_1)\ ;\end{split}\]

\(c\) is the speed of light in km/s.

The nominal use of this algorithm assumes \(R_1 \geq R_2\). However, in practice, convolution_variable_sigma() only uses a Gaussian kernel up to some minimum value of \(\epsilon_\sigma\); below this, the kernel is assumed to be a Delta function. Therefore, as long as

\[\sigma_{p,d} \equiv \sigma^2_{p,d}/\sqrt{|\sigma^2_{p,d}|} \geq -\epsilon_\sigma\ ,\]

the behavior of convolution_variable_sigma() should not be affected.

Even so, there may be spectral regions that do not have \(\sigma_{p,d} \geq -\epsilon_\sigma\); for such spectral regions there are three choices:

(Option 1) trim the spectral range to only those spectral regions where the existing resolution is better than the target resolution,

(Option 2) match the existing resolution to the target resolution up to some constant offset that must be accounted for in subsequent analyses, or

(Option 3) allow for a wavelength dependent difference in the spectral resolution that must be accounted for in subsequent analyses.

The choice of either Option 1 or 2 is selected by setting no_offset to, respectively, True or False; Option 1 is the default behavior. Currently, Option 3 is not allowed.

For Option 1, pixels with \(\sigma_{p,d} < -\epsilon_\sigma\) are masked (sigma_mask = 1); however, the returned values of \(\sigma_{p,d}\) are left unchanged.

For Option 2, we define

\[\sigma^2_{v,o} = -{\rm min}(\sigma^2_{v,d}) - {\rm max}(\epsilon_\sigma \delta v)^2\]

where \(\delta v\) is constant for the logarithmically binned spectrum and is wavelength dependent for the linearly binned spectra; in the latter case, the velocity step is determined for each pixel:

_wave = self.wave()
dv = c * (2.0*(_wave[1:] - _wave[0:-1]) / (_wave[1:] + _wave[0:-1]))

If \(\sigma^2_{v,o} > 0.0\), it must be that \({\rm min}(\sigma^2_{v,d}) < -{\rm max}(\epsilon_\sigma \delta v)^2\), such that an offset should be applied. In that case, the returned kernel parameters are

\[\sigma^\prime_{v,d} = \sqrt{\sigma^2_{v,d} + \sigma^2_{v,o}}\ .\]

with the units converted to pixels using the equations above, no pixels are masked, and \(\sqrt{\sigma^2_{v,o}}\) is returned for the offset. Otherwise, the offset is set to 0.

Args:
new_sres (spectral_resolution): Spectral resolution
to match to.
no_offset (bool): (Optional) Force :math:`sigma^2_{v,o} =
0` by masking regions with \(\sigma_{p,d} < -\epsilon_\sigma\); i.e., the value of this arguments selects Option 1 (True) or Option 2 (False).
min_sig_pix (float): (Optional) Minimum value of the
standard deviation allowed before assuming the kernel is a Delta function.
adjusted_resolution(indx=None)[source]

Return the resolution that should result from applying convolution_variable_sigma() to the spectrum associated with this spectral resolution object using sigma_pd. I.e., calculate:

\[R_2 = \left[ \left(\frac{f}{c}\right)^2 \sigma^2_{v,d} + R^{-2}_1\right]^{-1/2}\ . \]
Args:
indx (tuple): (Optional) Selection tuple used to return a
subset of the full resolution vector.
Returns:
np.ndarray: The (full or selected) vector with the
adjusted resolution.
match(new_sres, no_offset=True, min_sig_pix=0.0)[source]

Currently only an alias for GaussianKernelDifference().

offset_GaussianKernelDifference(offset)[source]

If the properties required to match the resolution of one spectrum to another has already been calculated (see GaussianKernelDifference()), this allows for one to apply an additional offset. The additional offset must be in km/s (not pixels).

The offset is applied in quadrature; however, the offset can be negative such that one can reduce sig_pd. Once converted to km/s, the offset is applied by calculating:

\[\sigma^{\prime\ 2}_{v,d} = \sigma^{2}_{v,d} + \sigma_{off}|\sigma_{off}|\ .\]

sig_vo is adjusted in the same way, and the change in \(\sigma^{\prime\ 2}_{v,d}\) is then propagated to sig_pd and sig_mask.

Args:
offset (float): Value of the standard deviation in km/s to
add in quadrature to a previously calculated sig_pd.
Raises:
Exception: Raised if the kernel properties have not yet been
defined.
sres()[source]

Return the resolution vector; held by interpolator.

wave()[source]

Return the wavelength vector; held by interpolator.

firefly_instrument.spectrum_velocity_scale(wave, log10=False)[source]

Determine the velocity sampling of an input wavelength coordinate vector. The wavelength vector is assumed to be logarithmically sampled, but its units are in angstroms.

Args:
wave (np.ndarray): Wavelength coordinates of each spectral
channel in angstroms. It is expected that the spectrum has been sampled geometrically.
log10 (bool): (Optional) Input spectrum has been sample
geometrically using the base 10 logarithm, instead of the natural logarithm.
Returns:
float : Velocity scale of the spectrum in km/s.
firefly_instrument.where_not(indx, size)[source]

Return a tuple with the indices of a vector that were not selected by a call to `np.where`_. The function currently only works for 1D vectors.

Args:
indx (tuple): Tuple returned by a call to `np.where`_ for a
1D vector.
size (int): Length of the original vector in the call to
`np.where`_.

Warning

Performs no checks of the input.