K-Space Hamiltonian and Band Structure
In this section, we will see how to generate the k-space Hamiltonian for a moiré system and compute its band structure. This is essential for understanding the electronic properties of the system, such as the presence of flat bands or Dirac cones.
1. Recap and Setup
To work in k-space, Periodic Boundary Conditions (pbc=True) are required, although not strictly enforced.
import numpy as np
import matplotlib.pyplot as plt
from moirepy import BilayerMoireLattice, HexagonalLayer
lattice = BilayerMoireLattice(
latticetype=HexagonalLayer,
ll1=1, ll2=2,
ul1=2, ul2=1,
n1=1, n2=1,
pbc=True # for k-space
)
# Output:
# twist angle = 0.3803 rad (21.7868 deg)
# 14 cells in upper lattice
# 14 cells in lower lattice
2. Single k-point Generation
To generate the k space hamiltonian, you first need to generate the real space hamiltonian.
lattice.generate_connections(inter_layer_radius=1)
real_ham = lattice.generate_hamiltonian(
tll=1, tuu=1,
tul=1, tlu=1,
tuself=1, tlself=1
).toarray().real
Note
we convert it to a dense array and take only the real part because the system size is small and the hamiltonian is purely real in this case. For larger systems, it is recommended to work with sparse matrices to save memory and computational resources. While the real-space hopping matrix is real, the k-space Hamiltonian is Hermitian and complex in general.
k = np.array([0.1, 0.2])
Hk = real_ham * lattice.get_phase(k).toarray()
print(f"Matrix shape: {Hk.shape}")
# Output:
# Matrix shape: (28, 28)
Manual Hermiticity
MoirePy gives you raw control over the hopping directions. To maintain a physical, Hermitian Hamiltonian, you must ensure your tlu (Lower to Upper) values are the complex conjugates of your tul (Upper to Lower) terms.
Phase Factors in k-space
We know that to go from real space to k-space, we need to multiply the hopping terms by a phase factor:
where \(H_{ij}\) is the real space hopping from site \(i\) to site \(j\). MoirePy keeps the \(R\) calculated when u run lattice.generate_connections and lattice.get_phase(k) will literally return the hamiltonian shaped grid of numbers (\(\cos\theta - i \sin\theta\)) where \(\theta = k \cdot R\).
However, this phase grid is not complete, it has 0s for all the places where hamiltonian is zero. To save memory and time, we don't calculate the phase for all pairs of sites.
3. Band Structure Calculation (The Physics Payoff)
The primary use of k-space is to visualize how energy levels evolve across the Moiré Brillouin Zone (MBZ).
Step 1: Define the Path
We sample the zone by interpolating between High-Symmetry points (\(\Gamma\), \(K\), \(M\)) in the MBZ. This path captures the essential physics of the system.
# Defined using the reciprocal lattice vectors of the moire system
G = np.array([0, 0])
K = (1/3) * (lattice.mlv1 + lattice.mlv2)
M = (1/2) * lattice.mlv1
def get_path(points, steps=50):
path = []
for i in range(len(points) - 1):
for t in np.linspace(0, 1, steps):
path.append(points[i] * (1 - t) + points[i + 1] * t)
return np.array(path)
k_path = get_path([G, K, M, G])
Step 2: Find eigenvalues of hamiltonian at each k-point
Since k-space Hamiltonians represent a single unit cell, they are much smaller than real-space matrices. We can use dense solvers like np.linalg.eigvalsh for high precision.
bands = []
for k in k_path:
Hk = real_ham * lattice.get_phase(k).toarray()
# Convert sparse k-matrix to dense for the eigensolver
eigvals = np.linalg.eigvalsh(Hk)
bands.append(eigvals)
bands = np.array(bands)
Step 3: Visualization
plt.figure(figsize=(6, 8))
# plt.plot(bands, color='black', linewidth=1)
plt.plot(bands, lw=0.8)
# Mark the high-symmetry transitions
plt.axvline(x=50, linestyle='--', color='gray', alpha=0.5)
plt.axvline(x=100, linestyle='--', color='gray', alpha=0.5)
plt.xticks([0, 50, 100, 150], [r'$\Gamma$', 'K', 'M', r'$\Gamma$'])
# plt.ylim(-4, 4)
plt.ylabel("Energy (t)")
plt.title("Moiré Band Structure")
plt.grid(axis='y', alpha=0.3)
plt.show()

Note
Using constant hopping (\(t=1\)) over a large supercell results in a dense, folded band structure. To see the isolated flat bands characteristic of the Magic Angle, refer to the Replicated Papers section where we implement distance-dependent (exponential) hopping.
Summary
- We constructed the k-space Hamiltonian from the same lattice definition
- We sampled momentum space along a high-symmetry path
- We computed and visualized the band structure
This is the standard workflow for studying electronic structure in periodic moiré systems.
Next Steps
- OBC vs PBC: Understand how boundary conditions affect the Hamiltonian.
- Defining Custom Layers: Go beyond built-in lattices.
- Designing Custom Hopping: Introduce realistic physics.
- Tutorials and Replicated Papers: See full physical results.