ch_util.timing

Tools for timing jitter and delay corrections.

This module contains tools for using noise sources to correct timing jitter and timing delay.

Example

The function construct_delay_template() generates a delay template from measurements of the visibility between noise source inputs, which can be used to remove the timing jitter in other data.

The user seldom needs to work with construct_delay_template() directly and can instead use several high-level functions and containers that load the timing data, derive the timing correction using construct_delay_template(), and then enable easy application of the timing correction to other data.

For example, to load the timing data and derive the timing correction from a list of timing acquisition files (i.e., YYYYMMSSTHHMMSSZ_chimetiming_corr), use the following:

`tdata = TimingData.from_acq_h5(timing_acq_filenames)`

This results in a andata.CorrData object that has additional methods avaiable for applying the timing correction to other data. For example, to obtain the complex gain for some freq, input, and time that upon multiplication will remove the timing jitter, use the following:

`tgain, tweight = tdata.get_gain(freq, input, time)`

To apply the timing correction to the visibilities in an andata.CorrData object called data, use the following:

`tdata.apply_timing_correction(data)`

The timing acquisitions must cover the span of time that you wish to correct. If you have a list of data acquisition files and would like to obtain the appropriate timing correction by searching the archive for the corresponding timing acquisitons files, then use:

`tdata = load_timing_correction(data_acq_filenames_full_path)`

To print a summary of the timing correction, use:

`print(tdata)`

Functions

construct_delay_template(data[, ...])

Construct a relative time delay template.

eigen_decomposition(vis, flag)

Eigenvalue decomposition of the visibility matrix.

fit_poly_to_phase(freq, resp, resp_error[, ...])

Fit complex data versus frequency to a model consisting of a polynomial in phase.

load_timing_correction(files[, start, stop, ...])

Find and load the timing correction for a list of corr acquisition files.

map_input_to_noise_source(inputs, noise_sources)

Find the appropriate noise source to use to correct the phase of each input.

model_poly_phase(freq, *param)

Evaluate a polynomial model for the phase.

Classes

TimingCorrection([h5_data])

Container that holds a timing correction.

TimingData([h5_data])

Subclass of andata.CorrData for timing data.

TimingInterpolator(x, y[, weight, flag, ...])

Interpolation that is aware of flagged data and weights.

class ch_util.timing.TimingCorrection(h5_data=None, **kwargs)[source]

Bases: BaseData

Container that holds a timing correction.

Provides methods for applying that correction to other datasets.

Used to pick which subclass to instantiate based on attributes in data.

property alpha

Provide convenience access to the alpha array.

property amp_to_delay

This is the conversion from noise source amplitude variations to delay variations.

apply_timing_correction(timestream, copy=False, **kwargs)[source]

Apply the timing correction to another visibility dataset.

This method uses the get_gain or get_stacked_tau method, depending on whether or not the visibilities have been stacked. It acccepts and passes along keyword arguments for those method.

Parameters:
  • timestream (andata.CorrData / equivalent or np.ndarray[nfreq, nprod, ntime]) – If timestream is an np.ndarray containing the visiblities, the you must also pass the corresponding freq, prod, input, and time axis as kwargs. Otherwise these quantities are obtained from the attributes of CorrData. If the visibilities have been stacked, then you must additionally pass the stack and reverse_stack axis as kwargs, and (optionally) the input flags.

  • copy (bool) – Create a copy of the input visibilities. Apply the timing correction to the copy and return it, leaving the original untouched. Default is False.

  • freq (np.ndarray[nfreq, ]) – Frequency in MHz. Must be passed as keyword argument if timestream is an np.ndarray.

  • prod (np.ndarray[nprod, ]) – Product map. Must be passed as keyword argument if timestream is an np.ndarray.

  • time (np.ndarray[ntime, ]) – Unix time. Must be passed as keyword argument if timestream is an np.ndarray.

  • input (np.ndarray[ninput, ] of dtype=('chan_id', 'correlator_input')) – Input axis. Must be passed as keyword argument if timestream is an np.ndarray.

  • stack (np.ndarray[nstack, ]) – Stack axis. Must be passed as keyword argument if timestream is an np.ndarray and the visibilities have been stacked.

  • reverse_stack (np.ndarray[nprod, ] of dtype=('stack', 'conjugate')) – The index of the stack axis that each product went into. Typically found in reverse_map[‘stack’] attribute. Must be passed as keyword argument if timestream is an np.ndarray and the visibilities have been stacked.

  • input_flags (np.ndarray [ninput, ntime]) – Array indicating which inputs were good at each time. Non-zero value indicates that an input was good. Optional. Only used for stacked visibilities.

