Qsymm Basics

See also

The complete source code of this example can be found in basics.py. A Jupyter notebook can be found in basics.ipynb.

Getting started with Qsymm is as simple as importing it:

import qsymm

To make effective use of Qsymm we’ll also need a few other utilities: numpy handling numeric arrays, and sympy for symbolic mathematics:

import numpy as np
import sympy

In all the following tutorials we will use these standard imports, and they won’t be explicitly shown

Defining a Qsymm model

Let’s start by defining a 3D Rashba Hamiltonian symbolically as a Python string:

ham = ("hbar^2 / (2 * m) * (k_x**2 + k_y**2 + k_z**2) * eye(2) +" +
        "alpha * sigma_x * k_x + alpha * sigma_y * k_y + alpha * sigma_z * k_z")

We can then create a Qsymm Model directly from this symbolic Hamiltonian:

H = qsymm.Model(ham)

We can then directly inspect the contents by printing the Model:

print(H)
{alpha*k_z:
[[ 1.+0.j  0.+0.j]
 [ 0.+0.j -1.+0.j]],

hbar**2*k_x**2/m:
[[0.5+0.j 0. +0.j]
 [0. +0.j 0.5+0.j]],

hbar**2*k_y**2/m:
[[0.5+0.j 0. +0.j]
 [0. +0.j 0.5+0.j]],

hbar**2*k_z**2/m:
[[0.5+0.j 0. +0.j]
 [0. +0.j 0.5+0.j]],

alpha*k_x:
[[0.+0.j 1.+0.j]
 [1.+0.j 0.+0.j]],

alpha*k_y:
[[0.+0.j 0.-1.j]
 [0.+1.j 0.+0.j]],

}

We can also extract a more readable representation by using the tosympy method, which converts the Model to a sympy expression:

H.tosympy(nsimplify=True)
\[\begin{split}\displaystyle \left[\begin{matrix}\alpha k_{z} + \frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{y}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m} & \alpha k_{x} - i \alpha k_{y}\\\alpha k_{x} + i \alpha k_{y} & - \alpha k_{z} + \frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{y}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m}\end{matrix}\right]\end{split}\]

The argument nsimplify=True makes the output more readable by forcing sympy to elide factors of 1.0 that multiply each term. Note that Qsymm automatically interprets the symbols sigma_x, sigma_y and sigma_z as the Pauli matrices, and eye(2) as the 2x2 identity matrix.

Model as a momenta attribute that specifies which symbols are considered the momentum variables:

H.momenta
(k_x, k_y, k_z)

By default Qsymm assumes that your model is written in 3D (even if it does not include all 3 momenta). To define a lower-dimensional model you must explicitly specify the momentum variables, e.g:

ham2D = ("hbar^2 / (2 * m) * (k_x**2 + k_z**2) * eye(2) +" +
         "alpha * sigma_x * k_x + alpha * sigma_y * k_z")
H2D = qsymm.Model(ham2D, momenta=['k_x', 'k_z'])
H2D.tosympy(nsimplify=True)
\[\begin{split}\displaystyle \left[\begin{matrix}\frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m} & \alpha k_{x} - i \alpha k_{z}\\\alpha k_{x} + i \alpha k_{z} & \frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m}\end{matrix}\right]\end{split}\]
H2D.momenta
(k_x, k_z)

Defining group elements

Qsymm is all about finding and generating symmetries of models, so it is unsurprising that it contains utilities for defining group elements.

Below are a few examples of the sorts of things you can define with Qsymm:

# Identity in 3D
E = qsymm.identity(3)
# Inversion in 3D
I = qsymm.inversion(3)
# 4-fold rotation around the x-axis
C4 = qsymm.rotation(1/4, [1, 0, 0])
# 3-fold rotation around the [1, 1, 1] axis
C3 = qsymm.rotation(1/3, [1, 1, 1])
# Time reversal
TR = qsymm.time_reversal(3)
# Particle-hole
PH = qsymm.particle_hole(3)

The documentation page of the qsymm.groups module contains an exhaustive list of what can be generated.

Alternatively, to define a point group element we initialize it directly by providing its real and orbital space action:

qsymm.PointGroupElement(R=np.diag([1, -1]), U='kron(sigma_y, sigma_z)')
\[M\left([ 0 -1]\right)\]

As with other Qsymm objects we can get a readable representation of these group elements:

C4
\[R\left(\frac{\pi}{2}, [1 0 0]\right)\]
TR
\[ \mathcal{T}\]

