Skip to content

Latest commit

 

History

History
301 lines (206 loc) · 8.47 KB

kdotp_generator.rst

File metadata and controls

301 lines (206 loc) · 8.47 KB

Generating k \cdot p models

.. seealso::
    The complete source code of this example can be found in :jupyter-download:script:`kdotp_generator`.
    A Jupyter notebook can be found in :jupyter-download:notebook:`kdotp_generator`.

.. jupyter-kernel::
    :id: kdotp_generator

.. jupyter-execute::
    :hide-code:

    import numpy as np
    import sympy

    import qsymm

In the :ref:`Bloch generator tutorial <tutorial_bloch_generator>` we saw how to generate tight-binding Hamiltonians using Qsymm. In this tutorial we will explore how to generate k \cdot p Hamiltonians.

Hexagonal warping

We reproduce the effective Hamiltonian derived in Phys. Rev. Lett. 103, 266801 (2009) to explain hexagonal warping effects in the surface dispersion of a topological insulator.

The symmetry group is generated by threefold rotation symmetry, mirror symmetry and time-reversal symmetry.

.. jupyter-execute::

    # C3 rotational symmetry - invariant under 2pi/3
    C3 = qsymm.rotation(1/3, spin=1/2)

    # Time reversal
    TR = qsymm.time_reversal(2, spin=1/2)

    # Mirror symmetry in x
    Mx = qsymm.mirror([1, 0], spin=1/2)

    symmetries = [C3, TR, Mx]

The Hamiltonian includes the momenta k_x and k_y, and contains terms up to third order in momenta.

.. jupyter-execute::

    dim = 2  # Momenta along x and y
    total_power = 3  # Maximum power of momenta

We use ~qsymm.hamiltonian_generator.continuum_hamiltonian to generate the Hamiltonian family:

.. jupyter-execute::

    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)

We can then verify that the Hamiltonian family satisfies the symmetries:

.. jupyter-execute::

    qsymm.hamiltonian_generator.check_symmetry(family, symmetries)


BHZ model

We reproduce the Hamiltonian for the quantum spin Hall effect derived in Science, 314, 1757 (2006).

The symmetry group is generated by spatial inversion symmetry, time-reversal symmetry and fourfold rotation symmetry.

.. jupyter-execute::

    # Spatial inversion
    pU = np.array([
        [1.0, 0.0, 0.0, 0.0],
        [0.0, -1.0, 0.0, 0.0],
        [0.0, 0.0, 1.0, 0.0],
        [0.0, 0.0, 0.0, -1.0],
    ])
    pS = qsymm.inversion(2, U=pU)

    # Time reversal
    trU = np.array([
        [0.0, 0.0, -1.0, 0.0],
        [0.0, 0.0, 0.0, -1.0],
        [1.0, 0.0, 0.0, 0.0],
        [0.0, 1.0, 0.0, 0.0],
    ])
    trS = qsymm.time_reversal(2, U=trU)

    # Rotation
    phi = 2.0 * np.pi / 4.0  # Impose 4-fold rotational symmetry
    rotU = np.array([
        [np.exp(-1j*phi/2), 0.0, 0.0, 0.0],
        [0.0, np.exp(-1j*3*phi/2), 0.0, 0.0],
        [0.0, 0.0, np.exp(1j*phi/2), 0.0],
        [0.0, 0.0, 0.0, np.exp(1j*3*phi/2)],
    ])
    rotS = qsymm.rotation(1/4, U=rotU)

    symmetries = [rotS, trS, pS]

The Hamiltonian includes the momenta k_x and k_y, with terms up to second order.

.. jupyter-execute::

    dim = 2
    total_power = 2

.. jupyter-execute::

    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)


Three-dimensional topological insulator

We reproduce the Hamiltonian for a three-dimensional topological insulator introduced in Nature Physics 5, 438–442 (2009).

The symmetry group is generated by threefold rotation symmetry around z, inversion symmetry and time-reversal symmetry.

.. jupyter-execute::

    # Spatial inversion
    pU = np.diag([1, -1, 1, -1])
    pS = qsymm.inversion(3, pU)

    # Time reversal
    trU = np.array([
        [0.0, 0.0, -1.0, 0.0],
        [0.0, 0.0, 0.0, -1.0],
        [1.0, 0.0, 0.0, 0.0],
        [0.0, 1.0, 0.0, 0.0],
    ])
    trS = qsymm.time_reversal(3, trU)

    # Rotation
    phi = 2.0 * np.pi / 3.0  # Impose 3-fold rotational symmetry
    rotU = np.array([
        [np.exp(1j*phi/2), 0.0, 0.0, 0.0],
        [0.0, np.exp(1j*phi/2), 0.0, 0.0],
        [0.0, 0.0, np.exp(-1j*phi/2), 0.0],
        [0.0, 0.0, 0.0, np.exp(-1j*phi/2)],
    ])
    rotS = qsymm.rotation(1/3, axis=[0, 0, 1], U=rotU)

    symmetries = [rotS, trS, pS]


