Tutorial: BMSTransformation/LorentzTransformation

The scri.bms_transformations.BMSTransformation() and scri.bms_transformations.LorentzTransformation() classes are a convenient way to store and work with BMS and Lorentz transformations. The LorentzTransformation class takes as input a frame rotation and a boost and then enables one to convert these transformations from their axis-angle representation to their $SL(2,\mathbb{C})$ representation. The BMSTransformation class is effectively the LorentzTransformation class, but also takes as input a supertranslation (which includes the usual spacetime translations).

The BMSTransformation and LorentzTransformation classes make it simple to reorder the individual components of a transformation, compute the inverse of a transformation, and compose two transformations.

Theory behind the Code

We consider a spacetime $(\hat{M},\hat{g}_{ab})$ to be asymptotically Minkowskian at null infinity if there exists a manifold $M$ equipped with a metric $g_{ab}$ and boundary $\mathcal{I}$, and a diffeomorphism from $\hat{M}$ to the interior $M\backslash I$ such that:

  1. there exists a smooth function $\Omega$ (the conformal factor) on $M$ such that $g_{ab}=\Omega^{2}\hat{g}_{ab}$ on $\hat{M}$, $\Omega=0$ on $\mathcal{I}$, and $n_{a}=\nabla_{a}\Omega$ is nowhere vanishing on $\mathcal{I}$;

  2. $\mathcal{I}$ is topologically $\mathbb{S}^{2}\times\mathbb{R}$;

  3. $\hat{g}_{ab}$ satisfies Einstein’s equations $\hat{R}_{ab}-\frac{1}{2}\hat{R}\hat{g}_{ab}=8\pi G\hat{T}_{ab}$, where $\Omega^{-2}\hat{T}_{ab}$ has a smooth limit to $\mathcal{I}$; and,

  4. the integral curves of $n^{a}$ are complete on $\mathcal{I}$ for any choice of the conformal factor such that $\nabla_{a}n^{a}=0$ on $\mathcal{I}$.

For these spacetimes, their asymptotic symmetry group is known as BMS group $\mathfrak{B}$. Formally, the BMS group is the semi-direct product of the Lorentz group $\mathcal{L}\simeq SL(2,\mathbb{C})$ and an infinite-dimensional (Abelian) group $\mathcal{S}$ of supertranslations $\mathcal{S}$. That is, $\mathfrak{B}\equiv SL(2,\mathbb{C})\ltimes\mathcal{S}$, where $\mathcal{S}$ is the group generated by vector fields $fn^{a}$ on $\mathcal{I}$ for $f$ a real-valued scalar function on $\mathbb{S}^{2}$ satisfying $\mathcal{L}_{n}f=0$, with $\mathcal{L}_{n}f$ being the Lie derivative of $f$ along the vector field $n$. (Note that $\mathcal{S}$ is an extension of the usual translation group $\mathcal{T}$; this can be realized by taking $f$ to be the $\ell\leq1$ spherical harmonics.)

Consequently, any element of $\mathfrak{B}$ can be written as a supertranslation followed by a rotation and a boost, seeing as any Lorentz transformation itself can be written as a rotation followed by a boost or vice versa.

Determining the $SL(2,\mathbb{C})$ representation of a Lorentz transformation

For a rotation $R$ by angle $\theta$ about the axis $\hat{r}=\left(r_{x},r_{y},r_{z}\right)$, one may write this as the quaternion

$$\mathbf{q}=\exp\left(\frac{1}{2}\theta\,\hat{r}\right)= \cos(\theta/2)\mathbf{I}+\left(r_{x}\mathbf{i}+r_{y}\mathbf{j}+r_{z}\mathbf{k}\right)\sin(\theta/2),$$

where $\mathbf{I}$, $\mathbf{i}$, $\mathbf{j}$, and $\mathbf{k}$ are the elementary quaternions obeying the usual multiplication rules. Using spin matrices, i.e., elements of $SL(2,\mathbb{C})$, one may also write these elementary quaternions as

