Multiple-Input Multiple-Output (MIMO)

Stream Management

Stream management determines which transmitter is sending which stream to which receiver. Transmitters and receivers can be user terminals or base stations, depending on whether uplink or downlink transmissions are considered. The StreamManagement class has various properties that are needed to recover desired or interfering channel coefficients for precoding and equalization. In order to understand how the various properties of StreamManagement can be used, we recommend to have a look at the source code of the LMMSEEqualizer or RZFPrecoder.

The following code snippet shows how to configure StreamManagement for a simple uplink scenario, where four transmitters send each one stream to a receiver. Note that StreamManagement is independent of the actual number of antennas at the transmitters and receivers.

num_tx = 4
num_rx = 1
num_streams_per_tx = 1

# Indicate which transmitter is associated with which receiver
# rx_tx_association[i,j] = 1 means that transmitter j sends one
# or mutiple streams to receiver i.
rx_tx_association = np.zeros([num_rx, num_tx])
rx_tx_association[0,0] = 1
rx_tx_association[0,1] = 1
rx_tx_association[0,2] = 1
rx_tx_association[0,3] = 1

sm = StreamManagement(rx_tx_association, num_streams_per_tx)
class sionna.phy.mimo.StreamManagement(rx_tx_association, num_streams_per_tx)[source]

Class for management of streams in multi-cell MIMO networks.

Parameters:
  • rx_tx_association ([num_rx, num_tx], np.int) – A binary NumPy array where rx_tx_association[i,j]=1 means that receiver i gets one or multiple streams from transmitter j.

  • num_streams_per_tx (int) – Indicates the number of streams that are transmitted by each transmitter.

Note

Several symmetry constraints on rx_tx_association are imposed to ensure efficient processing. All row sums and all column sums must be equal, i.e., all receivers have the same number of associated transmitters and all transmitters have the same number of associated receivers. It is also assumed that all transmitters send the same number of streams num_streams_per_tx.

property detection_desired_ind

Indices needed to gather desired channels for receive processing

A NumPy array of shape [num_rx*num_streams_per_rx] that can be used to gather desired channels from the flattened channel tensor of shape […,num_rx, num_tx, num_streams_per_tx,…]. The result of the gather operation can be reshaped to […,num_rx, num_streams_per_rx,…].

property detection_undesired_ind

Indices needed to gather undesired channels for receive processing

A NumPy array of shape [num_rx*num_streams_per_rx] that can be used to gather undesired channels from the flattened channel tensor of shape […,num_rx, num_tx, num_streams_per_tx,…]. The result of the gather operation can be reshaped to […,num_rx, num_interfering_streams_per_rx,…].

property num_interfering_streams_per_rx

Number of interfering streams received at each eceiver

property num_rx

Number of receivers

property num_rx_per_tx

Number of receivers communicating with a transmitter

property num_streams_per_rx

Number of streams transmitted to each receiver

property num_streams_per_tx

Number of streams per transmitter

property num_tx

Number of transmitters

property num_tx_per_rx

Number of transmitters communicating with a receiver

property precoding_ind

Indices needed to gather channels for precoding

A NumPy array of shape [num_tx, num_rx_per_tx], where precoding_ind[i,:] contains the indices of the receivers to which transmitter i is sending streams.

property rx_stream_ids

Mapping of streams to receivers

A Numpy array of shape [num_rx, num_streams_per_rx]. This array is obtained from tx_stream_ids together with the rx_tx_association. rx_stream_ids[i,:] contains the indices of streams that are supposed to be decoded by receiver i.

property rx_tx_association

Association between receivers and transmitters.

A binary NumPy array of shape [num_rx, num_tx], where rx_tx_association[i,j]=1 means that receiver i gets one ore multiple streams from transmitter j.

property stream_association

Association between receivers, transmitters, and streams

A binary NumPy array of shape [num_rx, num_tx, num_streams_per_tx], where stream_association[i,j,k]=1 means that receiver i gets the k th stream from transmitter j.

property stream_ind

Indices needed to gather received streams in the correct order

A NumPy array of shape [num_rx*num_streams_per_rx] that can be used to gather streams from the flattened tensor of received streams of shape […,num_rx, num_streams_per_rx,…]. The result of the gather operation is then reshaped to […,num_tx, num_streams_per_tx,…].

property tx_stream_ids

Mapping of streams to transmitters

A NumPy array of shape [num_tx, num_streams_per_tx]. Streams are numbered from 0,1,… and assiged to transmitters in increasing order, i.e., transmitter 0 gets the first num_streams_per_tx and so on.

Precoding

sionna.phy.mimo.cbf_precoding_matrix(h, precision=None)[source]

Computes the conjugate beamforming (CBF) Precoder

This function computes the CBF precoding matrix for a MIMO link, assuming the following model:

y=HGx+n

where yCK is the received signal vector, HCK×M is the known channel matrix, GCM×K is the precoding matrix, xCK is the symbol vector to be precoded, and nCK is a noise vector.

The precoding matrix G is defined as :

G=VD

where

V=HHD=diag(vk21,k=0,,K1).

The matrix V ensures that each stream is precoded with a unit-norm vector, i.e., tr(GGH)=K. The function returns the matrix G.

Input:
  • h ([…,K,M], tf.complex) – 2+D tensor containing the channel matrices

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used

Output:

g ([…,M,K], tf.complex) – 2+D tensor containing the precoding matrices

sionna.phy.mimo.rzf_precoding_matrix(h, alpha=0.0, precision=None)[source]

Computes the Regularized Zero-Forcing (RZF) Precoder

This function computes the RZF precoding matrix for a MIMO link, assuming the following model:

y=HGx+n

where yCK is the received signal vector, HCK×M is the known channel matrix, GCM×K is the precoding matrix, xCK is the symbol vector to be precoded, and nCK is a noise vector.

The precoding matrix G is defined as :

G=VD

where

V=HH(HHH+αI)1D=diag(vk21,k=0,,K1)

where α>0 is the regularization parameter. The matrix V ensures that each stream is precoded with a unit-norm vector, i.e., tr(GGH)=K. The function returns the matrix G.

Input:
  • h ([…,K,M], tf.complex) – 2+D tensor containing the channel matrices

  • alpha (0. (default) | […], tf.float32) – Regularization parameter

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used

Output:

g ([…,M,K], tf.complex) – 2+D tensor containing the precoding matrices

sionna.phy.mimo.rzf_precoder(x, h, alpha=0.0, return_precoding_matrix=False, precision=None)[source]

Regularized Zero-Forcing (RZF) Precoder

This function implements RZF precoding for a MIMO link, assuming the following model:

y=HGx+n

where yCK is the received signal vector, HCK×M is the known channel matrix, GCM×K is the precoding matrix, xCK is the symbol vector to be precoded, and nCK is a noise vector.

The precoding matrix G is defined as (Eq. 4.37) [BHS2017] :

G=VD

where

V=HH(HHH+αI)1D=diag(vk21,k=0,,K1)

where α>0 is the regularization parameter.

This ensures that each stream is precoded with a unit-norm vector, i.e., tr(GGH)=K. The function returns the precoded vector Gx.

Input:
  • x ([…,K], tf.complex) – Symbol vectors to be precoded

  • h ([…,K,M], tf.complex) – Channel matrices

  • alpha (0. (default) | […], tf.float32) – Regularization parameter

  • return_precoding_matrices (bool, (default, False)) – Indicates if the precoding matrices should be returned or not

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used

Output:
  • x_precoded ([…,M], tf.complex) – Precoded symbol vectors

  • g ([…,M,K], tf.complex) – Precoding matrices. Only returned if return_precoding_matrices=True.

sionna.phy.mimo.grid_of_beams_dft_ula(num_ant, oversmpl=1, precision=None)[source]

Computes the Discrete Fourier Transform (DFT) Grid of Beam (GoB) coefficients for a uniform linear array (ULA)

The coefficient applied to antenna n for beam m is expressed as:

cnm=e2πnmNO,n=0,,N1, m=0,,NO

where N is the number of antennas num_ant and O is the oversampling factor oversmpl.

Note that the main lobe of beam m points in the azimuth direction θ=arcsin(2mN) if mN/2 and θ=arcsin(2mNN) if mN/2, where θ=0 defines the perpendicular to the antenna array.

Input:
  • num_ant (int) – Number of antennas

  • oversmpl (int, (default 1)) – Oversampling factor

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:

gob ([num_ant x oversmpl, num_ant], tf.complex) – The m-th row contains the num_ant antenna coefficients for the m-th DFT beam.

sionna.phy.mimo.grid_of_beams_dft(num_ant_v, num_ant_h, oversmpl_v=1, oversmpl_h=1, precision=None)[source]

Computes the Discrete Fourier Transform (DFT) Grid of Beam (GoB) coefficients for a uniform rectangular array (URA)

GoB indices are arranged over a 2D grid indexed by (mv,mh). The coefficient of the beam with index (mv,mh) applied to the antenna located at row nv and column nh of the rectangular array is expressed as:

cnv,nhmv,mh=e2πnhmvNhOhe2πnhmhNvOv

where nv=0,,Nv1, nh=0,,Nh1, mv=0,,NvOv, mh=0,,NhOh, N is the number of antennas num_ant and Ov,Oh are the oversampling factor oversmpl_v, oversmpl_h in the vertical and horizontal direction, respectively.

We can rewrite more concisely the matrix coefficients cmv,mh as follows:

cmv,mh=cmvcmh

where denotes the Kronecker product and cmv,cmh are the ULA DFT beams computed as in grid_of_beams_dft_ula().

Such a DFT GoB is, e.g., defined in Section 5.2.2.2.1 [3GPP38214].

Input:
  • num_ant_v (int) – Number of antenna rows (i.e., in vertical direction) of the rectangular array

  • num_ant_h (int) – Number of antenna columns (i.e., in horizontal direction) of the rectangular array.

  • oversmpl_v (int, (default 1)) – Oversampling factor in vertical direction

  • oversmpl_h (int, (default 1)) – Oversampling factor in horizontal direction

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:

gob ([num_ant_v x oversmpl_v, num_ant_h x oversmpl_h, num_ant_v x num_ant_h], tf.complex) – The elements [mv,mh,:] contain the antenna coefficients of the DFT beam with index pair (mv,mh).

sionna.phy.mimo.flatten_precoding_mat(precoding_mat, by_column=True)[source]

Flattens a […, num_ant_v, num_ant_h] precoding matrix associated with a rectangular array by producing a […, num_ant_v x num_ant_h] precoding vector