Returns:

  • If copy == True

    visnp.ndarray[nfreq, nprod(nstack), ntime]

    New set of visibilities with timing correction applied.

  • else

    None

    Correction is applied to the input visibility data. Also, if timestream is an andata.CorrData instance and the gain dataset exists, then it will be updated with the complex gains that have been applied.

property coeff_alpha

Provide convenience access to the coeff_alpha array.

property coeff_tau

Provide convenience access to the coeff_tau array.

delete_coeff()[source]

Stop using coefficients to construct timing correction.

Calling this method will delete the coeff_tau, coeff_alpha, and reference_noise_source datasets if they exist.

property freq

Provide convenience access to the frequency bin centres.

classmethod from_dict(**kwargs)[source]

Instantiate a TimingCorrection object.

Parameters:
  • freq (np.ndarray[nfreq, ] of dtype=('centre', 'width')) – Frequencies in MHz that were used to construct the timing correction.

  • noise_source (np.ndarray[nsource,] of dtype=('chan_id', 'correlator_input')) – Correlator inputs that were used to construct the timing correction.

  • input (np.ndarray[ninput, ] of dtype=('chan_id', 'correlator_input')) – Correlator inputs to which the timing correction will be applied.

  • time (np.ndarray[ntime, ]) – Unix time.

  • param (np.ndarray[nparam, ]) – Parameters of the model fit to the static phase versus frequency.

  • tau (np.ndarray[nsource, ntime]) – The actual timing correction, which is the relative delay of each of the noise source inputs with respect to a reference input versus time.

  • weight_tau (np.ndarray[nsource, ntime]) – Estimate of the uncertainty (inverse variance) on the timing correction.

  • static_phi (np.ndarray[nfreq, nsource]) – The phase that was subtracted from each frequency and input prior to fitting for the timing correction. This is necessary to remove the approximately static ripple pattern caused by reflections.

  • weight_static_phi (np.ndarray[nfreq, nsource]) – Inverse variance on static_phi.

  • static_phi_fit (np.ndarray[nparam, nsource]) – Best-fit parameters of a fit to the static phase versus frequency for each of the noise source inputs.

  • alpha (np.ndarray[nsource, ntime]) – The coefficient of the spectral model of the amplitude variations of each of the noise source inputs versus time.

  • weight_alpha (np.ndarray[nsource, ntime]) – Estimate of the uncertainty (inverse variance) on the amplitude coefficients.

  • static_amp (np.ndarray[nfreq, nsource]) – The amplitude that was subtracted from each frequency and input prior to fitting for the amplitude variations. This is necessary to remove the approximately static ripple pattern caused by reflections.

  • weight_static_amp (np.ndarray[nfreq, nsource]) – Inverse variance on static_amp.

  • num_freq (np.ndarray[nsource, ntime]) – The number of frequencies used to determine the delay and alpha quantities. If num_freq is 0, then that time is ignored when deriving the timing correction.

  • coeff_tau (np.ndarray[ninput, nsource]) – If coeff is provided, then the timing correction applied to a particular input will be the linear combination of the tau correction from the noise source inputs, with the coefficients set by this array.

  • coeff_alpha (np.ndarray[ninput, nsource]) – If coeff is provided, then the timing correction applied to a particular input will be adjusted by the linear combination of the alpha correction from the noise source inputs, with the coefficients set by this array.

  • reference_noise_source (np.ndarray[ninput]) – The noise source input that was used as reference when fitting coeff_tau.

get_alpha(timestamp, interp='linear', extrap_limit=None)[source]

Return the amplitude variation for each noise source at the requested times.

Uses the TimingInterpolator to interpolate to the requested times.

Parameters:
  • timestamp (np.ndarray[ntime,]) – Unix timestamp.

  • interp (string) – Method to interpolate over time. Options include ‘linear’, ‘nearest’, ‘zero’, ‘slinear’, ‘quadratic’, ‘cubic’, ‘previous’, and ‘next’.

  • extrap_limit (float) – Do not extrapolate the underlying data beyond its boundaries by this amount in seconds. Default is 2 integrations.