$$\mathbf{\sigma}\equiv\left(\mathbf{I},\mathbf{i},\mathbf{j},\mathbf{k}\right)= \left(\begin{pmatrix}1&0\\0&1\end{pmatrix},\begin{pmatrix}0&i\\i&0\end{pmatrix}, \begin{pmatrix}0&-1\\1&0\end{pmatrix},\begin{pmatrix}i&0\\0&-i\end{pmatrix}\right).$$

Therefore, a rotation by angle $\theta$ about the axis $\hat{r}=\left(r_{x},r_{y},r_{z}\right)$ has the $SL(2,\mathbb{C})$ representation

$$\tilde{R}=\left(\cos(\theta/2),\hat{r}\sin(\theta/2)\right)\mathbf{\sigma}.$$

Note that in this equation $\tilde{R}$ is a unitary matrix. Furthermore, since a boost $B$ by rapidity $w$ along the axis $\hat{v}=\left(v_{x},v_{y},v_{z}\right)$ is nothing more than a rotation by angle $iw$ about the axis $\hat{v}=\left(v_{x},v_{y},v_{z}\right)$, a boost has the $SL(2,\mathbb{C})$ representation

$$\tilde{B}=\left(\cosh(w/2),i\hat{v}\sinh(w/2)\right)\cdot\mathbf{\sigma}.$$

Note that in this equation $\tilde{B}$ is a Hermitian matrix; by a similarity transformation, it can be brought to the form $\text{diag}(\exp(w/2),\exp(-w/2))$. Consequently, a Lorentz transformation $L=B\cdot R$ has the $SL(2,\mathbb{C})$ representation

$$\tilde{L}=\tilde{B}\cdot\tilde{R}=\left[\left(\cosh(w/2),i\hat{v}\sinh(w/2)\right)\cdot\mathbf{\sigma}\right] \cdot\left[\left(\cos(\theta/2),\hat{r}\sin(\theta/2)\right)\cdot\mathbf{\sigma}\right].$$

Let us make the following remark. What is shown above is for a passive Lorentz transformation. Put differently, the transformation is acting on the coordinates themselves. Under $L$, our original Minkowski coordinates $X$ are transformed to an intermediate coordinate system $X’$ by a rotation, and then to a final coordinate system $X’’$ by a boost. The parameters of the boost vector $\vec{v}$ can be interpreted as a boost vector in the $X’$ coordinate system.

Determining the rotation and boost from a $SL(2,\mathbb{C})$ representation

Given the $SL(2,\mathbb{C})$ representation of a Lorentz transformation $L$, say

$$\tilde{L}=\begin{pmatrix}a&b\\c&d\end{pmatrix}$$

for $a,b,c,d\in\mathbb{C}$ with $ad-bc=1$, using singular value decomposition we may write $\tilde{L}$ as

$$\tilde{L}=U\cdot\Sigma\cdot V^{\dagger},$$

where $U$ is a unitary matrix, $\Sigma$ is a diagonal Hermitian matrix, and $V$ is a complex unitary matrix with $V^{\dagger}$ being its conjugate transpose. Because of this, we may also write $\tilde{L}$ as

$$\tilde{L}=\left(U\cdot\Sigma\cdot U^{\dagger}\right)\cdot\left(U\cdot V^{\dagger}\right)$$

or

$$\tilde{L}=\left(U\cdot V^{\dagger}\right)\cdot\left(V\cdot\Sigma\cdot V^{\dagger}\right),$$

which may then be interpreted in terms of a pure rotation and a pure boost as

$$\tilde{L}=\tilde{B}\cdot\tilde{R}\quad\text{with}\quad \tilde{B}\equiv U\cdot\Sigma\cdot U^{\dagger}\quad\text{and} \quad\tilde{R}\equiv U\cdot V^{\dagger}$$

or

$$\tilde{L}=\tilde{R}’\cdot\tilde{B}’\quad\text{with}\quad \tilde{R}’\equiv U\cdot V^{\dagger}\quad\text{and} \quad\tilde{B}’\equiv V\cdot\Sigma\cdot V^{\dagger},$$

