ShiftAndResample

class baseband.tasks.ShiftAndResample(ih, shift, offset=None, whence='start', *, lo=None, pad=64, samples_per_frame=None)[source] [edit on github]

Bases: baseband_tasks.convolution.Convolve

Shift and optionally resample a stream in time.

The shift is added to the sample times, and the stream is optionally resampled to ensure a sample falls on the given offset.

If the shift corresponds to a time delay, and the signal was recorded after mixing, the phases should be adjusted, which can be done by passing in the local oscillator frequency. For an upper sideband, the phase correction is,

\[\phi = - \tau f_{lo}.\]

For the lower sideband, the rotation is in the opposite direction.

Parameters
ihtask or baseband stream reader

Input data stream, with time as the first axis.

shiftfloat, float array-like, or Quantity

Amount by which to shift samples in time, as (float) samples, or a quantity with units of time. Should broadcast to the sample shape.

offsetfloat, Quantity, or Time

Offset that the output stream should include. Can be an absolute time, or a (float) number of samples or time offset relative to the start of the underlying stream. The default of None implies that the output stream is free to adjust. Hence, if a single shift is given, all that will happen is a change in start_time. To ensure the grid stays fixed, pass in 0.

whence{0, 1, 2, ‘start’, ‘current’, or ‘end’}, optional

Like regular seek, the offset is taken to be from the start if whence=0 (default), from the current position if 1, and from the end if 2. One can alternativey use ‘start’, ‘current’, or ‘end’ for 0, 1, or 2, respectively. Ignored if offset is a time.

loQuantity, or None, optional

Local oscillator frequency. Should be passed in when the signal was mixed down before sampling and the shift should be interpreted as a time delay (rather than a clock correction). For raw data, this can usually just be if.frequency. But for channelized data, the actual frequency needs to be passed in. If data were recorded without mixing (like for CHIME), no phase correction is necessary and one should pass in None. Default: None.

padint, optional

Padding to apply on each side when shifting data. This sets the size of the sinc function which which the data is convolved (see below). Numerical tests suggest that with the default of 64, the accuracy is better than 0.1%.

samples_per_frameint, optional

Number of resampled samples which should be produced in one go. The number of input samples used will be larger by 2*pad+1. If not given, works on the larger of the samples per frame from If not given, the larger of the sampler per frame in the underlying stream or 14 times the padding (to ensure ~87.5% efficiency).

See also

Resample

resample a stream to a new grid, without time shifts

ShiftSamples

shift channels by integer number of samples

TimeDelay

delay a complex stream, changing phases but no resampling

baseband_tasks.base.SetAttribute

change start time without resampling

Notes

The precision of the shifting and resampling is controlled by pad, as it sets the length of the response, which consists of a sinc function combined with a Hann window,

\[\begin{split}R(x) &= {\rm sinc}(x-s) \cos^2(\pi(x-s)/L),\\ -{\rm pad} &\le x \le {\rm pad},\quad L = 2{\rm pad} + 2\end{split}\]

Here, \(s\) is the fractional pixel offset. Note that \(L\) is chosen such that \(R(\pm{\rm pad})\) is still non-zero.

The convolution is done in the Fourier domain, and the Hann window, combined with with the padding done when reading, ensures there is no aliasing between the front and the back of each frame. There is, however, aliasing near the Nyquist frequency. For most recorded signals, this will not be important, as the band-pass filter will likely have ensured there is very little signal near the band edges. For channelized data, however, it may be more problematic and some pre-filtering stage may be necessary.

Attributes Summary

closed

complex_data

dtype

Data type of the output.

meta

ndim

Number of dimensions of the output.

offset

sample_rate

Number of complete samples per second.

sample_shape

Shape of a complete sample.

samples_per_frame

Number of samples per frame of data.

shape

Shape of the output.

size

Number of component samples in the output.

start_time

Start time of the output.

stop_time

Time at the end of the output, just after the last sample.

time

Time of the sample pointer's current offset in the output.

Methods Summary

close()

Close task.

read([count, out])

Read a number of complete samples.

seek(offset[, whence])

Change the sample pointer position.

task(data)

tell([unit])

Current offset in the file.

Attributes Documentation

closed = False
complex_data
dtype

Data type of the output.

meta
ndim

Number of dimensions of the output.

offset = 0
sample_rate

Number of complete samples per second.

sample_shape

Shape of a complete sample.

samples_per_frame

Number of samples per frame of data.

For compatibility with file readers, to help indicate what a nominal chunk of data is.

shape

Shape of the output.

size

Number of component samples in the output.

start_time

Start time of the output.

See also time and stop_time.

stop_time

Time at the end of the output, just after the last sample.

See also start_time and time.

time

Time of the sample pointer’s current offset in the output.

See also start_time and stop_time.

Methods Documentation

close() [edit on github]

Close task.

Note that this does not explicitly close the underlying source; instead, it just deletes the reference to it.

read(count=None, out=None) [edit on github]

Read a number of complete samples.

Parameters
countint or None, optional

Number of complete samples to read. If None (default) or negative, the number of samples left. Ignored if out is given.

outNone or array, optional

Array to store the samples in. If given, count will be inferred from the first dimension; the remaining dimensions should equal sample_shape.

Returns
outndarray of float or complex

The first dimension is sample-time, and the remaining ones are as given by sample_shape.

seek(offset, whence=0) [edit on github]

Change the sample pointer position.

This works like a normal filehandle seek, but the offset is in samples (or a relative or absolute time).

Parameters
offsetint, Quantity, or Time

Offset to move to. Can be an (integer) number of samples, an offset in time units, or an absolute time. For the latter two, the pointer will be moved to the nearest integer sample.

whence{0, 1, 2, ‘start’, ‘current’, or ‘end’}, optional

Like regular seek, the offset is taken to be from the start if whence=0 (default), from the current position if 1, and from the end if 2. One can alternativey use ‘start’, ‘current’, or ‘end’ for 0, 1, or 2, respectively. Ignored if offset is a time.

task(data) [edit on github]
tell(unit=None) [edit on github]

Current offset in the file.

Parameters
unitUnit or str, optional

Time unit the offset should be returned in. By default, no unit is used, i.e., an integer enumerating samples is returned. For the special string ‘time’, the absolute time is calculated.

Returns
offsetint, Quantity, or Time

Offset in current file (or time at current position).