Returns:

  • alpha (np.ndarray[nsource, ntime]) – Amplitude coefficient as a function of time for each of the noise sources.

  • weight (np.ndarray[nsource, ntime]) – The uncertainty on the amplitude coefficient, expressed as an inverse variance.

get_gain(freq, inputs, timestamp, **kwargs)[source]

Return the complex gain for the requested frequencies, inputs, and times.

Multiplying the visibilities by the outer product of these gains will remove the fluctuations in phase due to timing jitter. This method uses the get_tau method. It acccepts and passes along keyword arguments for that method.

Parameters:
  • freq (np.ndarray[nfreq, ]) – Frequency in MHz.

  • inputs (np.ndarray[ninput, ]) – Must contain ‘correlator_input’ field.

  • timestamp (np.ndarray[ntime, ]) – Unix timestamps.

Returns:

  • gain (np.ndarray[nfreq, ninput, ntime]) – Complex gain. Multiplying the visibilities by the outer product of this vector at a given time and frequency will correct for the timing jitter.

  • weight (np.ndarray[nfreq, ninput, ntime]) – Uncerainty on the gain expressed as an inverse variance.

get_stacked_alpha(timestamp, inputs, prod, reverse_stack, input_flags=None, **kwargs)[source]

Return the stacked alphas for the noise source amplitude variations.

Averages the alphas from the noise source inputs that map to the set of redundant baseline included in each stacked visibility. If input_flags is provided, then the bad inputs that were excluded from the stack are also excluded from the alpha template averaging. This method can be used to generate a stacked alpha template that can be used to correct a stacked tau template for variations in the noise source distribution system. However, it is recommended that the tau template be corrected before stacking. This is accomplished by setting the amp_to_delay property prior to calling get_stacked_tau.

Parameters:
  • timestamp (np.ndarray[ntime,]) – Unix timestamp.

  • inputs (np.ndarray[ninput,]) – Must contain ‘correlator_input’ field.

  • prod (np.ndarray[nprod,]) – The products that were included in the stack. Typically found in the index_map[‘prod’] attribute of the andata.CorrData object.

  • reverse_stack (np.ndarray[nprod,] of dtype=('stack', 'conjugate')) – The index of the stack axis that each product went into. Typically found in reverse_map[‘stack’] attribute of the andata.CorrData.

  • input_flags (np.ndarray [ninput, ntime]) – Array indicating which inputs were good at each time. Non-zero value indicates that an input was good.

Returns:

alpha – Noise source amplitude variation as a function of time for each stacked visibility.

Return type:

np.ndarray[nstack, ntime]

get_stacked_tau(timestamp, inputs, prod, reverse_stack, input_flags=None, **kwargs)[source]

Return the delay for each stacked visibility at the requested time.

Averages the delays from the noise source inputs that map to the set of redundant baseline included in each stacked visibility. This yields the appropriate common-mode delay correction. If input_flags is provided, then the bad inputs that were excluded from the stack are also excluded from the delay template averaging.

Parameters:
  • timestamp (np.ndarray[ntime,]) – Unix timestamp.

  • inputs (np.ndarray[ninput,]) – Must contain ‘correlator_input’ field.

  • prod (np.ndarray[nprod,]) – The products that were included in the stack. Typically found in the index_map[‘prod’] attribute of the andata.CorrData object.

  • reverse_stack (np.ndarray[nprod,] of dtype=('stack', 'conjugate')) – The index of the stack axis that each product went into. Typically found in reverse_map[‘stack’] attribute of the andata.CorrData.

  • input_flags (np.ndarray [ninput, ntime]) – Array indicating which inputs were good at each time. Non-zero value indicates that an input was good.

Returns:

tau – Delay as a function of time for each stacked visibility.

Return type:

np.ndarray[nstack, ntime]

get_tau(timestamp, ignore_amp=False, interp='linear', extrap_limit=None)[source]

Return the delay for each noise source at the requested times.

Uses the TimingInterpolator to interpolate to the requested times.

Parameters:
  • timestamp (np.ndarray[ntime,]) – Unix timestamp.

  • ignore_amp (bool) – Do not apply a noise source based amplitude correction, even if one exists.

  • interp (string) – Method to interpolate over time. Options include ‘linear’, ‘nearest’, ‘zero’, ‘slinear’, ‘quadratic’, ‘cubic’, ‘previous’, and ‘next’.

  • extrap_limit (float) – Do not extrapolate the underlying data beyond its boundaries by this amount in seconds. Default is 2 integrations.

