Hamiltonian set up

HQS Qolossal supports different ways of defining Hamiltonians: via abstract definition, via directly providing the Hamiltonian in matrix form or as a linear operator. Additionally, you can add disorder to your systems.

Hamiltonian as a structured input

Setting up a Hamiltonian can be done via structured input using the Lattice Builder. Please check out its documentation for more information. In short, to define a lattice system, you need to specify which atoms/sites make up the unit cell, which bonds connect them, and other important information regarding the structure. The corresponding input dictionary can be directly handed to HQS Qolossal for initialization, for example:

from qolossal import HamiltonianSparse

unitcell = {
    'atoms': [
        {'id': 0, 'name': 'A', 'position': [0.0, 0.0, 0.0], 'e0': 0.0}
    ],
    'lattice_vectors': [[1.0, 0.0, 0.0]],
    'bonds': [
        {'id_from': 0, 'id_to': 0, 'translation': [1, 0, 0], 't': -1},
    ],
}

system = {
    'system_size': [1000, 1, 1],
    'site_type': 'spinless',
    'system_boundary_condition': ['PBC', 'HW', 'HW']
}

lb_input = {
    'unitcell': unitcell,
    'system': system
}

H = HamiltonianSparse(lb_input)

The script above initializes a 1D chain of length 1000 as a sparse Hamiltonian.

Hamiltonian as a matrix

Setting up a Hamiltonian with its matrix form is as simple as expected: we simply need to initialize the corresponding Hamiltonian class with it. Let's make the example of a 1D chain, whose non-zero entries in the Hamiltonian are just the two off-diagonals with value \(-t\)

import numpy as np
from scipy.sparse import diags
from qolossal import HamiltonianSparseDirect

# parameters
t = 1  # hopping amplitude
dimension = 10  # dimension of the chain
lattice_vector_length = 0.5  # in Angstroms

# initialize the hamiltonian matrix.
# NOTE: we use the csr format as this is what qolossal uses
diagonals = -t * np.ones(dimension)
Hamiltonian_matrix = diags(
    diagonals=[diagonals, diagonals],
    offsets=[1, -1],
    shape=(dimension, dimension),
    format='csr'
)

# initialize HamiltonianSparseDirect.
H = HamiltonianSparseDirect(
    H = Hamiltonian_matrix,
    total_volume = dimension * lattice_vector_length
)

This will initialize our Hamiltonian object and we will be able to start computing the systems properties. The above code will result in a RuntimeWarning: as we have not provided the velocity operators, the calculation of quantities that require these, such as the DC conductivity, will not work. To do so, we follow the same rationale

sites_position = np.arange(dimension) * lattice_vector_length
position_operator = diags(
    [sites_position],
    offsets=[0],
    shape=(dimension, dimension),
    format='csr'
)

# compute the velocity as V_x=i [H, r_x] (hbar=1)
velocity_x = 1.j * (
    Hamiltonian_matrix @ position_operator - position_operator @ Hamiltonian_matrix
)
# the velocity matrices for the y and z directions are just zeros.
velocity_y = velocity_z = diags(
    np.zeros(dimension),
    shape=(dimension, dimension),
    format='csr'
)

# initialize HamiltonianSparseDirect.
H = HamiltonianSparseDirect(
    H = Hamiltonian_matrix,
    total_volume = dimension * lattice_vector_length,
    velocity_times_hbar_matrices=[velocity_x, velocity_y, velocity_z]
)

Now we can proceed to compute velocity-dependent quantities as well.

Hamiltonian as a linear operator

To initialize the Hamiltonian with a linear operator we proceed exactly like for the matrix version.

import numpy as np
from scipy.sparse import diags
from qolossal import HamiltonianLinearOperator

# parameters
t = 1  # hopping amplitude
dimension = 10  # dimension of the chain
lattice_vector_length = 0.5  # in Angstroms

# initialize the hamiltonian matrix.
# NOTE: we use the csr format as this is what qolossal uses
diagonals = -t * np.ones(dimension)
Hamiltonian_matrix = diags(
    diagonals=[diagonals, diagonals],
    offsets=[1, -1],
    shape=(dimension, dimension),
    format='csr'
)


# we define the linear operator starting from the matrix for the example, but this can be
# done however one likes
def linear_operator(v: np.ndarray, out: np.ndarray) -> np.ndarray:
    """Apply Hamiltonian to vector.

    Args:
        v (np.ndarray): input.
        out (np.ndarray): preallocated output.

    Returns:
        np.ndarray: H @ v
    """
    # syntax used to make sure that the result is written in the already allocated space of "out"
    out[:] = Hamiltonian_matrix @ v
    return out


H = HamiltonianLinearOperator(
    linear_operator,
    dim=dimension,
    dtype=float,
    total_volume=dimension * lattice_vector_length
)

For the velocity we proceed in exactly the same way, ensuring the linear operators have the same signature as the Hamiltonian linear operator.

Disorder

HQS Qolossal can also add disorder to the system. Currently only local disorder is implemented: in this case a random potential is applied to the different sites. To add local disorder to a system we can simply specify its properties during the Hamiltonian initialization:

from qolossal import HamiltonianSparseDIrect

# ... initialization of all relevant arguments

H = HamiltonianSparseDirect(
    H = Hamiltonian_matrix,
    total_volume = dimension * lattice_vector_length,
    disDict={"type": "onsite", "distr": "uniform", "sigma": 0.1, "avg": 0}
)

where the disorder dictionary specifies the type of disorder (currently only "onsite"), the distribution type of the random potential (allowed: uniform, gaussian, laplace) and average value and standard deviation of the random potential.