Example 1: Pulsebroadening by chromatic dispersion

A simple example to show the effect of chromatic dispersion is to show the pulse broadening of a single pulse, which has field $E$ has a gaussian shape with a $T_0$ 1/e-pulse width:
$$ E(0, t) = e^{-\frac{1}{2} \left( \frac{t^2}{T_0^2} \right)^2 } $$

According to Agrawal the pulse width increases with the transmission length $z$: $$ T_1(z) = T_0 \underbrace{ \sqrt{1 + \left( \frac{z}{L_D} \right)^2} }_{\text{Broadening factor}} $$ where $L_D = T_0^2 / |\beta_2|$ represents the dispersion length. This formula can be used to calculate the FWHM-pulse width, as well. For this the input FWHM-pulse width has to be multiplicated simply by the broadening factor.

In the following a signal with $T_{FWHM}=20ps$ pulse width is transmitted over $z=100km$ fiber with a dispersion coefficent $D=17 ps/nm/km$. The coefficent $\alpha$, $\beta_{1x}$, $\beta_{1y}$, $\beta_3$ and the non linear coeffiecent $\gamma$ is set to zero, so that only the effect of first order chromatic dispersion can be studied.

From the plot can be seen, that the FWHM pulse width is $T_{FWHM}=47.1 ps$. After the transmission the FWHM pulse width is $T_{FWHM}=259.72 ps$. This is in good accordance to the theorie, as the calculated FWHM pulse width is 259.59ps. The difference gets smaller with a higher number of samples per symbol.

Pulsebroadening by chromatic dispersion
Fig. 1 : Input and broadened output signal

Good to know!
Definition of pulse width

The FWHM pulse width of a gaussian pulse of the electrical field as defined above can be calculated by $T_{FWHM} = 2\sqrt{\ln4} T_0 \approx 2.35482 T_0$.
If the power of the pulse $|E(0, t)|^2$ is analysed the FWHM pulse width is given by $T_{FWHM,power} = \sqrt{0.5} T_{FWHM} = 2\sqrt{\ln2} T_0 \approx 1.665109 T_0$. The 1/e width is $T_{0,power} = \sqrt{0.5} T_{0}$ accordingly.

#Import functions and libraries
from pypho_setup import pypho_setup
from pypho_bits import pypho_bits
from pypho_signalsrc import pypho_signalsrc
from pypho_lasmod import pypho_lasmod
from pypho_fiber import pypho_fiber
from pypho_functions import *
import numpy as np
import copy
import matplotlib.pyplot as plt
from scipy.interpolate import UnivariateSpline

# Define network elements
gp       = pypho_setup(nos = 16, sps = 256, symbolrate = 10e9)
bitsrc   = pypho_bits(glova = gp, nob = gp.nos, pattern = 'singlepulse')
esigsrc  = pypho_signalsrc(glova = gp, pulseshape = 'gauss_rz' , fwhm = 0.30)
sig_1550 = pypho_lasmod(glova = gp, power = 0, Df = 0, teta = 0)
SSMF     = pypho_fiber(glova = gp, l = 80.0e3,  D = 17,   S = 0, alpha = 0.2e-12, gamma = 0, phi_max=10.0)

# Simulation
bits = bitsrc()
esig = esigsrc(bitsequence = bits)
E_Tx = sig_1550(esig = esig)                                                      


# Define your parameters here
T_0 = 20.0e-12
z = 100.0e3
D = 17.0
beta_2, beta_3 = DS2beta(17.0, 0, gp.lambda0)


# Create a single pulse with gaussian shape (not power!)
E_Tx[0]['E'][0] = E_Tx[0]['E'][0]*0 + np.exp(-(gp.timeax-gp.timeax[-1]/2)**2 / (2.0*T_0**2) )
E = copy.deepcopy(E_Tx)

# Fiber trannsmission
E = SSMF(E = E, D = D, l = z) 

# Plot Input and Output signal 
plt.figure(1)
plt.plot(gp.timeax*1.0e12, np.abs(E_Tx[0]['E'][0]), 'r', label='$E(0, t)$')
plt.plot(gp.timeax*1.0e12, np.abs(E[0]['E'][0]), 'g', label='$E(z=100$km$, t)$')


# Get FWHM of the input signal E_Tx
spline = UnivariateSpline(gp.timeax*1.0e12, np.abs(E_Tx[0]['E'][0])-1*np.max(np.abs(E_Tx[0]['E'][0]))/2, s=0)
r1, r2 = spline.roots() # find the roots
plt.annotate(s='', xy=(r1,np.max(np.abs(E_Tx[0]['E'][0]))/2), xytext=(r2,np.max(np.abs(E_Tx[0]['E'][0]))/2), arrowprops=dict(arrowstyle='<->'))
plt.text(r1+(r2-r1)/2.0, 0.01 +np.max(np.abs(E_Tx[0]['E'][0]))/2, '$T_{FWHM,0}$ = ' + str(np.round(r2-r1,2)) + ' ps', fontsize=12, horizontalalignment='center')
T_FWHM_0 = (r2-r1) * 1e-12
T_0_plot = T_FWHM_0 / 2.35482


# Get FWHM of the output signal E
spline = UnivariateSpline(gp.timeax*1.0e12, np.abs(E[0]['E'][0])-1*np.max(np.abs(E[0]['E'][0]))/2, s=0)
r1, r2 = spline.roots() # find the roots
plt.annotate(s='', xy=(r1,np.max(np.abs(E[0]['E'][0]))/2), xytext=(r2,np.max(np.abs(E[0]['E'][0]))/2), arrowprops=dict(arrowstyle='<->'))
plt.text(r1+(r2-r1)/2.0, 0.01 + np.max(np.abs(E[0]['E'][0]))/2, '$T_{FWHM,1}$ = ' + str(np.round(r2-r1,2)) + ' ps', fontsize=12, horizontalalignment='center')
T_FWHM_1 = (r2-r1) * 1e-12
plt.ylabel('$|E|$ a.u.'); plt.xlabel('Time $t$ [ps]'); legend = plt.legend(loc='upper right')

L_D = (T_0_plot)**2  / np.abs(beta_2)
# Print the results
print 'Input signal 1/e-pulse width by definition : T_0 = ' , T_0*1e12, ' ps'
print 'Input signal 1/e-pulse width from plot : T_0 = ' , T_0_plot*1e12, ' ps'
print 'Input signal FWHM-pulse width from plot : T_FWHM,0 = ' , T_FWHM_0*1e12, ' ps'
print 'Output signal FWHM-pulse width from plot : T_FWHM,1 = ' , T_FWHM_1*1e12, ' ps'
print 'Calculated output FWHM-pulse width : T_FWHM,1 = ' , T_0 * np.sqrt(1 + (z/L_D)**2)*1e12, ' ps'