Returns:

  • tau (np.ndarray[nsource, ntime]) – Delay as a function of time for each of the noise sources.

  • weight (np.ndarray[nsource, ntime]) – The uncertainty on the delay, expressed as an inverse variance.

get_timing_correction(freq, timestamp, **kwargs)[source]

Return the phase correction from each noise source.

Assumes the phase correction scales with frequency nu as

phi = 2 pi nu tau

and uses the get_tau method to interpolate over time. It acccepts and passes along keyword arguments for that method.

Parameters:
  • freq (np.ndarray[nfreq, ]) – Frequency in MHz.

  • timestamp (np.ndarray[ntime, ]) – Unix timestamp.

Returns:

  • gain (np.ndarray[nfreq, nsource, ntime]) – Complex gain containing a pure phase correction for each of the noise sources.

  • weight (np.ndarray[nfreq, nsource, ntime]) – Uncerainty on the gain for each of the noise sources, expressed as an inverse variance.

property has_amplitude

Determine if this timing correction contains amplitude data.

property has_coeff_alpha

Indicates if there are valid coeff that map noise source alpha to inputs.

property has_coeff_tau

Indicates if there are valid coeff that map noise source tau to inputs.

property has_num_freq

Inidicates if there is a num_freq flag that identifies missing data.

property input

Provide convenience access to the correlator inputs.

property noise_source

Provide convenience access to the noise source inputs.

Note that in older versions of the timing correction, the noise_source axis does not exist. Instead, the equivalent quantity is labeled as input. Since the addition of the coeff dataset it has become necessary to distinguish between the noise source inputs from which the timing correction is derived and the correlator inputs to which the timing correction is applied.

property nsource

Provide convenience access to the number of noise source inputs.

property num_freq

Provide convenience access to the num_freq array.

property reference_noise_source

Return the index of the reference noise source.

search_input(inputs)[source]

Find inputs in the input axis.

Parameters:

inputs (np.ndarray[ninput,] of dtype=('chan_id', 'correlator_input'))

Returns:

index – Indices of the input axis that yield the requested inputs.

Return type:

np.ndarray[ninput,] of .int

set_coeff(coeff_tau, inputs, noise_source, coeff_alpha=None, reference_noise_source=None)[source]

Use coefficients to construct timing correction.

Setting the coefficients changes how the timing corretion for a particular correlator input is derived. Without coefficients, each input is matched to the timing correction from a single noise source input through the map_input_to_noise_source method. With coefficients, each input is a linear combination of the timing correction from all noise source inputs.

Parameters:
  • coeff_tau (np.ndarray[ninput, nsource]) – The timing correction applied to a particular input will be the linear combination of the tau correction from the noise source inputs, with the coefficients set by this array.

  • inputs (np.ndarray[ninput, ] of dtype=('chan_id', 'correlator_input')) – Correlator inputs to which the timing correction will be applied.

  • noise_source (np.ndarray[nsource,] of dtype=('chan_id', 'correlator_input')) – Correlator inputs that were used to construct the timing correction.

  • coeff_alpha (np.ndarray[ninput, nsource]) – The timing correction applied to a particular input will be adjusted by the linear combination of the alpha correction from the noise source inputs, with the coefficients set by this array.

  • reference_noise_source (np.ndarray[ninput,]) – For each input, the index into noise_source that was used as reference in the fit for coeff_tau.

set_global_reference_time(tref, window=0.0, interpolate=False, **kwargs)[source]

Normalize the delay and alpha template to the value at a single time.

Useful for referencing the template to the value at the time that you plan to calibrate.

Parameters:
  • tref (unix time) – Reference the templates to the values at this time.

  • window (float) – Reference the templates to the median value over a window (in seconds) around tref. If nonzero, this will override the interpolate keyword.

  • interpolate (bool) – Interpolate the delay template to time tref. Otherwise take the measured time nearest to tref. The get_tau method is use to perform the interpolation, and kwargs for that method will be passed along.

set_reference_time(tref, tstart, tend=None, tinit=None, tau_init=None, alpha_init=None, interpolate=False, **kwargs)[source]

Normalize the delay and alpha template to specific times.

Required if applying the timing correction to data that has already been calibrated.