The model includes the momenta k_x, k_y and k_z, up to second order.

.. jupyter-execute::

    dim = 3
    total_power = 2

.. jupyter-execute::

    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)


Continuous rotations

Rotation in real and spin space with spin 1/2. Should have linear Rashba-term as there is no inversion symmetry.

.. jupyter-execute::

    # 3D real space rotation generators
    L = qsymm.groups.L_matrices(3)
    # Spin-1/2 matrices
    S = qsymm.groups.spin_matrices(1/2)
    # Spin-3/2 spin matrices
    J = qsymm.groups.spin_matrices(3/2)
    # Continuous rotation generators
    symmetries = [qsymm.ContinuousGroupGenerator(l, s) for l, s in zip(L, S)]

    dim = 3
    total_power = 2

    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)

Also impose inversion symmetry, which removes the linear Rashba term.

.. jupyter-execute::

    symmetries = [qsymm.ContinuousGroupGenerator(l, s) for l, s in zip(L, S)]
    # Add inversion
    symmetries.append(qsymm.PointGroupElement(-np.eye(3), False, False, np.eye(2)))
    dim = 3
    total_power = 2
    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)

Rotations in real and spin space with spin 3/2.

.. jupyter-execute::

    symmetries = [qsymm.ContinuousGroupGenerator(l, s) for l, s in zip(L, J)]
    dim = 3
    total_power = 2
    family = qsymm.continuum_hamiltonian(symmetries, dim, total_power, prettify=True)
    qsymm.display_family(family)


Distorted SnTe

Reproduce the k \cdot p model for SnTe used in Phys. Rev. Lett. 122, 186801 (2019)

The symmetry group is generated by three-fold rotation symmetry, mirror symmetry in x, inversion symmetry, and time-reversal symmetry.

.. jupyter-execute::

    # Double spin-1/2 representation
    spin = [np.kron(s, np.eye(2)) for s in qsymm.groups.spin_matrices(1/2)]
    # C3 rotational symmetry
    C3 = qsymm.rotation(1/3, axis=[0, 0, 1], spin=spin)

    # Time reversal
    TR = qsymm.time_reversal(3, spin=spin)

    # Mirror x
    Mx = qsymm.mirror([1, 0, 0], spin=spin)

    # Inversion
    IU = np.kron(np.eye(2), np.diag([1, -1]))
    I = qsymm.inversion(3, IU)

.. jupyter-execute::

    dim = 3
    total_power = 2
    family = qsymm.continuum_hamiltonian([C3, TR, Mx, I], dim, total_power, prettify=True)
    qsymm.display_family(family)

There are 3 combinations proportional to the identity. Generate all terms proportional to the identity and subtract them.

.. jupyter-execute::

    identity_terms = qsymm.continuum_hamiltonian([qsymm.groups.identity(dim, 1)], dim, total_power)
    identity_terms = [
        qsymm.Model({
            key: np.kron(val, np.eye(4))
            for key, val in term.items()
        })
        for term in identity_terms
    ]

    family = qsymm.hamiltonian_generator.subtract_family(family, identity_terms, prettify=True)
    qsymm.display_family(family)


Break C3 symmetry

Removing the C3 symmetry, we find 8 additional terms, or 11 after removing terms proportional to the identity.

.. jupyter-execute::

    dim = 3
    total_power = 2
    no_c3_family = qsymm.continuum_hamiltonian([TR, Mx, I], dim, total_power, prettify=True)
    no_c3_family = qsymm.hamiltonian_generator.subtract_family(no_c3_family, identity_terms, prettify=True)

The new terms are given in the following.

.. jupyter-execute::

    new_terms = qsymm.hamiltonian_generator.subtract_family(no_c3_family, family, prettify=True)
    qsymm.display_family(new_terms)


Break inversion symmetry

Breaking inversion symmetry as well yields 22 additional terms.

.. jupyter-execute::

    dim = 3
    total_power = 2
    broken_family = qsymm.continuum_hamiltonian([TR, Mx], dim, total_power, prettify=True)
    broken_family = qsymm.hamiltonian_generator.subtract_family(broken_family, identity_terms, prettify=True)
    new_terms = qsymm.hamiltonian_generator.subtract_family(broken_family, no_c3_family, prettify=True)
    qsymm.display_family(new_terms)