SSFM#

class sionna.phy.channel.SSFM(alpha: float = 0.046, beta_2: float = -21.67, f_c: float = 193550000000000.0, gamma: float = 1.27, half_window_length: int = 0, length: float = 80, n_ssfm: int | str = 1, n_sp: float = 1.0, sample_duration: float = 1.0, t_norm: float = 1e-12, with_amplification: bool = False, with_attenuation: bool = True, with_dispersion: bool = True, with_manakov: bool = False, with_nonlinearity: bool = True, phase_inc: float = 0.0001, precision: str | None = None, device: str | None = None, **kwargs)[source]#

Bases: sionna.phy.block.Block

Block implementing the split-step Fourier method (SSFM)

The SSFM (first mentioned in [HT1973]) numerically solves the generalized nonlinear Schrödinger equation (NLSE)

\[\frac{\partial E(t,z)}{\partial z}=-\frac{\alpha}{2} E(t,z)+j\frac{\beta_2}{2}\frac{\partial^2 E(t,z)}{\partial t^2}-j\gamma |E(t,z)|^2 E(t,z) + n(n_{\text{sp}};\,t,\,z)\]

for an unpolarized (or single polarized) optical signal; or the Manakov equation (according to [WMC1991])

\[\frac{\partial \mathbf{E}(t,z)}{\partial z}=-\frac{\alpha}{2} \mathbf{E}(t,z)+j\frac{\beta_2}{2}\frac{\partial^2 \mathbf{E}(t,z)}{\partial t^2}-j\gamma \frac{8}{9}||\mathbf{E}(t,z)||_2^2 \mathbf{E}(t,z) + \mathbf{n}(n_{\text{sp}};\,t,\,z)\]

for dual polarization, with attenuation coefficient \(\alpha\), group velocity dispersion parameters \(\beta_2\), and nonlinearity coefficient \(\gamma\). The noise terms \(n(n_{\text{sp}};\,t,\,z)\) and \(\mathbf{n}(n_{\text{sp}};\,t,\,z)\), respectively, stem from an (optional) ideally distributed Raman amplification with spontaneous emission factor \(n_\text{sp}\). The optical signal \(E(t,\,z)\) has the unit \(\sqrt{\text{W}}\). For the dual polarized case, \(\mathbf{E}(t,\,z)=(E_x(t,\,z), E_y(t,\,z))\) is a vector consisting of the signal components of both polarizations.

The symmetrized SSFM is applied according to Eq. (7) of [FMF1976] that can be written as

\[E(z+\Delta_z,t) \approx \exp\left(\frac{\Delta_z}{2}\hat{D}\right)\exp\left(\int^{z+\Delta_z}_z \hat{N}(z')dz'\right)\exp\left(\frac{\Delta_z}{2}\hat{D}\right)E(z,\,t)\]

where only the single-polarized case is shown. The integral is approximated by \(\Delta_z\hat{N}\) with \(\hat{D}\) and \(\hat{N}\) denoting the linear and nonlinear SSFM operator, respectively [A2012].

Additionally, ideally distributed Raman amplification may be applied, which is implemented as in [MFFP2009]. Please note that the implemented Raman amplification currently results in a transparent fiber link. Hence, the introduced gain cannot be parametrized.

The SSFM operates on normalized time \(T_\text{norm}\) (e.g., \(T_\text{norm}=1\,\text{ps}=1\cdot 10^{-12}\,\text{s}\)) and distance units \(L_\text{norm}\) (e.g., \(L_\text{norm}=1\,\text{km}=1\cdot 10^{3}\,\text{m}\)). Hence, all parameters as well as the signal itself have to be given with the same unit prefix for the same unit (e.g., always pico for time, or kilo for distance). Despite the normalization, the SSFM is implemented with physical units, which is different from the normalization, e.g., used for the nonlinear Fourier transform. For simulations, only \(T_\text{norm}\) has to be provided.

To avoid reflections at the signal boundaries during simulation, a Hamming window can be applied in each SSFM-step, whose length can be defined by half_window_length.

Parameters:
  • alpha (float) – Attenuation coefficient \(\alpha\) in \((1/L_\text{norm})\). Defaults to 0.046.

  • beta_2 (float) – Group velocity dispersion coefficient \(\beta_2\) in \((T_\text{norm}^2/L_\text{norm})\). Defaults to -21.67.

  • f_c (float) – Carrier frequency \(f_\mathrm{c}\) in \((\text{Hz})\). Defaults to 193.55e12.

  • gamma (float) – Nonlinearity coefficient \(\gamma\) in \((1/L_\text{norm}/\text{W})\). Defaults to 1.27.

  • half_window_length (int) – Half of the Hamming window length. Defaults to 0.

  • length (float) – Fiber length \(\ell\) in \((L_\text{norm})\). Defaults to 80.0.

  • n_ssfm (int | str) – Number of steps \(N_\mathrm{SSFM}\). Set to “adaptive” to use nonlinear-phase rotation to calculate the step widths adaptively (maximum rotation can be set in phase_inc). Defaults to 1.

  • n_sp (float) – Spontaneous emission factor \(n_\mathrm{sp}\) of Raman amplification. Defaults to 1.0.

  • sample_duration (float) – Normalized time step \(\Delta_t\) in \((T_\text{norm})\). Defaults to 1.0.

  • t_norm (float) – Time normalization \(T_\text{norm}\) in \((\text{s})\). Defaults to 1e-12.

  • with_amplification (bool) – If True, enables ideal inline amplification and corresponding noise. Defaults to False.

  • with_attenuation (bool) – If True, enables attenuation. Defaults to True.

  • with_dispersion (bool) – If True, applies chromatic dispersion. Defaults to True.

  • with_manakov (bool) – If True, considers axis [-2] as x- and y-polarization and calculates the nonlinear step as given by the Manakov equation. Defaults to False.

  • with_nonlinearity (bool) – If True, applies Kerr nonlinearity. Defaults to True.

  • phase_inc (float) – Maximum nonlinear-phase rotation in rad allowed during simulation. To be used with n_ssfm = “adaptive”. Defaults to 1e-4.

  • precision (str | None) – Precision used for internal calculations and outputs. If set to None, precision is used.

  • device (str | None) – Device for computation. If None, device is used.

Inputs:

x – […, n] or […, 2, n], torch.complex. Input signal in \((\sqrt{\text{W}})\). If with_manakov is True, the second last dimension is interpreted as x- and y-polarization, respectively.

Outputs:

y – Tensor (same shape as x), torch.complex. Channel output.

Examples

import torch
from sionna.phy.channel.optical import SSFM

ssfm = SSFM(
    alpha=0.046,
    beta_2=-21.67,
    f_c=193.55e12,
    gamma=1.27,
    half_window_length=100,
    length=80,
    n_ssfm=200,
    n_sp=1.0,
    t_norm=1e-12,
    with_amplification=False,
    with_attenuation=True,
    with_dispersion=True,
    with_manakov=False,
    with_nonlinearity=True)

# x is the optical input signal
x = torch.randn(10, 1024, dtype=torch.complex64)
y = ssfm(x)
print(y.shape)
# torch.Size([10, 1024])