Parameters:
  • tref (np.ndarray[nref]) – Reference the delays to the values at this unix time.

  • tstart (np.ndarray[nref]) – Begin transition to the reference delay at this unix time.

  • tend (np.ndarray[nref]) – Complete transition to the reference delay at this unix time.

  • tinit (float) – Use the delay at this time for the period before the first tstart. Takes prescendent over tau_init.

  • tau_init (np.ndarray[nsource]) – Use this delay for times before the first tstart. Must provide a value for each noise source input. If None, then will reference with respect to the average delay over the full time series.

  • alpha_init (np.ndarray[nsource]) – Use this alpha for times before the first tstart. Must provide a value for each noise source input. If None, then will reference with respect to the average alpha over the full time series.

  • interpolate (bool) – Interpolate the delay template to times tref. Otherwise take the measured times nearest to tref. The get_tau method is use to perform the interpolation, and kwargs for that method will be passed along.

property static_amp

Provide convenience access to the static_amp array.

property static_phi

Provide convenience access to the static_phi array.

property static_phi_fit

Provide convenience access to the static_phi_fit array.

summary()[source]

Provide a summary of the timing correction.

Returns:

summary – Contains useful information about the timing correction. Specifically contains for each noise source input the time averaged phase offset and delay. Also contains estimates of the variance in the timing for both the shortest and longest timescale probed by the underlying dataset. Meant to be joined with new lines and printed.

Return type:

list of strings

property tau

Provide convenience access to the tau array.

property weight_alpha

Provide convenience access to the weight_alpha array.

property weight_static_amp

Provide convenience access to the weight_static_amp array.

property weight_static_phi

Provide convenience access to the weight_static_phi array.

property weight_tau

Provide convenience access to the weight_tau array.

property zero_delay_noise_source

Return the index of the noise source with zero delay.

class ch_util.timing.TimingData(h5_data=None, **kwargs)[source]

Bases: CorrData, TimingCorrection

Subclass of andata.CorrData for timing data.

Automatically computes the timing correction when data is loaded and inherits the methods of TimingCorrection that enable the application of that correction to other datasets.

Used to pick which subclass to instantiate based on attributes in data.

classmethod from_acq_h5(acq_files, only_correction=False, **kwargs)[source]

Load a list of acquisition files and computes the timing correction.

Accepts and passes on all keyword arguments for andata.CorrData.from_acq_h5 and the construct_delay_template function.

Parameters:
  • acq_files (str or list of str) – Path to file(s) containing the timing data.

  • only_correction (bool) – Only return the timing correction. Do not return the underlying data from which that correction was derived.

Returns:

data

Return type:

TimingData or TimingCorrection

summary()[source]

Provide a summary of the timing data and correction.

Returns:

summary – Contains useful information about the timing correction and data. Includes the reduction in the standard deviation of the phase after applying the timing correction. This is presented as quantiles over frequency for each of the noise source products.

Return type:

list of strings

class ch_util.timing.TimingInterpolator(x, y, weight=None, flag=None, kind='linear', extrap_limit=None)[source]

Bases: object

Interpolation that is aware of flagged data and weights.

Flagged data is ignored during the interpolation. The weights from the data are propagated to obtain weights for the interpolated points.

Instantiate a callable TimingInterpolator object.

Parameters:
  • x (np.ndarray[nsample,]) – The points where the data was sampled. Must be monotonically increasing.

  • y (np.ndarray[..., nsample]) – The data to interpolate.

  • weight (np.ndarray[..., nsample]) – The uncertainty on the data, expressed as an inverse variance.

  • flag (np.ndarray[..., nsample]) – Boolean indicating if the data is to be included in the interpolation.

  • kind (str) – String that specifies the kind of interpolation. The value nearest, previous, next, and linear will use custom methods that propagate uncertainty to obtain the interpolated weights. The value zero, slinear, quadratic, and cubic will use spline interpolation from scipy.interpolation.interp1d and use the weight from the nearest point.

Returns:

interpolator – Callable that will interpolate the data that was provided to a new set of x values.

Return type:

TimingInterpolator

ch_util.timing.construct_delay_template(data, min_frac_kept=0.0, threshold=0.5, min_freq=420.0, max_freq=780.0, mask_rfi=False, max_iter_weight=None, check_amp=False, nsigma_amp=None, check_phi=True, nsigma_phi=None, nparam=2, static_phi=None, weight_static_phi=None, static_phi_fit=None, static_amp=None, weight_static_amp=None)[source]

Construct a relative time delay template.