Input:
  • precoding_mat ([…, num_antennas_vertical, num_antennas_horizontal], tf.complex) – Precoding matrix. The element (i,j) contains the precoding coefficient of the antenna element located at row i and column j of a rectangular antenna array.

  • by_column (bool, (default True)) – If True, then flattening occurs on a per-column basis, i.e., the first column is appended to the second, and so on. Else, flattening is performed on a per-row basis.

Output:

[…, num_antennas_vertical x num_antennas_horizontal], tf.complex – Flattened precoding matrix

sionna.phy.mimo.normalize_precoding_power(precoding_vec, tx_power_list=None, precision=None)[source]

Normalizes the beam coefficient power to 1 by default, or to tx_power_list if provided as input

Input:
  • precoding_vec ([N,M], tf.complex) – Each row contains a set of antenna coefficients whose power is to be normalized.

  • tx_power_list ([N], float) – The i-th element defines the power of the i-th precoding vector.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:

[N,M] tf.complex – Normalized antenna coefficients

Equalization

sionna.phy.mimo.lmmse_matrix(h, s=None, precision=None)[source]

MIMO LMMSE Equalization matrix

This function computes the LMMSE equalization matrix for a MIMO link, assuming the following model:

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector. It is assumed that E[x]=E[n]=0, E[xxH]=IK and E[nnH]=S.

This function returns the LLMSE equalization matrix:

G=HH(HHH+S)1.

If S=IM, a numerically more stable version of the equalization matrix is computed:

G=(HHH+I)1HH.
Input:
  • h ([…,M,K], tf.complex) – Channel matrices

  • s (None (default) | […,M,M], tf.complex) – Noise covariance matrices. If None, the noise is assumed to be white with unit variance.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:

g ([…,K,M], tf.complex) – LLMSE equalization matrices

sionna.phy.mimo.lmmse_equalizer(y, h, s, whiten_interference=True, precision=None)[source]

MIMO LMMSE Equalizer

This function implements LMMSE equalization for a MIMO link, assuming the following model:

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector. It is assumed that E[x]=E[n]=0, E[xxH]=IK and E[nnH]=S.

The estimated symbol vector x^CK is given as (Lemma B.19) [BHS2017] :

x^=diag(GH)1Gy

where

G=HH(HHH+S)1.

This leads to the post-equalized per-symbol model:

x^k=xk+ek,k=0,,K1

where the variances σk2 of the effective residual noise terms ek are given by the diagonal elements of

diag(E[eeH])=diag(GH)1I.

Note that the scaling by diag(GH)1 is important for the Demapper although it does not change the signal-to-noise ratio.

The function returns x^ and σ2=[σ02,,σK12]T.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,K], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

  • whiten_interference (bool, (default True)) – If True, the interference is first whitened before equalization. In this case, an alternative expression for the receive filter is used that can be numerically more stable.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:
  • x_hat ([…,K], tf.complex) – Estimated symbol vectors

  • no_eff (tf.float) – Effective noise variance estimates

sionna.phy.mimo.mf_equalizer(y, h, s, precision=None)[source]

MIMO Matched Filter (MF) Equalizer

This function implements matched filter (MF) equalization for a MIMO link, assuming the following model:

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector. It is assumed that E[x]=E[n]=0, E[xxH]=IK and E[nnH]=S.

The estimated symbol vector x^CK is given as (Eq. 4.11) [BHS2017] :

x^=Gy

where

G=diag(HHH)1HH.

This leads to the post-equalized per-symbol model:

x^k=xk+ek,k=0,,K1

where the variances σk2 of the effective residual noise terms ek are given by the diagonal elements of the matrix

E[eeH]=(IGH)(IGH)H+GSGH.

Note that the scaling by diag(HHH)1 in the definition of G is important for the Demapper although it does not change the signal-to-noise ratio.

The function returns x^ and σ2=[σ02,,σK12]T.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,K], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:
  • x_hat ([…,K], tf.complex) – Estimated symbol vectors

  • no_eff (tf.float) – Effective noise variance estimates

sionna.phy.mimo.zf_equalizer(y, h, s, precision=None)[source]

Applies MIMO ZF Equalizer

This function implements zero-forcing (ZF) equalization for a MIMO link, assuming the following model:

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector. It is assumed that E[x]=E[n]=0, E[xxH]=IK and E[nnH]=S.

The estimated symbol vector x^CK is given as (Eq. 4.10) [BHS2017] :

x^=Gy

where

G=(HHH)1HH.

This leads to the post-equalized per-symbol model:

x^k=xk+ek,k=0,,K1

where the variances σk2 of the effective residual noise terms ek are given by the diagonal elements of the matrix

E[eeH]=GSGH.

The function returns x^ and σ2=[σ02,,σK12]T.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,K], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Output:
  • x_hat ([…,K], tf.complex) – Estimated symbol vectors

  • no_eff (tf.float) – Effective noise variance estimates

Detection

class sionna.phy.mimo.EPDetector(output, num_bits_per_symbol, hard_out=False, l=10, beta=0.9, precision=None, **kwargs)[source]

MIMO Expectation Propagation (EP) detector

This block implements Expectation Propagation (EP) MIMO detection as described in [EP2014]. It can generate hard- or soft-decisions for symbols or bits.

This block assumes the following channel model:

y=Hx+n