since pure rotations and pure boosts have unitary spin-matrix and positive definite Hermitian spin-matrix representations. These are examples of polar decompositions, which for invertible matrices are unique. In what follows we drop the primes on $\tilde{R}’$ and $\tilde{B}’$ in the final equation if using the form of $\tilde{L}$ in said equation. With these $SL(2,\mathbb{C})$ representations, one may then obtain the angle-axis representation of these pure Lorentz transformations via

$$ R=\left(\cos(\theta/2),\hat{r}\sin(\theta/2)\right)_{k}= \frac{1}{2}\eta_{k}\text{Tr}\left[\tilde{R}\cdot\mathbf{\sigma}_{k}\right];\\ B=\left(\cosh(w/2),i\hat{v}\sinh(w/2)\right)_{k}= \frac{1}{2}\eta_{k}\text{Tr}\left[\tilde{B}\cdot\mathbf{\sigma}_{k}\right],$$

where $\eta_{k}=(1,-1,-1,-1)_{k}$. With an understanding of how to map an arbitrary Lorentz transformation to its unique $SL(2,\mathbb{C})$ representation and back, reordering, inverting, and composing Lorentz transformations now becomes trivial. Thus, working with BMS transformations is also trivial, since supertranslations are normal in the BMS group.

Compose two elements of the BMS group

We will write a BMS element as

$$g=l\circ s$$

and avoid using $b$ or $B$, which could ambiguously stand for either BMS or boost. The first thing to note is that since supertranslations are normal in BMS, this decomposition is well-defined. Specifically, consider the quotient homomorphism $\varphi:\mathfrak{B}\to\mathfrak{B}/\mathcal{S}\simeq\mathcal{L}$. We can say that a BMS element is a ``pure supertranslation’’ if it is in $\ker\varphi$, i.e. $\varphi(s)=e\in\mathcal{L}$. The Lorentz transformation of a BMS transformation $L=\varphi(g)$ is well-defined. All the elements in the preimage $\varphi^{-1}(L)$ differ by pure supertranslations (elements of the kernel). But the element $l\in \varphi^{-1}(L)$ above is unique whether we choose to do a left or right decomposition: since supertranslations are a normal subgroup, we can write $s = l^{-1} \circ \hat{s} \circ l$, where $\hat{s}=l \circ s \circ l^{-1}$ is another supertranslation that is distinct unless $l=e$. This provides the other decomposition of $g$ in terms of a pure supertranlsation and a Lorentz transformation, $g=\hat{s}\circ l$.

With this, consider taking the composition of $g_{1},g_{2}\in\mathfrak{B}$,

$$g=g_2 \circ g_1$$

$$g=l_2 \circ s_2\circ l_1\circ s_1.$$

As before, since $s_2$ is in the normal subgroup, we can write it as conjugate to another pure supertranslation, namely

$$s_2=l_1\circ s’\circ l_1^{-1},$$

$$l_1^{-1}\circ s_2\circ l_1 = s’.$$

Combining and adding parentheses we get

$$g = (l_2\circ l_1)\circ(s’\circ s_1).$$

Creating a BMSTransformation

A BMSTransformation object consisting of a supertranslation, frame rotation, and boost (in that order) can be created via the following:

>>> S = np.array([1, 2 + 4j, 3, -2 + 4j, 7 - 5j, -3 - 2j, 4, 3 - 2j, 7 + 5j]) * 1e-3
>>> q = np.quaternion(1, 2, 3, 4).normalized()
>>> v = np.array([1, 2, 3]) * 1e-4

>>> BMS1 = bms_transformations.BMSTransformation(supertranslation=S,
 ...     frame_rotation=q.components,
 ...     boost_velocity=v)

This transformation can have its components reordered via, e.g.,

>>> BMS1.reorder(["boost_velocity", "supertranslation", "frame_rotation"])

where the component listed first is the transformation to be performed first. This transformation’s inverse can be computed via

>>> BMS1_inv = BMS1.inverse()

Finally, this transformation can be composed with another, BMS2, via

>>> BMS2 * BMS1

where here the order implies that we perform BMS1 first, followed by BMS2.