Fits the phase of the cross-correlation between noise source inputs to a model that increases linearly with frequency.

Parameters:
  • data (andata.CorrData) –

    Correlation data. Must contain the following attributes:
    freq: np.ndarray[nfreq, ]

    Frequency in MHz.

    vis: np.ndarray[nfreq, nprod, ntime]

    Upper-triangle, product packed visibility matrix containing ONLY the noise source inputs.

    weight: np.ndarray[nfreq, nprod, ntime]

    Flag indicating the data points to fit.

    flags/frac_lost: np.ndarray[nfreq, ntime]

    Flag indicating the fraction of data lost. If provided, then data will be weighted by the fraction of data that remains when solving for the delay template.

  • min_frac_kept (float) – Do not include frequencies and times where the fraction of data that remains is less than this threshold. Default is 0.0.

  • threshold (float) – A (frequency, input) must pass the checks specified above more than this fraction of the time, otherwise it will be flaged as bad for all times. Default is 0.50.

  • min_freq (float) – Minimum frequency in MHz to include in the fit. Default is 420.

  • max_freq (float) – Maximum frequency in MHz to include in the fit. Default is 780.

  • mask_rfi (bool) – Mask frequencies that occur within known RFI bands. Note that the noise source data does not contain RFI, however the real-time pipeline does not distinguish between noise source inputs and sky inputs, and as a result will discard large amounts of data in these bands.

  • max_iter_weight (int) – The weight for each frequency is estimated from the variance of the residuals of the template fit from the previous iteration. Outliers are also flagged at each iteration with an increasingly aggresive threshold. This is the total number of times to iterate. Setting to 1 corresponds to linear least squares. Default is 1, unless check_amp or check_phi is True, in which case this defaults to the maximum number of thresholds provided.

  • check_amp (bool) – Do not fit frequencies and times where the residual amplitude is an outlier. Default is False.

  • nsigma_amp (list of float) – If check_amp is True, then residuals greater than this number of sigma will be considered an outlier. Provide a list containing the value to be used at each iteration. If the length of the list is less than max_iter_weight, then the last value in the list will be repeated for the remaining iterations. Default is [1000, 500, 200, 100, 50, 20, 10, 5].

  • check_phi (bool) – Do not fit frequencies and times where the residual phase is an outlier. Default is True.

  • nsigma_phi (list of float) – If check_phi is True, then residuals greater than this number of sigma will be considered an outlier. Provide a list containing the value to be used at each iteration. If the length of the list is less than max_iter_weight, then the last value in the list will be repeated for the remaining iterations. Default is [1000, 500, 200, 100, 50, 20, 10, 5].

  • nparam (int) – Number of parameters for polynomial fit to the time averaged phase versus frequency. Default is 2.

  • static_phi (np.ndarray[nfreq, nsource]) – Subtract this quantity from the noise source phase prior to fitting for the timing correction. If None, then this will be estimated from the median of the noise source phase over time.

  • weight_static_phi (np.ndarray[nfreq, nsource]) – Inverse variance of the time averaged phased. Set to zero for frequencies and inputs that are missing or should be ignored. If None, then this will be estimated from the residuals of the fit.

  • static_phi_fit (np.ndarray[nparam, nsource]) – Polynomial fit to static_phi versus frequency.

  • static_amp (np.ndarray[nfreq, nsource]) – Subtract this quantity from the noise source amplitude prior to fitting for the amplitude variations. If None, then this will be estimated from the median of the noise source amplitude over time.

  • weight_static_amp (np.ndarray[nfreq, nsource]) – Inverse variance of the time averaged amplitude. Set to zero for frequencies and inputs that are missing or should be ignored. If None, then this will be estimated from the residuals of the fit.