Given a set of group generators we can also generate a group:

cubic_gens = {I, C4, C3, TR, PH}
cubic_group = qsymm.groups.generate_group(cubic_gens)
print(len(cubic_group))
192

Group elements can be multiplied and inverted, as we would expect:

C3 * C4
\[R\left(\pi, [1 1 0]\right)\]
C3**-1
\[R\left(-\frac{2\pi}{3}, [1 1 1]\right)\]

We can also apply group elements to the Model that we defined in the previous section:

H_with_TR = TR.apply(H)
H_with_TR.tosympy(nsimplify=True)
\[\begin{split}\displaystyle \left[\begin{matrix}- \alpha k_{z} + \frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{y}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m} & - \alpha k_{x} - i \alpha k_{y}\\- \alpha k_{x} + i \alpha k_{y} & \alpha k_{z} + \frac{\hbar^{2} k_{x}^{2}}{2 m} + \frac{\hbar^{2} k_{y}^{2}}{2 m} + \frac{\hbar^{2} k_{z}^{2}}{2 m}\end{matrix}\right]\end{split}\]

Defining continuous group generators

In addition to the group elements we can also define generators of continuous groups using qsymm.groups.ContinuousGroupGenerator:

sz = qsymm.ContinuousGroupGenerator(None, np.array([[1, 0], [0, -1]]))

The first argument to ContinuousGroupGenerator is the realspace rotation generator; by specifying None we indicate that we want the rotation part to be zero. The second argument is the unitary action of the generator on the Hilbert space as a Hermitian matrix.

Applying a ContinuousGroupGenerator to a Model calculates the commutator:

sz.apply(H).tosympy(nsimplify=True)
\[\begin{split}\displaystyle \left[\begin{matrix}0 & - 2 i \alpha k_{x} - 2 \alpha k_{y}\\2 i \alpha k_{x} - 2 \alpha k_{y} & 0\end{matrix}\right]\end{split}\]

For the 3D Rashba Hamiltonian we defined at the start of the tutorial spin-z is not conserved, hence the commutator is non-zero.

Finding symmetries

The function symmetries allows us to find the symmetries of a Model. Let us find out whether the 3D Rashba Hamiltonian defined earlier has cubic group symmetry:

discrete_symm, continuous_symm = qsymm.symmetries(H, cubic_group)
print(len(discrete_symm), len(continuous_symm))
48 0

It has 48 discrete symmetries (cubic group without inversion and time-reversal) and no continuous symmetries (conserved quantities).

For more detailed examples see Finding Symmetries and Finding symmetries of the Kekule-Y continuum model.

Generating Hamiltonians from symmetry constraints

The qsymm.hamiltonian_generator module contains algorithms for generating Hamiltonians from symmetry constraints.

For example let us generate all 2-band \(k \cdot p\) Hamiltonians with the same discrete symmetries as the Rashba Hamiltonian that we found in the previous section:

family = qsymm.continuum_hamiltonian(discrete_symm, dim=3, total_power=2, prettify=True)
qsymm.display_family(family)
\[\begin{split}\displaystyle \left[\begin{matrix}1 & 0\\0 & 1\end{matrix}\right]\end{split}\]
\[\begin{split}\displaystyle \left[\begin{matrix}k_{z} & k_{x} - i k_{y}\\k_{x} + i k_{y} & - k_{z}\end{matrix}\right]\end{split}\]
\[\begin{split}\displaystyle \left[\begin{matrix}k_{x}^{2} + k_{y}^{2} + k_{z}^{2} & 0\\0 & k_{x}^{2} + k_{y}^{2} + k_{z}^{2}\end{matrix}\right]\end{split}\]

It is exactly the Hamiltonian family we started with.

For more detailed examples see Generating k \cdot p models, Generating tight-binding models and Finding symmetries of the Kekule-Y continuum model.

Saving and loading Qsymm models

Qsymm models and identified symmetries don’t guarantee consistent ordering and basis selection across multiple runs. To avoid irrerproducible results you may use the Model.tosympy method and serialize the resulting sympy expression as shown below.

To save we do:

H2D_sympy = H2D.tosympy()

with open("H2D.txt", "w") as f:
    f.write(repr(H2D_sympy))

To load we do:

with open("H2D.txt") as f:
    data = f.read()

loaded_H2D = qsymm.Model(
    sympy.parsing.sympy_parser.parse_expr(data),
    momenta=['k_x', 'k_z']
)