mmwave.dsp.angle_estimation

mmwave.dsp.angle_estimation.azimuth_processing(radar_cube, det_obj_2d, config, window_type_2d=None)

Calculate the X/Y coordinates for all detected objects.

The following procedures will be performed in this function:

  1. Filter radarCube based on the range indices from detObj2D and optional clutter removal.
  2. Re-do windowing and 2D FFT, select associated doppler indices to form the azimuth input.
  3. Doppler compensation on the virtual antennas related to tx2. Save optional copy for near field compensation and velocity disambiguation.
  4. Perform azimuth FFT.
  5. Optional near field correction and velocity disambiguation. Currently mutual exclusive.
  6. Magnitude squared.
  7. Calculate X/Y coordinates.
Parameters:
  • radar_cube – (numChirpsPerFrame, numRxAntennas, numRangeBins). Because of the space limitation, TI Demo starts again from the 1D FFT, recalculate the 2D FFT on the selected range bins and pick the associated doppler bins for the azimuth FFT.
  • det_obj_2d – (numDetObj, 3)
  • config – [TBD]
  • window_type_2d – windowing function for the 2D FFT. Default is None.
Returns:

(numDetObj, 5). Copy of detObj2D but with populated X/Y coordinates.

Return type:

azimuthOut

mmwave.dsp.angle_estimation.aoa_bartlett(steering_vec, sig_in, axis)

Perform AOA estimation using Bartlett Beamforming on a given input signal (sig_in). Make sure to specify the correct axis in (axis) to ensure correct matrix multiplication. The power spectrum is calculated using the following equation:

\[P_{ca} (\theta) = a^{H}(\theta) R_{xx}^{-1} a(\theta)\]

This steers the beam using the steering vector as weights:

\[w_{ca} (\theta) = a(\theta)\]
Parameters:
  • steering_vec (ndarray) – A 2D-array of size (numTheta, num_ant) generated from gen_steering_vec
  • sig_in (ndarray) – Either a 2D-array or 3D-array of size (num_ant, numChirps) or (numChirps, num_vrx, num_adc_samples) respectively, containing ADC sample data sliced as described
  • axis (int) – Specifies the axis where the Vrx data in contained.
Returns:

A 3D-array of size (numChirps, numThetas, numSamples)

Return type:

doa_spectrum (ndarray)

Example

>>> # In this example, dataIn is the input data organized as numFrames by RDC
>>> frame = 0
>>> dataIn = np.random.rand((num_frames, num_chirps, num_vrx, num_adc_samples))
>>> aoa_bartlett(steering_vec,dataIn[frame],axis=1)
mmwave.dsp.angle_estimation.aoa_capon(x, steering_vector, magnitude=False)

Perform AOA estimation using Capon (MVDR) Beamforming on a rx by chirp slice

Calculate the aoa spectrum via capon beamforming method using one full frame as input. This should be performed for each range bin to achieve AOA estimation for a full frame This function will calculate both the angle spectrum and corresponding Capon weights using the equations prescribed below.

\[ \begin{align}\begin{aligned}P_{ca} (\theta) = \frac{1}{a^{H}(\theta) R_{xx}^{-1} a(\theta)}\\w_{ca} (\theta) = \frac{R_{xx}^{-1} a(\theta)}{a^{H}(\theta) R_{xx}^{-1} a(\theta)}\end{aligned}\end{align} \]
Parameters:
  • x (ndarray) – Output of the 1d range fft with shape (num_ant, numChirps)
  • steering_vector (ndarray) – A 2D-array of size (numTheta, num_ant) generated from gen_steering_vec
  • magnitude (bool) – Azimuth theta bins should return complex data (False) or magnitude data (True). Default=False
Raises:

ValueError – steering_vector and or x are not the correct shape

Returns:

A list containing numVec and steeringVectors den (ndarray: A 1D-Array of size (numTheta) containing azimuth angle estimations for the given range weights (ndarray): A 1D-Array of size (num_ant) containing the Capon weights for the given input data

Example

>>> # In this example, dataIn is the input data organized as numFrames by RDC
>>> Frame = 0
>>> dataIn = np.random.rand((num_frames, num_chirps, num_vrx, num_adc_samples))
>>> for i in range(256):
>>>     scan_aoa_capon[i,:], _ = dss.aoa_capon(dataIn[Frame,:,:,i].T, steering_vector, magnitude=True)
mmwave.dsp.angle_estimation.cov_matrix(x)
Calculates the spatial covariance matrix (Rxx) for a given set of input data (x=inputData).
Assumes rows denote Vrx axis.
Parameters:x (ndarray) – A 2D-Array with shape (rx, adc_samples) slice of the output of the 1D range fft
Returns:A 2D-Array with shape (rx, rx)
Return type:Rxx (ndarray)
mmwave.dsp.angle_estimation.forward_backward_avg(Rxx)

Performs forward backward averaging on the given input square matrix

Parameters:Rxx (ndarray) – A 2D-Array square matrix containing the covariance matrix for the given input data
Returns:The 2D-Array square matrix containing the forward backward averaged covariance matrix
Return type:R_fb (ndarray)
Wrapper function to perform scipy.signal’s prescribed peak search algorithm
Tested Runtime: 45 µs ± 2.61 µs per loop (mean ± std. dev. of 7 runs, 10000 loops each)
Parameters:
  • doa_spectrum (ndarray) – A 1D-Array of size (numTheta, 1) containing the theta spectrum at a given range bin
  • peak_threshold_weight (float) – A float specifying the desired peak declaration threshold weight to be applied
Returns:

The number of max points found by the algorithm peaks (list): List of indexes where peaks are located total_power (float): Total power in the current spectrum slice

Return type:

num_max (int)

mmwave.dsp.angle_estimation.peak_search_full(doa_spectrum, gamma=1.2, peak_threshold_weight=0.251188643150958)

Perform TI prescribed peak search algorithm Tested Runtime: 147 µs ± 4.27 µs per loop (mean ± std. dev. of 7 runs, 10000 loops each)

Parameters:
  • doa_spectrum (ndarray) – A 1D-Array of size (numTheta, 1) containing the theta spectrum at a given range bin
  • gamma (float) – A float specifying the maximum/minimum wiggle necessary to qualify as a peak
  • peak_threshold_weight (float) – A float specifying the desired peak declaration threshold weight to be applied
Returns:

The number of max points found by the algorithm ang_est (list): List of indexes where the peaks are located

Return type:

num_max (int)

mmwave.dsp.angle_estimation.peak_search_full_variance(doa_spectrum, steering_vec_size, sidelobe_level=0.251188643150958, gamma=1.2)

Performs peak search (TI’s full search) will retaining details about each peak including each peak’s width, location, and value.

Parameters:
  • doa_spectrum (ndarray) – a 1D numpy array containing the power spectrum generated via some aoa method (naive,
  • or capon) (bartlett,) –
  • steering_vec_size (int) – Size of the steering vector in terms of number of theta bins
  • sidelobe_level (float) – A low value threshold used to avoid sidelobe detections as peaks
  • gamma (float) – Weight to determine when a peak will pass as a true peak
Returns:

A 1D numpy array of custom data types with length numberOfPeaksDetected. Each detected peak is organized as [peak_location, peak_value, peak_width] total_power (float): The total power of the spectrum. Used for variance calculations

Return type:

peak_data (ndarray)

mmwave.dsp.angle_estimation.variance_estimation(num_max, est_resolution, peak_data, total_power, width_adjust_3d_b=2.5, input_snr=10000)
This function will calculate an estimated variance value for each detected peak. This should
be run after running peak_search_full_variance
Parameters:
  • num_max (int) – The number of detected peaks
  • est_resolution (float) – The desired resolution in terms of theta
  • peak_data (ndarray) – A numpy array of dictionaries, where each dictionary is of the form: {“peakLoc”: , “peakVal”: , “peakWid”: }
  • total_power (float) – The total power of the spectrum
  • width_adjust_3d_b (float) – Constant to adjust the gamma bandwidth to 3dB level
  • input_snr (int) – the linear snr for the input signal samples
Returns:

A 1D array of variances (of the peaks). The order of the peaks is preserved from peak_data

Return type:

est_var (ndarray)

mmwave.dsp.angle_estimation.gen_steering_vec(ang_est_range, ang_est_resolution, num_ant)

Generate a steering vector for AOA estimation given the theta range, theta resolution, and number of antennas

Defines a method for generating steering vector data input –Python optimized Matrix format The generated steering vector will span from -angEstRange to angEstRange with increments of ang_est_resolution The generated steering vector should be used for all further AOA estimations (bartlett/capon)

Parameters:
  • ang_est_range (int) – The desired span of thetas for the angle spectrum.
  • ang_est_resolution (float) – The desired resolution in terms of theta
  • num_ant (int) – The number of Vrx antenna signals captured in the RDC
Returns:

Number of vectors generated (integer divide angEstRange/ang_est_resolution) steering_vectors (ndarray): The generated 2D-array steering vector of size (num_vec,num_ant)

Return type:

num_vec (int)

Example

>>> #This will generate a numpy array containing the steering vector with
>>> #angular span from -90 to 90 in increments of 1 degree for a 4 Vrx platform
>>> _, steering_vec = gen_steering_vec(90,1,4)
mmwave.dsp.angle_estimation.aoa_estimation_bf_one_point(num_ant, sig_in, steering_vec)

Calculates the total power of the given spectrum

Parameters:
  • num_ant (int) – The number of virtual antennas (Vrx) being used
  • sig_in (ndarray) – A 2D-array of size (num_ant, numChirps) containing ADC sample data sliced for used Vrx
  • steering_vec (ndarray) – A 2D-array of size (numTheta, num_ant) generated from gen_steering_vec
Returns:

The total power of the given input spectrum

Return type:

out_value (complex)

mmwave.dsp.angle_estimation.aoa_est_bf_single_peak_det(sig_in, steering_vec)
Beamforming Estimate Angle of Arrival for single peak (single peak should be known a priori)
Function call does not include variance calculations Function does not generate a spectrum. Rather, it only returns the array index (theta) to the highest peak
Parameters:
  • sig_in (ndarray) – A 2D-array of size (num_ant, numChirps) containing ADC sample data sliced as described
  • steering_vec (ndarray) – A generated 2D-array steering vector of size (numVec,num_ant)
Returns:

Index of the theta spectrum at a given range bin that contains the max peak

Return type:

max_index (int)

mmwave.dsp.angle_estimation.aoa_est_bf_single_peak(num_ant, noise, est_resolution, sig_in, steering_vec_size, steering_vec)
Beamforming Estimate Angle of Arrival for single peak (single peak should be known a priori)
Function call includes variance calculations Function does generate a spectrum.
Parameters:
  • num_ant (int) – The number of virtual receivers in the current radar setup
  • noise (float) – Input noise figure
  • est_resolution (float) – Desired theta spectrum resolution used when generating steering_vec
  • sig_in (ndarray) – A 2D-array of size (num_ant, numChirps) containing ADC sample data sliced as described
  • steering_vec_size (int) – Length of the steering vector array
  • steering_vec (ndarray) – A generated 2D-array steering vector of size (numVec,num_ant)
Returns:

The estimated variance of the doa_spectrum max_index (int): Index of the theta spectrum at a given range bin that contains the max peak doa_spectrum (ndarray): A 1D-Array of size (numTheta, 1) containing the theta spectrum at a given range bin

Return type:

est_var (float)

mmwave.dsp.angle_estimation.aoa_est_bf_multi_peak_det(gamma, sidelobe_level, sig_in, steering_vec, steering_vec_size, ang_est, search=False)

Use Bartlett beamforming to estimate AOA for multi peak situation (a priori), no variance calculation

Parameters:
  • gamma (float) – Weight to determine when a peak will pass as a true peak
  • sidelobe_level (float) – A low value threshold used to avoid sidelobe detections as peaks
  • sig_in (ndarray) – A 2D-array of size (num_ant, numChirps) containing ADC sample data sliced as described
  • steering_vec (ndarray) – A generated 2D-array steering vector of size (numVec,num_ant)
  • steering_vec_size (int) – Length of the steering vector array
  • ang_est (ndarray) – An empty 1D numpy array that gets populated with max indexes
  • search (bool) – Flag that determines whether search is done to find max points
Returns:

The number of max points found across the theta bins at this particular range bin doa_spectrum (ndarray): A 1D-Array of size (numTheta, 1) containing the theta spectrum at a given range bin

Return type:

num_max (int)

mmwave.dsp.angle_estimation.aoa_est_bf_multi_peak(gamma, sidelobe_level, width_adjust_3d_b, input_snr, est_resolution, sig_in, steering_vec, steering_vec_size, peak_data, ang_est)

This function performs all sections of the angle of arrival process in one function.

  1. Performs bartlett beamforming
  2. Performs multi-peak search
  3. Calculates an estimated variance
Parameters:
  • gamma (float) – Weight to determine when a peak will pass as a true peak
  • sidelobe_level (float) – A low value threshold used to avoid sidelobe detections as peaks
  • width_adjust_3d_b (float) – Constant to adjust gamma bandwidth to 3dB bandwidth
  • input_snr (float) – Input data SNR value
  • est_resolution (float) – User defined target resolution
  • sig_in (ndarray) – A 2D-array of size (num_ant, numChirps) containing ADC sample data sliced as described
  • steering_vec (ndarray) – A generated 2D-array steering vector of size (numVec,num_ant)
  • steering_vec_size (int) – Length of the steering vector array
  • peak_data (ndarray) – A 2D ndarray with custom data-type that contains information on each detected point
  • ang_est (ndarray) – An empty 1D numpy array that gets populated with max indexes
Returns:

Tuple [ndarray, ndarray]
  1. num_max (int): The number of max values detected by search algorithm
  2. est_var (ndarray): The estimated variance of this range of thetas at this range bin

mmwave.dsp.angle_estimation.naive_xyz(virtual_ant, num_tx=3, num_rx=4, fft_size=64)

Estimate the phase introduced from the elevation of the elevation antennas

Parameters:
  • virtual_ant – Signal received by the rx antennas, shape = [#angleBins, #detectedObjs], zero-pad #virtualAnts to #angleBins
  • num_tx – Number of transmitter antennas used
  • num_rx – Number of receiver antennas used
  • fft_size – Size of the fft performed on the signals
Returns:

Estimated x axis coordinate in meters (m) y_vector (float): Estimated y axis coordinate in meters (m) z_vector (float): Estimated z axis coordinate in meters (m)

Return type:

x_vector (float)

mmwave.dsp.angle_estimation.beamforming_naive_mixed_xyz(azimuth_input, input_ranges, range_resolution, method='Capon', num_vrx=12, est_range=90, est_resolution=1)

This function estimates the XYZ location of a series of input detections by performing beamforming on the azimuth axis and naive AOA on the vertical axis.

TI xWR1843 virtual antenna map Row 1 8 9 10 11 Row 2 0 1 2 3 4 5 6 7

phi (ndarray): theta (ndarray): ranges (ndarray): xyz_vec (ndarray):

Parameters:
  • azimuth_input (ndarray) – Must be a numpy array of shape (numDetections, numVrx)
  • input_ranges (ndarray) – Numpy array containing the rangeBins that have detections (will determine x, y, z for
  • detection) (each) –
  • range_resolution (float) – The range_resolution in meters per rangeBin for rangeBin->meter conversion
  • method (string) – Determines which beamforming method to use for azimuth aoa estimation.
  • num_vrx (int) – Number of virtual antennas in the radar platform. Default set to 12 for 1843
  • est_range (int) – The desired span of thetas for the angle spectrum. Used for gen_steering_vec
  • est_resolution (float) – The desired angular resolution for gen_steering_vec
Raises:
  • ValueError – If method is not one of two AOA implementations (‘Capon’, ‘Bartlett’)
  • ValueError – azimuthInput’s second axis should have same shape as the number of Vrx
Returns:

  1. A numpy array of shape (numDetections, ) where each element represents the elevation angle in degrees
  2. A numpy array of shape (numDetections, ) where each element represents the azimuth in degrees
  3. A numpy array of shape (numDetections, ) where each element represents the polar range in rangeBins
  4. A numpy array of shape (3, numDetections) and format: [x, y, z] where x, y, z are 1D arrays. x, y, z should be in meters

Return type:

tuple [ndarray, ndarray, ndarray, ndarray, list]