Returns:

  • phi (np.ndarray[nfreq, nsource, ntime]) – Phase of the signal from the noise source.

  • weight_phi (np.ndarray[nfreq, nsource, ntime]) – Inverse variance of the phase of the signal from the noise source.

  • tau (np.ndarray[nsource, ntime]) – Delay template for each noise source input.

  • weight_tau (np.ndarray[nfreq, nsource]) – Estimate of the uncertainty on the delay template (inverse variance).

  • static_phi (np.ndarray[nfreq, nsource]) – Time averaged phase versus frequency.

  • weight_static_phi (np.ndarray[nfreq, nsource]) – Inverse variance of the time averaged phase.

  • static_phi_fit (np.ndarray[nparam, nsource]) – Best-fit parameters of the polynomial fit to the time averaged phase versus frequency.

  • amp (np.ndarray[nfreq, nsource, ntime]) – Amplitude of the signal from the noise source.

  • weight_amp (np.ndarray[nfreq, nsource, ntime]) – Inverse variance of the amplitude of the signal from the noise source.

  • alpha (np.ndarray[nsource, ntime]) – Amplitude coefficient for each noise source input.

  • weight_alpha (np.ndarray[nfreq, nsource]) – Estimate of the uncertainty on the amplitude coefficient (inverse variance).

  • static_amp (np.ndarray[nfreq, nsource]) – Time averaged amplitude versus frequency.

  • weight_static_amp (np.ndarray[nfreq, nsource]) – Inverse variance of the time averaged amplitude.

  • num_freq (np.ndarray[nsource, ntime]) – Number of frequencies used to construct the delay and amplitude templates.

ch_util.timing.eigen_decomposition(vis, flag)[source]

Eigenvalue decomposition of the visibility matrix.

Parameters:
  • vis (np.ndarray[nfreq, nprod, ntime]) – Upper-triangle, product packed visibility matrix.

  • flag (np.ndarray[nfreq, nsource, ntime] (optional)) – Array of 1 or 0 indicating the inputs that should be included in the eigenvalue decomposition for each frequency and time.

Returns:

resp – Eigenvector corresponding to the largest eigenvalue for each frequency and time.

Return type:

np.ndarray[nfreq, nsource, ntime]

ch_util.timing.fit_poly_to_phase(freq, resp, resp_error, nparam=2)[source]

Fit complex data versus frequency to a model consisting of a polynomial in phase.

Nonlinear least squares algorithm is applied to the complex data to avoid problems caused by phase wrapping.

Parameters:
  • freq (np.ndarray[nfreq, ]) – Frequency in MHz.

  • resp (np.ndarray[nfreq, ]) – Complex data with magnitude equal to 1.0.

  • resp_error (np.ndarray[nfreq, ]) – Uncertainty on the complex data.

  • nparam (int) – Number of parameters in the polynomial. Default is 2 (i.e, linear).

Returns:

  • popt (np.ndarray[nparam, ]) – Best-fit parameters.

  • pcov (np.ndarray[nparam, nparam]) – Covariance of the best-fit parameters. Assumes that it obtained a good fit and returns the errors necessary to achieve that.

ch_util.timing.load_timing_correction(files, start=None, stop=None, window=43200.0, instrument='chime', **kwargs)[source]

Find and load the timing correction for a list of corr acquisition files.

For example, if the instrument keyword is set to ‘chime’, then this function will accept all types of chime corr acquisition files, such as ‘chimetiming’, ‘chimepb’, ‘chimeN2’, ‘chimecal’, and then find the relevant set of ‘chimetiming’ files to load.

Accepts and passes on all keyword arguments for the functions andata.CorrData.from_acq_h5 and construct_delay_template.

Should consider modifying this method to use Finder at some point in future.

Parameters:
  • files (string or list of strings) – Absolute path to corr acquisition file(s).

  • start (integer, optional) – What frame to start at in the full set of files.

  • stop (integer, optional) – What frame to stop at in the full set of files.

  • window (float) – Use the timing data -window from start and +window from stop. Default is 12 hours.

  • instrument (string) – Name of the instrument. Default is ‘chime’.

Returns:

data

Return type:

TimingData

ch_util.timing.map_input_to_noise_source(inputs, noise_sources)[source]

Find the appropriate noise source to use to correct the phase of each input.

Searches for a noise source connected to the same slot, then crate, then hut, then correlator.

Parameters:
  • inputs (np.ndarray[ninput, ] of dtype=('chan_id', 'correlator_input')) – The input axis from a data acquisition file.

  • noise_sources (np.ndarray[nsource, ] of dtype=('chan_id', 'correlator_input')) – The noise sources.

ch_util.timing.model_poly_phase(freq, *param)[source]

Evaluate a polynomial model for the phase.

To be used with the parameters output from fit_poly_to_phase.

Parameters:
  • freq (np.ndarray[nfreq, ]) – Frequency in MHz.

  • *param (float) – Coefficients of the polynomial.

Returns:

phi – Phase in radians between -pi and +pi.

Return type:

np.ndarray[nfreq, ]