setigen API reference
Signal parameter functions
setigen.frame module
- class setigen.frame.Frame(waterfall=None, fchans=None, tchans=None, df=<Quantity 2.79396772 Hz>, dt=<Quantity 18.25361101 s>, fch1=<Quantity 6. GHz>, ascending=False, data=None, seed=None, **kwargs)[source]
Bases:
object
A class to facilitate the creation of entirely synthetic radio data (narrowband signals + background noise) as well as signal injection into existing observations.
- __init__(waterfall=None, fchans=None, tchans=None, df=<Quantity 2.79396772 Hz>, dt=<Quantity 18.25361101 s>, fch1=<Quantity 6. GHz>, ascending=False, data=None, seed=None, **kwargs)[source]
Initialize a Frame object either from an existing .fil/.h5 file or from frame resolution / size.
If you are initializing based on a .fil or .h5, pass in either the filename or the Waterfall object into the waterfall keyword.
Otherwise, you can initialize a frame by specifying the parameters
fchans
,tchans
,df
,dt
, and evenfch1
, if it’s important to specify frequencies (8 GHz is an arbitrary but reasonable choice otherwise). Note that the frame resolutionsdf
anddt
are given defaults based on the Breakthrough Listen high frequency resolution data product – be sure to change these if you are working with different kinds of data.The data keyword is only necessary if you are also preloading data that matches your specified frame dimensions and resolutions.
- Parameters:
waterfall (str or Waterfall, optional) – Name of filterbank file or Waterfall object for preloading data
fchans (int, optional) – Number of frequency samples
tchans (int, optional) – Number of time samples
df (astropy.Quantity, optional) – Frequency resolution (e.g. in u.Hz)
dt (astropy.Quantity, optional) – Time resolution (e.g. in u.s)
fch1 (astropy.Quantity, optional) – Central frequency of first channel. If
ascending=True
,fch1
is the minimum frequency; ifascending=False
(default),fch1
is the maximum frequency.ascending (bool, optional) – Specify whether frequencies should be in ascending order, so that
fch1
is the minimum frequency. Default is False, for whichfch1
is the maximum frequency. This is overwritten if a waterfall object is provided, whereascending
will be automatically determined by observational parameters.data (ndarray, optional) – 2D array of intensities to preload into frame
seed (None, int, Generator, optional) – Random seed or seed generator
**kwargs – For convenience, the
shape
keyword can be used in place of individually settingfchans
andtchans
, so thatshape=(tchans, fchans)
.
- classmethod from_data(df, dt, fch1, ascending, data, metadata={}, waterfall=None)[source]
Initialize Frame more directly from 2D numpy array of data.
- Parameters:
df (astropy.Quantity) – Frequency resolution (e.g. in u.Hz)
dt (astropy.Quantity) – Time resolution (e.g. in u.s)
fch1 (astropy.Quantity) – Central frequency of first channel. If
ascending=True
,fch1
is the minimum frequency; ifascending=False
(default),fch1
is the maximum frequency.ascending (bool) – Specify whether frequencies should be in ascending order, so that
fch1
is the minimum frequency. Default is False, for whichfch1
is the maximum frequency. This is overwritten if a waterfall object is provided, whereascending
will be automatically determined by observational parameters.data (ndarray) – 2D array of intensities to preload into frame
metadata (dict, optional) – Dictionary of features associated with the frame
waterfall (Waterfall, optional) – Associated Waterfall object if data is derived from another frame object (accessed via
frame.get_waterfall()
) or a blimpy waterfall object
- Returns:
frame – Frame object with preloaded data
- Return type:
- classmethod from_waterfall(waterfall)[source]
Instantiate Frame using a filterbank file or blimpy Waterfall object.
- classmethod from_backend_params(fchans=None, obs_length=300, sample_rate=3000000000.0, num_branches=1024, fftlength=1048576, int_factor=51, fch1=<Quantity 6. GHz>, ascending=False, data=None)[source]
Create frame based on backend / software related parameters. Either
fchans
ordata
must be provided to get number of frequency channels to create. If a 2D numpy array fordata
is provided,fchans
will be inferred. The parameterint_factor
must still be provided to determinetchans
; there is a check that the data dimensions also match. Since multipleint_factor
values may correspond to the sametchans
, for clarity we do not inferint_factor
just from the dimensions of the data.- Parameters:
fchans (int, optional) – Number of frequency samples. Should be provided if
data
is None.obs_length (float, optional) – Length of observation in seconds
sample_rate (float, optional) – Physical sample rate, in Hz, for collecting real voltage data
num_branches (int, optional) – Number of PFB branches. Note that this corresponds to
num_branches / 2
coarse channels.fftlength (int, optional) – FFT length to be used in fine channelization
int_factor (int, optional) – Integration factor used in fine channelization. Determines
tchans
.fch1 (astropy.Quantity, optional) – Central frequency of first channel. If
ascending=True
,fch1
is the minimum frequency; ifascending=False
(default),fch1
is the maximum frequency.ascending (bool, optional) – Specify whether frequencies should be in ascending order, so that
fch1
is the minimum frequency. Default is False, for whichfch1
is the maximum frequency.data (ndarray, optional) – 2D array of intensities to preload into frame. If provided,
fchans
will be inferred from this.
- Returns:
frame – Frame object with appropriate dimensions.
- Return type:
- property fmid
- property mjd
- property t_stop
- property obs_length
- property ts_ext
Extended time array of length
tchans + 1
, including the ending timestamp.
- property mean
- property std
- add_noise(x_mean, x_std=None, x_min=None, noise_type='chi2')[source]
By default, synthesize radiometer noise based on a chi-squared distribution. Alternately, can generate pure Gaussian noise.
Specifying
noise_type='chi2'
will only usex_mean
, and ignore other parameters. Specifyingnoise_type='gaussian'
will use all arguments (if provided).When adding Gaussian noise to the frame, the minimum is simply a lower bound for intensities in the data (e.g. it may make sense to cap intensities at 0), but this is optional.
- Parameters:
x_mean (float) – Target mean
x_std (float, optional) – Target standard deviation
x_min (float, optional) – Lower bound for Gaussian noise
noise_type ({"chi2", "gaussian", "normal"}, default: "chi2") – Distribution to use for synthetic noise
- Returns:
noise – Array of synthetic noise
- Return type:
ndarray
- add_noise_from_obs(x_mean_array=None, x_std_array=None, x_min_array=None, share_index=True, noise_type='chi2')[source]
By default, synthesize radiometer noise based on a chi-squared distribution. Alternately, can generate pure Gaussian noise.
If no arrays are specified from which to sample, noise samples will be drawn from saved GBT C-Band observations at (dt, df) = (1.4 s, 1.4 Hz) resolution, from frames of shape (tchans, fchans) = (32, 1024). These sample noise parameters consist of 126500 samples for mean, std, and min of each observation.
Specifying
noise_type='chi2'
will only usex_mean_array
(if provided), and ignore other parameters. Specifying noise_type=’gaussian’ will use all arrays (if provided).Note: this method will attempt to scale the noise parameters to match self.dt and self.df. This assumes that the observation data products are not normalized by the FFT length used to construct them.
- Parameters:
x_mean_array (ndarray, optional) – Array of potential means
x_std_array (ndarray, optional) – Array of potential standard deviations
x_min_array (ndarray, optional) – Array of potential minimum values
share_index (bool, optional, default: True) – Whether to select noise parameters from the same index across each provided array. If True, then each array must be the same length.
noise_type ({"chi2", "gaussian", "normal"}, default: "chi2") – Distribution to use for synthetic noise
- Returns:
noise – Array of synthetic noise
- Return type:
ndarray
- add_signal(path, t_profile, f_profile, bp_profile=None, bounding_f_range=None, integrate_path=False, integrate_t_profile=False, integrate_f_profile=False, doppler_smearing=False, t_subsamples=10, f_subsamples=10, smearing_subsamples=10)[source]
Generate synthetic signal.
Add a synethic signal using given path in time-frequency domain and brightness profiles in time and frequency directions.
- Parameters:
path (function, np.ndarray, list, float) – Function in time that returns frequencies, or provided array or single value of frequencies for the center of the signal at each time sample
t_profile (function, np.ndarray, list, float) – Time profile: function in time that returns an intensity (scalar), or provided array or single value of intensities at each time sample
f_profile (function) – Frequency profile: function in frequency that returns an intensity (scalar), relative to the signal frequency within a time sample. Note that unlike the other parameters, this must be a function
bp_profile (function, np.ndarray, list, float, optional) – Bandpass profile: function in frequency that returns a relative intensity (scalar, between 0 and 1), or provided array or single value of relative intensities at each frequency sample
bounding_f_range (tuple) – Tuple (bounding_min, bounding_max) that constrains the computation of the signal to only a range in frequencies
integrate_path (bool, optional) – Option to average path along time to get a more accurate frequency position in t-f space. Note that this option only makes sense if the provided path can be evaluated at the sub frequency sample level (e.g. as opposed to returning a pre-computed array of frequencies of length
tchans
). Makest_subsamples
calculations per time sample.integrate_t_profile (bool, optional) – Option to integrate
t_profile
in the time direction. Note that this option only makes sense if the providedt_profile
can be evaluated at the sub time sample level (e.g. as opposed to returning an array of intensities of lengthtchans
). Makest_subsamples
calculations per time sample.integrate_f_profile (bool, optional) – Option to integrate
f_profile
in the frequency direction. Makesf_subsamples
calculations per time sample.doppler_smearing (bool, optional) – Option to numerically “Doppler smear” spectral power over frequency bins. At time t, averages
smearing_subsamples
copies of the signal centered at evenly spaced center frequencies between times t and t+1. This causes the effective drop in power when the signal crosses multiple bins.t_subsamples (int, optional) – Number of bins for integration in the time direction, using Riemann sums. Default is 10.
f_subsamples (int, optional) – Number of bins for integration in the frequency direction, using Riemann sums. Default is 10.
smearing_subsamples (int, optional) – Number of steps for averaging evenly spaced copies of the signal between center frequencies at times t and t+1. Default is 10.
- Returns:
signal – Two-dimensional NumPy array containing synthetic signal data
- Return type:
ndarray
Examples
Here’s an example that creates a linear Doppler-drifted signal with chi-squared noise with sampled parameters:
>>> from astropy import units as u >>> import setigen as stg >>> fchans = 1024 >>> tchans = 32 >>> df = 2.7939677238464355*u.Hz >>> dt = tsamp = 18.253611008*u.s >>> fch1 = 6095.214842353016*u.MHz >>> frame = stg.Frame(fchans=fchans, tchans=tchans, df=df, dt=dt, fch1=fch1) >>> noise = frame.add_noise(x_mean=10) >>> signal = frame.add_signal(stg.constant_path(f_start=frame.get_frequency(200), drift_rate=2*u.Hz/u.s), stg.constant_t_profile(level=frame.get_intensity(snr=30)), stg.gaussian_f_profile(width=40*u.Hz), stg.constant_bp_profile(level=1))
Saving the noise and signals individually may be useful depending on the application, but the combined data can be accessed via frame.get_data(). The synthetic signal can then be visualized and saved within a Jupyter notebook using:
>>> %matplotlib inline >>> import matplotlib.pyplot as plt >>> fig = plt.figure(figsize=(10, 6)) >>> frame.plot() >>> plt.savefig('image.png', bbox_inches='tight') >>> plt.show()
To run within a script, simply exclude the first line:
%matplotlib inline
.
- add_constant_signal(f_start, drift_rate, level, width, f_profile_type='sinc2', doppler_smearing=False)[source]
A wrapper around add_signal() that injects a constant intensity, constant drift_rate signal into the frame.
- Parameters:
f_start (astropy.Quantity) – Starting signal frequency
drift_rate (astropy.Quantity) – Signal drift rate, in units of frequency per time
level (float) – Signal intensity
width (astropy.Quantity) – Signal width in frequency units
f_profile_type ({"sinc2", "box", "gaussian", "lorentzian", "voigt}, default: "sinc2") – Signal spectral profile
doppler_smearing (bool, optional, default: False) – Option to numerically “Doppler smear” spectral power over frequency bins. At time t, averages
drift_rate / frame.unit_drift_rate
copies of the signal centered at evenly spaced center frequencies between times t and t+1. This causes the effective drop in power when the signal crosses multiple bins.
- Returns:
signal – Two-dimensional NumPy array containing synthetic signal data
- Return type:
ndarray
- get_intensity(snr)[source]
Calculate intensity from SNR, based on estimates of the noise in the frame.
Note that there must be noise present in the frame for this to make sense.
- get_snr(intensity)[source]
Calculate SNR from intensity.
Note that there must be noise present in the frame for this to make sense.
- plot(*args, **kwargs)[source]
Plot frame spectrogram data.
- Parameters:
frame (Frame) – Frame to plot
xtype ({"fmid", "fmin", "f", "px"}, default: "fmid") – Types of axis labels, particularly the x-axis. “px” puts axes in units of pixels. The others are all in frequency: “fmid” shows frequencies relative to the central frequency, “fmin” is relative to the minimum frequency, and “f” is absolute frequency.
db (bool, default: True) – Option to convert intensities to dB
colorbar (bool, default: True) – Whether to display colorbar
label (bool, default: False) – Option to place target name as a label in plot
minor_ticks (bool, default: False) – Option to include minor ticks on both axes
grid (bool, default: False) – Option to overplot grid from major ticks
- Returns:
p – Spectrogram axes object
- Return type:
- get_slice(l, r)[source]
Slice frame data with left and right index bounds.
- Parameters:
l (int) – Left bound
r (int) – Right bound
- Returns:
s_fr – Sliced frame
- Return type:
- integrate(axis='t', mode='mean', normalize=False)[source]
Integrate along either time (‘t’, 0) or frequency (‘f’, 1) axes, to create spectra or time series data. Mode is either ‘mean’ or ‘sum’.
- Parameters:
data (Frame, or 2D ndarray) – Input frame or Numpy array
axis (int or str) – Axis over which to integrate; time (‘t’, 0) or frequency (‘f’, 1)
mode ({"mean", "sum"}, default: "mean") – Integration mode
normalize (bool) – Whether to normalize integrated array to mean 0, std 1
- Returns:
output – Integrated product
- Return type:
ndarray
- get_waterfall()[source]
Return current frame as a Waterfall object. Note: some filterbank metadata may not be accurate anymore, depending on prior frame manipulations.
- check_waterfall()[source]
If an associated Waterfall object exists, update and return it. Otherwise, return None. Useful to chain with
setigen.Frame.from_data()
if manipulating completely synthetic data.
- classmethod load_pickle(filename)[source]
Load Frame object from a pickled file (.pickle), created with
save_pickle()
.
- setigen.frame.params_from_backend(obs_length=300, sample_rate=3000000000.0, num_branches=1024, fftlength=1048576, int_factor=51)[source]
Return frame parameters calculated from data backend characteristics.
- Parameters:
obs_length (float, optional) – Length of observation in seconds
sample_rate (float, optional) – Physical sample rate, in Hz, for collecting real voltage data
num_branches (int, optional) – Number of PFB branches. Note that this corresponds to
num_branches / 2
coarse channels.fftlength (int, optional) – FFT length to be used in fine channelization
int_factor (int, optional) – Integration factor used in fine channelization. Determines
tchans
.
- Returns:
param_dict – Dictionary of parameters
- Return type:
dict
setigen.cadence module
- class setigen.cadence.Cadence(frame_list=None, t_slew=0, t_overwrite=False)[source]
Bases:
MutableSequence
A class for organizing cadences of Frame objects.
- Parameters:
frame_list (list of Frames, optional) – List of frames to be included in the cadence
t_slew (float, optional) – Slew time between frames
t_overwrite (bool, optional) – Option to overwrite time data in frames to enforce slew time spacing
- property t_start
- property fch1
- property ascending
- property fmin
- property fmax
- property fmid
- property df
- property dt
- property fchans
- property tchans
- property obs_range
- overwrite_times()[source]
Overwrite the starting time and time arrays used to compute signals in each frame of the cadence, using
slew_time
(s) to space between frames.
- property slew_times
Compute slew times in between frames.
- add_signal(*args, **kwargs)[source]
Add signal to each frame in the cadence. Arguments are passed through to
add_signal()
.
- plot(*args, **kwargs)[source]
Plot cadence as a multi-panel figure.
- Parameters:
cadence (Cadence) – Cadence to plot
xtype ({"fmid", "fmin", "f", "px"}, default: "fmid") – Types of axis labels, particularly the x-axis. “px” puts axes in units of pixels. The others are all in frequency: “fmid” shows frequencies relative to the central frequency, “fmin” is relative to the minimum frequency, and “f” is absolute frequency.
db (bool, default: True) – Option to convert intensities to dB
slew_times (bool, default: False) – Option to space subplots vertically proportional to slew times
colorbar (bool, default: True) – Whether to display colorbar
labels (bool, default: True) – Option to place target name as a label in each subplot
title (bool, default: False) – Option to place first source name as the figure title
minor_ticks (bool, default: False) – Option to include minor ticks on both axes
grid (bool, default: False) – Option to overplot grid from major ticks
- Returns:
axs (matplotlib.axes.Axes) – Axes subplots
cax (matplotlib.axes.Axes) – Colorbar axes, if created
- class setigen.cadence.OrderedCadence(frame_list=None, order='ABACAD', t_slew=0, t_overwrite=False)[source]
Bases:
Cadence
A class that extends Cadence for imposing a cadence order, such as “ABACAD” or “ABABAB”. Order labels are added to each frame’s metadata with the
order_label
keyword.- Parameters:
frame_list (list of Frames, optional) – List of frames to be included in the cadence
order (str, optional) – Cadence order, expressed as a string of letters (e.g. “ABACAD”)
t_slew (float, optional) – Slew time between frames
t_overwrite (bool, optional) – Option to overwrite time data in frames to enforce slew time spacing
setigen.frame_utils module
- setigen.frame_utils.array(fr)[source]
Return a Numpy array for input frame or Numpy array.
- Parameters:
fr (Frame, or 2D ndarray) – Input frame or Numpy array
- Returns:
data – Data array
- Return type:
ndarray
- setigen.frame_utils.integrate(fr, axis='t', mode='mean', normalize=False)[source]
Integrate along either time (‘t’, 0) or frequency (‘f’, 1) axes, to create spectra or time series data. Mode is either ‘mean’ or ‘sum’.
- Parameters:
fr (Frame, or 2D ndarray) – Input frame or Numpy array
axis (int or str) – Axis over which to integrate; time (‘t’, 0) or frequency (‘f’, 1)
mode ({"mean", "sum"}, default: "mean") – Integration mode
normalize (bool) – Option to normalize integrated array to mean 0, std 1
- Returns:
output – Integrated product
- Return type:
ndarray
setigen.dedrift module
setigen.split_utils module
- setigen.split_utils.split_waterfall_generator(waterfall_fn, fchans, tchans=None, f_shift=None)[source]
Create a generator that returns smaller Waterfall objects by ‘splitting’ an input filterbank file according to the number of frequency samples.
Since this function only loads in data in chunks according to fchans, it handles very large observations well. Specifically, it will not attempt to load all the data into memory before splitting, which won’t work when the data is very large anyway.
- Parameters:
waterfall_fn (str) – Filterbank filename with .fil extension
fchans (int) – Number of frequency samples per new filterbank file
tchans (int, optional) – Number of time samples to select - will default from start of observation. If None, just uses the entire integration time
f_shift (int, optional) – Number of samples to shift when splitting filterbank. If None, defaults to
f_shift=fchans
so that there is no overlap between new filterbank files
- Returns:
waterfall – A blimpy Waterfall object containing a smaller section of the data
- Return type:
Waterfall
- setigen.split_utils.split_fil(waterfall_fn, output_dir, fchans, tchans=None, f_shift=None)[source]
Create a set of new filterbank files by ‘splitting’ an input filterbank file according to the number of frequency samples.
- Parameters:
waterfall_fn (str) – Filterbank filename with .fil extension
output_dir (str) – Directory for new filterbank files
fchans (int) – Number of frequency samples per new filterbank file
tchans (int, optional) – Number of time samples to select - will default from start of observation. If None, just uses the entire integration time
f_shift (int, optional) – Number of samples to shift when splitting filterbank. If None, defaults to
f_shift=fchans
so that there is no overlap between new filterbank files
- Returns:
split_fns – List of new filenames
- Return type:
list of str
- setigen.split_utils.split_array(data, f_sample_num=None, t_sample_num=None, f_shift=None, t_shift=None, f_trim=False, t_trim=False)[source]
Split NumPy arrays into a list of smaller arrays according to limits in frequency and time. This doesn’t reduce/combine data, it simply cuts the data into smaller chunks.
- Parameters:
data (ndarray) – Time-frequency data
- Returns:
split_data – List of new time-frequency data frames
- Return type:
list of ndarray
setigen.normalize module
- setigen.normalize.sigma_clip_norm(fr, axis=None, as_data=None)[source]
Normalize data by subtracting out noise background, determined by sigma clipping.
- Parameters:
- Returns:
Returns normalized data in the same type as f
- Return type:
n_data
- setigen.normalize.sliding_norm(data, cols=0, exclude=0.0, db=False, use_median=False)[source]
Normalize data per frequency channel so that the noise level in data is controlled; using mean or median filter.
Uses a sliding window to calculate mean and standard deviation to preserve non-drifted signals. Excludes a fraction of brightest pixels to better isolate noise.
- Parameters:
data (ndarray) – Time-frequency data
cols (int) – Number of columns on either side of the current frequency bin. The width of the sliding window is thus 2 * cols + 1
exclude (float, optional) – Fraction of brightest samples in each frequency bin to exclude in calculating mean and standard deviation
db (bool, optional) – Convert values to decibel equivalents before normalization
use_median (bool, optional) – Use median and median absolute deviation instead of mean and standard deviation
- Returns:
normalized_data – Normalized data
- Return type:
ndarray
- setigen.normalize.blimpy_clip(data, exclude=0)[source]
Naive clipping method by excluding top fraction of intensities.
- Parameters:
data (ndarray) – Time-frequency data
exclude (float, optional) – Fraction of brightest samples to exclude in calculating mean and standard deviation
- Returns:
data – Clipped data
- Return type:
ndarray
setigen.distributions module
- setigen.distributions.fwhm(sigma)[source]
Get full width at half maximum (FWHM) for a provided sigma / standard deviation, assuming a Gaussian distribution.
- setigen.distributions.truncated_gaussian(x_mean, x_std, x_min, shape, seed=None)[source]
Sample from a normal distribution, but enforces a minimum value.
- setigen.distributions.chi2(x_mean, chi2_df, shape, seed=None)[source]
Chi-squared distribution centered at a specific mean.
- Parameters:
x_mean (float) –
chi2_df (int) – Degrees of freedom for chi-squared
shape (list) – Shape of output noise array
- Returns:
dist – Array of chi-squared noise
- Return type:
ndarray
setigen.waterfall_utils module
- setigen.waterfall_utils.max_freq(waterfall)[source]
Return central frequency of the highest-frequency bin in a .fil file.
- Parameters:
waterfall (str or Waterfall) – Name of filterbank file or Waterfall object
- Returns:
fmax – Maximum frequency in data
- Return type:
float
- setigen.waterfall_utils.min_freq(waterfall)[source]
Return central frequency of the lowest-frequency bin in a .fil file.
- Parameters:
waterfall (str or Waterfall) – Name of filterbank file or Waterfall object
- Returns:
fmin – Minimum frequency in data
- Return type:
float
- setigen.waterfall_utils.get_data(waterfall, db=False)[source]
Get time-frequency data from filterbank file as a 2d NumPy array.
Note: when multiple Stokes parameters are supported, this will have to be expanded.
- Parameters:
waterfall (str or Waterfall) – Name of filterbank file or Waterfall object
- Returns:
data – Time-frequency data
- Return type:
ndarray
setigen.sample_from_obs module
- setigen.sample_from_obs.sample_gaussian_params(x_mean_array, x_std_array, x_min_array=None, seed=None)[source]
Sample Gaussian parameters (mean, std, min) from provided arrays.
Typical usage would be for select Gaussian noise properties for injection into data frames.
- Parameters:
x_mean_array (ndarray) – Array of potential means
x_std_array (ndarray) – Array of potential standard deviations
x_min_array (ndarray, optional) – Array of potential minimum values
seed (None, int, Generator, optional) – Random seed or seed generator
- Returns:
x_mean – Selected mean of distribution
x_std – Selected standard deviation of distribution
x_min – If x_min_array provided, selected minimum of distribution
- setigen.sample_from_obs.get_parameter_distributions(waterfall_fn, fchans, tchans=None, f_shift=None)[source]
Calculate parameter distributions for the mean, standard deviation, and minimum of split filterbank data from real observations.
- Parameters:
waterfall_fn (str) – Filterbank filename with .fil extension
fchans (int) – Number of frequency samples per new filterbank file
tchans (int, optional) – Number of time samples to select - will default from start of observation. If None, just uses the entire integration time
f_shift (int, optional) – Number of samples to shift when splitting filterbank. If None, defaults to
f_shift=f_window
so that there is no overlap between new filterbank files
- Returns:
x_mean_array – Distribution of means calculated from observations
x_std_array – Distribution of standard deviations calculated from observations
x_min_array – Distribution of minimums calculated from observations
- setigen.sample_from_obs.get_mean_distribution(waterfall_fn, fchans, tchans=None, f_shift=None)[source]
Calculate parameter distributions for the mean of split filterbank frames from real observations.
- Parameters:
waterfall_fn (str) – Filterbank filename with .fil extension
fchans (int) – Number of frequency samples per new filterbank file
tchans (int, optional) – Number of time samples to select - will default from start of observation. If None, just uses the entire integration time
f_shift (int, optional) – Number of samples to shift when splitting filterbank. If None, defaults to
f_shift=f_window
so that there is no overlap between new filterbank files
- Returns:
Distribution of means calculated from observations
- Return type:
x_mean_array
setigen.unit_utils module
This module contains a couple unit conversion utilities used in frame.Frame.
In general, we rely on astropy units for conversions, and note that float values are assumed to be in SI units (e.g. Hz, s).
- setigen.unit_utils.cast_value(value, unit)[source]
If value is already an astropy Quantity, then cast it into the desired unit. Otherwise, value is assumed to be a float and converted directly to the desired unit.
- setigen.unit_utils.get_value(value, unit=None)[source]
This function converts a value, which may be a float or astropy Quantity, into a float (in terms of a desired unit).
If we know that value is an astropy Quantity, then grabbing the value is simple (and we can cast this to a desired unit, if we need to change this.
If value is already a float, it simply returns value.