where yCM is the received signal vector, xCS is the vector of transmitted symbols which are uniformly and independently drawn from the constellation C, HCM×S is the known channel matrix, and nCM is a complex Gaussian noise vector. It is assumed that E[n]=0 and E[nnH]=S, where S has full rank.

The channel model is first whitened using whiten_channel() and then converted to its real-valued equivalent, see complex2real_channel(), prior to MIMO detection.

The computation of LLRs is done by converting the symbol logits that naturally arise in the algorithm to LLRs using PAM2QAM(). Custom conversions of symbol logits to LLRs can be implemented by using the soft-symbol output.

The detector is currently restricted to QAM constellations.

Parameters:
  • output ("bit" | "symbol") – Type of output, either LLRs on bits or logits on constellation symbols

  • num_bits_per_symbol (int) – Number of bits per QAM constellation symbol, e.g., 4 for QAM16.

  • hard_out (bool, (default False)) – If True, the detector computes hard-decided bit values or constellation point indices instead of soft-values.

  • l (int, (default 10)) – Number of iterations

  • beta (float, (default 0.9)) – Parameter β[0,1] for update smoothing

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,num_streams], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

Output:
  • One of

  • […, num_streams, num_bits_per_symbol], tf.float – LLRs or hard-decisions for every bit of every stream, if output equals “bit”.

  • […, num_streams, num_points], tf.float or […, num_streams], tf.int32 – Logits or hard-decisions for constellation symbols for every stream, if output equals “symbol”. Hard-decisions correspond to the symbol indices.

class sionna.phy.mimo.KBestDetector(output, num_streams, k, constellation_type=None, num_bits_per_symbol=None, constellation=None, hard_out=False, use_real_rep=False, list2llr=None, precision=None, **kwargs)[source]

MIMO K-Best detector

This block implements K-Best MIMO detection as described in (Eq. 4-5) [FT2015]. It can either generate hard decisions (for symbols or bits) or compute LLRs.

The algorithm operates in either the complex or real-valued domain. Although both options produce identical results, the former has the advantage that it can be applied to arbitrary non-QAM constellations. It also reduces the number of streams (or depth) by a factor of two.

The way soft-outputs (i.e., LLRs) are computed is determined by the list2llr function. The default solution List2LLRSimple assigns a predetermined value to all LLRs without counter-hypothesis.

This block assumes the following channel model:

y=Hx+n

where yCM is the received signal vector, xCS is the vector of transmitted symbols which are uniformly and independently drawn from the constellation C, HCM×S is the known channel matrix, and nCM is a complex Gaussian noise vector. It is assumed that E[n]=0 and E[nnH]=S, where S has full rank.

In a first optional step, the channel model is converted to its real-valued equivalent, see complex2real_channel(). We assume in the sequel the complex-valued representation. Then, the channel is whitened using whiten_channel():

y~=S12y=S12Hx+S12n=H~x+n~.

Next, the columns of H~ are sorted according to their norm in descending order. Then, the QR decomposition of the resulting channel matrix is computed:

H~=QR

where QCM×S is unitary and RCS×S is upper-triangular. The channel outputs are then pre-multiplied by QH. This leads to the final channel model on which the K-Best detection algorithm operates:

y¯=Rx¯+n¯

where y¯CS, x¯CS, and n¯CS with E[n¯]=0 and E[n¯n¯H]=I.

LLR Computation

The K-Best algorithm produces K candidate solutions x¯kCS and their associated distance metrics dk=y¯Rx¯k2 for k=1,,K. If the real-valued channel representation is used, the distance metrics are scaled by 0.5 to account for the reduced noise power in each complex dimension. A hard-decision is simply the candidate with the shortest distance. Various ways to compute LLRs from this list (and possibly additional side-information) are possible. The (sub-optimal) default solution is List2LLRSimple. Custom solutions can be provided.

Parameters:
  • output ("bit" | "symbol") – Type of output, either LLRs on bits or logits on constellation symbols

  • num_streams (int) – Number of transmitted streams

  • k (int`) – Number of paths to keep. Cannot be larger than the number of constellation points to the power of the number of streams.

  • constellation_type ("qam" | "pam" | "custom") – For “custom”, an instance of Constellation must be provided.

  • num_bits_per_symbol (int) – Number of bits per constellation symbol, e.g., 4 for QAM16. Only required for constellation_type in [“qam”, “pam”].

  • constellation (None (default) | Constellation) – If None, constellation_type and num_bits_per_symbol must be provided.

  • hard_out (bool, (default False)) – If True, the detector computes hard-decided bit values or constellation point indices instead of soft-values.

  • use_real_rep (bool, (default False)) – If True, the detector use the real-valued equivalent representation of the channel. Note that this only works with a QAM constellation.

  • list2llr (None (default) | List2LLR) – The function to be used to compute LLRs from a list of candidate solutions. If None, the default solution List2LLRSimple is used.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,num_streams], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

Output:
  • One of

  • […, num_streams, num_bits_per_symbol], tf.float – LLRs or hard-decisions for every bit of every stream, if output equals “bit”.

  • […, num_streams, num_points], tf.float or […, num_streams], tf.int32 – Logits or hard-decisions for constellation symbols for every stream, if output equals “symbol”. Hard-decisions correspond to the symbol indices.

property list2llr

Set/get the function to be used to compute LLRs from a list of candidate solutions

Type:

List2LLR

class sionna.phy.mimo.LinearDetector(equalizer, output, demapping_method, constellation_type=None, num_bits_per_symbol=None, constellation=None, hard_out=False, precision=None, **kwargs)[source]

Convenience class that combines an equalizer, such as lmmse_equalizer(), and a Demapper.

Parameters:
  • equalizer ("lmmse" | "zf" | "mf" | equalizer function) – The equalizer to be used. Either one of the existing equalizers lmmse_equalizer(), zf_equalizer(), or mf_equalizer() can be used, or a custom equalizer callable provided that has the same input/output specification.

  • output ("bit" | "symbol") – Type of output, either LLRs on bits or logits on constellation symbols

  • demapping_method ("app" | "maxlog") – Demapping method to be used

  • constellation_type ("qam" | "pam" | "custom") – For “custom”, an instance of Constellation must be provided.

  • num_bits_per_symbol (int) – Number of bits per constellation symbol, e.g., 4 for QAM16. Only required for constellation_type in [“qam”, “pam”].

  • constellation (None (default) | Constellation) – If None, constellation_type and num_bits_per_symbol must be provided.

  • hard_out (bool, (default False)) – If True, the detector computes hard-decided bit values or constellation point indices instead of soft-values.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,num_streams], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

Output:
  • One of

  • […, num_streams, num_bits_per_symbol], tf.float – LLRs or hard-decisions for every bit of every stream, if output equals “bit”

  • […, num_streams, num_points], tf.float or […, num_streams], tf.int32 – Logits or hard-decisions for constellation symbols for every stream, if output equals “symbol” Hard-decisions correspond to the symbol indices.

class sionna.phy.mimo.MaximumLikelihoodDetector(output, demapping_method, num_streams, constellation_type=None, num_bits_per_symbol=None, constellation=None, hard_out=False, precision=None, **kwargs)[source]

MIMO maximum-likelihood (ML) detector

This block implements MIMO maximum-likelihood (ML) detection assuming the following channel model:

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols which are uniformly and independently drawn from the constellation C, HCM×K is the known channel matrix, and nCM is a complex Gaussian noise vector. It is assumed that E[n]=0 and E[nnH]=S, where S has full rank. Optionally, prior information of the transmitted signal x can be provided, either as LLRs on the bits mapped onto x or as logits on the individual constellation points forming x.

Prior to demapping, the received signal is whitened:

y~=S12y=S12Hx+S12n=H~x+n~

The block can compute ML detection of symbols or bits with either soft- or hard-decisions. Note that decisions are computed symbol-/bit-wise and not jointly for the entire vector x (or the underlying vector of bits).

ML detection of bits:

Soft-decisions on bits are called log-likelihood ratios (LLR). With the “app” demapping method, the LLR for the ith bit of the kth user is then computed according to

LLR(k,i)=ln(Pr(bk,i=1|y,H)Pr(bk,i=0|y,H))=ln(xCk,i,1exp(y~H~x2)Pr(x)xCk,i,0exp(y~H~x2)Pr(x))

where Ck,i,1 and Ck,i,0 are the sets of vectors of constellation points for which the ith bit of the kth user is equal to 1 and 0, respectively. Pr(x) is the prior distribution of the vector of constellation points x. Assuming that the constellation points and bit levels are independent, it is computed from the prior of the bits according to

Pr(x)=k=1Ki=1Iσ(LLRp(k,i))

where LLRp(k,i) is the prior knowledge of the ith bit of the kth user given as an LLR and which is set to 0 if no prior knowledge is assumed to be available, and σ() is the sigmoid function. The definition of the LLR has been chosen such that it is equivalent with that of logit. This is different from many textbooks in communications, where the LLR is defined as LLR(k,i)=ln(Pr(bk,i=0|y,H)Pr(bk,i=1|y,H)).

With the “maxlog” demapping method, the LLR for the ith bit of the kth user is approximated like

LLR(k,i)ln(maxxCk,i,1(exp(y~H~x2)Pr(x))maxxCk,i,0(exp(y~H~x2)Pr(x)))=minxCk,i,0(y~H~x2ln(Pr(x)))minxCk,i,1(y~H~x2ln(Pr(x))).

ML detection of symbols:

Soft-decisions on symbols are called logits (i.e., unnormalized log-probability).

With the “app” demapping method, the logit for the constellation point cC of the kth user is computed according to

logit(k,c)=ln(x:xk=cexp(y~H~x2)Pr(x)).

With the “maxlog” demapping method, the logit for the constellation point cC of the kth user is approximated like

logit(k,c)maxx:xk=c(y~H~x2+ln(Pr(x))).

When hard decisions are requested, this block returns for the k th stream

c^k=argmaxcC(x:xk=cexp(y~H~x2)Pr(x))

where C is the set of constellation points.

Parameters:
  • output ("bit" | "symbol") – Type of output, either LLRs on bits or logits on constellation symbols

  • demapping_method ("app" | "maxlog") – Demapping method to be used

  • num_streams (int) – Number of transmitted streams

  • constellation_type ("qam" | "pam" | "custom") – For “custom”, an instance of Constellation must be provided.

  • num_bits_per_symbol (int) – Number of bits per constellation symbol, e.g., 4 for QAM16. Only required for constellation_type in [“qam”, “pam”].

  • constellation (None (default) | Constellation) – If None, constellation_type and num_bits_per_symbol must be provided.

  • hard_out (bool, (default False)) – If True, the detector computes hard-decided bit values or constellation point indices instead of soft-values.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,num_streams], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

  • prior (None (default) | […,num_streams,num_bits_per_symbol] or […,num_streams,num_points], tf.float) – Prior of the transmitted signals. If output equals “bit”, then LLRs of the transmitted bits are expected. If output equals “symbol”, then logits of the transmitted constellation points are expected.

Output:
  • One of

  • […, num_streams, num_bits_per_symbol], tf.float – LLRs or hard-decisions for every bit of every stream, if output equals “bit”.

  • […, num_streams, num_points], tf.float or […, num_streams], tf.int32 – Logits or hard-decisions for constellation symbols for every stream, if output equals “symbol”. Hard-decisions correspond to the symbol indices.

class sionna.phy.mimo.MMSEPICDetector(output, demapping_method='maxlog', num_iter=1, constellation_type=None, num_bits_per_symbol=None, constellation=None, hard_out=False, precision=None, **kwargs)[source]

Minimum mean square error (MMSE) with parallel interference cancellation (PIC) detector

This block implements the MMSE PIC detector, as proposed in [CST2011]. For num_iter>1, this implementation performs MMSE PIC self-iterations. MMSE PIC self-iterations can be understood as a concatenation of MMSE PIC detectors from [CST2011], which forward intrinsic LLRs to the next self-iteration.

Compared to [CST2011], this implementation also accepts priors on the constellation symbols as an alternative to priors on the bits.

This block assumes the following channel model:

y=Hx+n

where yCM is the received signal vector, xCS is the vector of transmitted symbols which are uniformly and independently drawn from the constellation C, HCM×S is the known channel matrix, and nCM is a complex Gaussian noise vector. It is assumed that E[n]=0 and E[nnH]=S, where S has full rank.

The algorithm starts by computing the soft symbols x¯s=E[xs] and variances vs=E[|es|2] from the priors, where es=xsx¯s, for all s=1,,S.

Next, for each stream, the interference caused by all other streams is cancelled from the observation y, leading to

y^s=yjshjxj=hsxs+n~s,s=1,,S

where n~s=jshjej+n.

Then, a linear MMSE filter ws is computed to reduce the resdiual noise for each observation y^s, which is given as

ws=hsH(HDsHH+S)1

where DsCS×S is diagonal with entries

[Ds]i,i={viis1i=s.

The filtered observations

z~s=wsHy^s=μ~sxs+wsHn~s

where μ~s=wsHhs, are then demapped to either symbol logits or LLRs, assuming that the remaining noise is Gaussian with variance

νs2=Var[z~s]=wsH(jshjhjHvj+S)ws.

The resulting soft-symbols can then be used for the next self-iteration of the algorithm.

Note that this algorithm can be substantially simplified as described in [CST2011] to avoid the computation of different matrix inverses for each stream. This is the version which is implemented.

Parameters:
  • output ("bit" | "symbol") – Type of output, either LLRs on bits or logits on constellation symbols

  • demapping_method ("maxlog" (default) | "app") – Demapping method to be used

  • num_iter (int, (default 1)) – Number of MMSE PIC iterations

  • constellation_type ("qam" | "pam" | "custom") – For “custom”, an instance of Constellation must be provided.

  • num_bits_per_symbol (int) – Number of bits per constellation symbol, e.g., 4 for QAM16. Only required for constellation_type in [“qam”, “pam”].

  • constellation (None (default) | Constellation) – If None, constellation_type and num_bits_per_symbol must be provided.

  • hard_out (bool, (default False)) – If True, the detector computes hard-decided bit values or constellation point indices instead of soft-values.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex) – Received signals

  • h ([…,M,num_streams], tf.complex) – Channel matrices

  • s ([…,M,M], tf.complex) – Noise covariance matrices

  • prior ([…,num_streams,num_bits_per_symbol] or […,num_streams,num_points], tf.float) – Prior of the transmitted signals. If output equals “bit”, then LLRs of the transmitted bits are expected. If output equals “symbol”, then logits of the transmitted constellation points are expected.

Output:
  • One of

  • […, num_streams, num_bits_per_symbol], tf.float – LLRs or hard-decisions for every bit of every stream, if output equals “bit”.

  • […, num_streams, num_points], tf.float or […, num_streams], tf.int32 – Logits or hard-decisions for constellation symbols for every stream, if output equals “symbol”. Hard-decisions correspond to the symbol indices.

Utility Functions

class sionna.phy.mimo.List2LLR(precision=None, **kwargs)[source]

Abstract class defining a callable to compute LLRs from a list of candidate vectors (or paths) provided by a MIMO detector

The following channel model is assumed

y¯=Rx¯+n¯

where y¯CS are the channel outputs, RCS×S is an upper-triangular matrix, x¯CS is the transmitted vector whose entries are uniformly and independently drawn from the constellation C, and n¯CS is white noise with E[n¯]=0 and E[n¯n¯H]=I.

It is assumed that a MIMO detector such as KBestDetector produces K candidate solutions x¯kCS and their associated distance metrics dk=y¯Rx¯k2 for k=1,,K. This layer can also be used with the real-valued representation of the channel.

Parameters:

precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex or tf.float) – Channel outputs of the whitened channel

  • r ([…,num_streams, num_streams], same dtype as y) – Upper triangular channel matrix of the whitened channel

  • dists ([…,num_paths], tf.float) – Distance metric for each path (or candidate)

  • path_inds ([…,num_paths,num_streams], tf.int32) – Symbol indices for every stream of every path (or candidate)

  • path_syms ([…,num_path,num_streams], same dtype as y) – Constellation symbol for every stream of every path (or candidate)

Output:

llr ([…num_streams,num_bits_per_symbol], tf.float) – LLRs for all bits of every stream

Note

An implementation of this class does not need to make use of all of the provided inputs which enable various different implementations.

class sionna.phy.mimo.List2LLRSimple(num_bits_per_symbol, llr_clip_val=20.0, precision=None, **kwargs)[source]

Computes LLRs from a list of candidate vectors (or paths) provided by a MIMO detector

The following channel model is assumed:

y¯=Rx¯+n¯

where y¯CS are the channel outputs, RCS×S is an upper-triangular matrix, x¯CS is the transmitted vector whose entries are uniformly and independently drawn from the constellation C, and n¯CS is white noise with E[n¯]=0 and E[n¯n¯H]=I.

It is assumed that a MIMO detector such as KBestDetector produces K candidate solutions x¯kCS and their associated distance metrics dk=y¯Rx¯k2 for k=1,,K. This layer can also be used with the real-valued representation of the channel.

The LLR for the ith bit of the kth stream is computed as

LLR(k,i)=log(Pr(bk,i=1|y¯,R)Pr(bk,i=0|y¯,R))minjCk,i,0djminjCk,i,1dj

where Ck,i,1 and Ck,i,0 are the set of indices in the list of candidates for which the ith bit of the kth stream is equal to 1 and 0, respectively. The LLRs are clipped to ±LLRclip which can be configured through the parameter llr_clip_val.

If Ck,i,0 is empty, LLR(k,i)=LLRclip; if Ck,i,1 is empty, LLR(k,i)=LLRclip.

Parameters:
  • num_bits_per_symbol (int) – Number of bits per constellation symbol

  • llr_clip_val (float, (default 20.0)) – The absolute values of LLRs are clipped to this value.

  • precision (None (default) | “single” | “double”) – Precision used for internal calculations and outputs. If set to None, precision is used.

Input:
  • y ([…,M], tf.complex or tf.float) – Channel outputs of the whitened channel

  • r ([…,num_streams, num_streams], same dtype as y) – Upper triangular channel matrix of the whitened channel

  • dists ([…,num_paths], tf.float) – Distance metric for each path (or candidate)

  • path_inds ([…,num_paths,num_streams], tf.int32) – Symbol indices for every stream of every path (or candidate)

  • path_syms ([…,num_path,num_streams], same dtype as y) – Constellation symbol for every stream of every path (or candidate)

Output:

llr ([…num_streams,num_bits_per_symbol], tf.float) – LLRs for all bits of every stream

property llr_clip_val

Get/set the value to which the absolute values of LLRs are clipped

Type:

float

sionna.phy.mimo.complex2real_vector(z)[source]

Transforms a complex-valued vector into its real-valued equivalent

Transforms the last dimension of a complex-valued tensor into its real-valued equivalent by stacking the real and imaginary parts on top of each other.

For a vector zCM with real and imaginary parts xRM and yRM, respectively, this function returns the vector [xT,yT]TR2M.

Input:

[…,M], tf.complex

Output:

[…,2M], tf.complex.real_dtype

sionna.phy.mimo.real2complex_vector(z)[source]

Transforms a real-valued vector into its complex-valued equivalent

Transforms the last dimension of a real-valued tensor into its complex-valued equivalent by interpreting the first half as the real and the second half as the imaginary part.

For a vector z=[xT,yT]TR2M with xRM and yRM, this function returns the vector x+jyCM.

Input:

[…,2M], tf.float

Output:

[…,M], tf.complex

sionna.phy.mimo.complex2real_matrix(z)[source]

Transforms a complex-valued matrix into its real-valued equivalent

Transforms the last two dimensions of a complex-valued tensor into their real-valued matrix equivalent representation.

For a matrix ZCM×K with real and imaginary parts XRM×K and YRM×K, respectively, this function returns the matrix Z~R2M×2K, given as

Z~=(XYYX).
Input:

[…,M,K], tf.complex

Output:

[…,2M, 2K], tf.complex.real_dtype

sionna.phy.mimo.real2complex_matrix(z)[source]

Transforms a real-valued matrix into its complex-valued equivalent

Transforms the last two dimensions of a real-valued tensor into their complex-valued matrix equivalent representation.

For a matrix Z~R2M×2K, satisfying

Z~=(XYYX)

with XRM×K and YRM×K, this function returns the matrix Z=X+jYCM×K.

Input:

[…,2M,2K], tf.float

Output:

[…,M, 2], tf.complex

sionna.phy.mimo.complex2real_covariance(r)[source]

Transforms a complex-valued covariance matrix to its real-valued equivalent

Assume a proper complex random variable zCM [ProperRV] with covariance matrix R=∈CM×M and real and imaginary parts xRM and yRM, respectively. This function transforms the given R into the covariance matrix of the real-valued equivalent vector z~=[xT,yT]TR2M, which is computed as [CovProperRV]

E[z~z~H]=(12{R}12{R}12{R}12{R}).
Input:

[…,M,M], tf.complex

Output:

[…,2M, 2M], tf.complex.real_dtype

sionna.phy.mimo.real2complex_covariance(q)[source]

Transforms a real-valued covariance matrix to its complex-valued equivalent

Assume a proper complex random variable zCM [ProperRV] with covariance matrix R=∈CM×M and real and imaginary parts xRM and yRM, respectively. This function transforms the given covariance matrix of the real-valued equivalent vector z~=[xT,yT]TR2M, which is given as [CovProperRV]

E[z~z~H]=(12{R}12{R}12{R}12{R}),

into is complex-valued equivalent R.

Input:

[…,2M,2M], tf.float

Output:

[…,M, M], tf.complex

sionna.phy.mimo.complex2real_channel(y, h, s)[source]

Transforms a complex-valued MIMO channel into its real-valued equivalent

Assume the canonical MIMO channel model

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector with covariance matrix SCM×M.

This function returns the real-valued equivalent representations of y, H, and S, which are used by a wide variety of MIMO detection algorithms (Section VII) [YH2015]. These are obtained by applying complex2real_vector() to y, complex2real_matrix() to H, and complex2real_covariance() to S.

Input:
  • y ([…,M], tf.complex) – Complex-valued received signals

  • h ([…,M,K], tf.complex) – Complex-valued channel matrices

  • s ([…,M,M], tf.complex) – Complex-valued noise covariance matrices

Output:
  • […,2M], tf.complex.real_dtype – Real-valued equivalent received signals

  • […,2M,2K], tf.complex.real_dtype – Real-valued equivalent channel matrices

  • […,2M,2M], tf.complex.real_dtype – Real-valued equivalent noise covariance matrices

sionna.phy.mimo.real2complex_channel(y, h, s)[source]

Transforms a real-valued MIMO channel into its complex-valued equivalent

Assume the canonical MIMO channel model

y=Hx+n

where yCM is the received signal vector, xCK is the vector of transmitted symbols, HCM×K is the known channel matrix, and nCM is a noise vector with covariance matrix SCM×M.

This function transforms the real-valued equivalent representations of y, H, and S, as, e.g., obtained with the function complex2real_channel(), back to their complex-valued equivalents (Section VII) [YH2015].

Input:
  • y ([…,2M], tf.float) – Real-valued received signals

  • h ([…,2M,2K], tf.float) – Real-valued channel matrices

  • s ([…,2M,2M], tf.float) – Real-valued noise covariance matrices

Output:
  • […,M], tf.complex – Complex-valued equivalent received signals

  • […,M,K], tf.complex – Complex-valued equivalent channel matrices

  • […,M,M], tf.complex – Complex-valued equivalent noise covariance matrices

sionna.phy.mimo.whiten_channel(y, h, s, return_s=True)[source]

Whitens a canonical MIMO channel

Assume the canonical MIMO channel model

y=Hx+n

where yCM(RM) is the received signal vector, xCK(RK) is the vector of transmitted symbols, HCM×K(RM×K) is the known channel matrix, and nCM(RM) is a noise vector with covariance matrix SCM×M(RM×M).

This function whitens this channel by multiplying y and H from the left by S12=L1, where LCM×M is the Cholesky decomposition of S. Optionally, the whitened noise covariance matrix IM can be returned.

Input:
  • y ([…,M], tf.float or tf.complex) – Received signals

  • h ([…,M,K], tf.float or tf.complex) – Channel matrices

  • s ([…,M,M], tf.float or tf.complex) – Noise covariance matrices

  • return_s (bool, (default True)) – If True, the whitened covariance matrix is returned.

Output:
  • […,M], tf.float or tf.complex – Whitened received signals

  • […,M,K], tf.float or tf.complex – Whitened channel matrices

  • […,M,M], tf.float or tf.complex – Whitened noise covariance matrices. Only returned if return_s is True.

References:
[ProperRV] (1,2)

Proper complex random variables, Wikipedia, accessed 11 September, 2022.

[CovProperRV] (1,2)

Covariance matrices of real and imaginary parts, Wikipedia, accessed 11 September, 2022.

[YH2015] (1,2)

S. Yang and L. Hanzo, “Fifty Years of MIMO Detection: The Road to Large-Scale MIMOs”, IEEE Communications Surveys & Tutorials, vol. 17, no. 4, pp. 1941-1988, 2015.

[FT2015]

W. Fu and J. S. Thompson, “Performance analysis of K-best detection with adaptive modulation”, IEEE Int. Symp. Wirel. Commun. Sys. (ISWCS), 2015.

[EP2014]

J. Céspedes, P. M. Olmos, M. Sánchez-Fernández, and F. Perez-Cruz, “Expectation Propagation Detection for High-Order High-Dimensional MIMO Systems”, IEEE Trans. Commun., vol. 62, no. 8, pp. 2840-2849, Aug. 2014.

[CST2011] (1,2,3,4)

C. Studer, S. Fateh, and D. Seethaler, “ASIC Implementation of Soft-Input Soft-Output MIMO Detection Using MMSE Parallel Interference Cancellation”, IEEE Journal of Solid-State Circuits, vol. 46, no. 7, pp. 1754–1765, July 2011.