# Guan is an open-source python package developed and maintained by https://www.guanjihuan.com/about. The primary location of this package is on website https://py.guanjihuan.com.
# Modules
# # Module 1: basic functions
# # Module 2: Fourier transform
# # Module 3: Hamiltonian of finite size systems
# # Module 4: Hamiltonian of models in the reciprocal space
# # Module 5: band structures and wave functions
# # Module 6: Green functions
# # Module 7: density of states
# # Module 8: quantum transport
# # Module 9: topological invariant
# # Module 10: read and write
# # Module 11: plot figures
# # Module 12: data processing
# # Module 13: others
# import packages
import numpy as np
import math
import cmath
import copy
import guan
# Module 1: basic functions
## test
def test():
print('\nSuccess in the installation of Guan package!\n')
## Pauli matrices
def sigma_0():
return np.eye(2)
def sigma_x():
return np.array([[0, 1],[1, 0]])
def sigma_y():
return np.array([[0, -1j],[1j, 0]])
def sigma_z():
return np.array([[1, 0],[0, -1]])
## Kronecker product of Pauli matrices
def sigma_00():
return np.kron(sigma_0(), sigma_0())
def sigma_0x():
return np.kron(sigma_0(), sigma_x())
def sigma_0y():
return np.kron(sigma_0(), sigma_y())
def sigma_0z():
return np.kron(sigma_0(), sigma_z())
def sigma_x0():
return np.kron(sigma_x(), sigma_0())
def sigma_xx():
return np.kron(sigma_x(), sigma_x())
def sigma_xy():
return np.kron(sigma_x(), sigma_y())
def sigma_xz():
return np.kron(sigma_x(), sigma_z())
def sigma_y0():
return np.kron(sigma_y(), sigma_0())
def sigma_yx():
return np.kron(sigma_y(), sigma_x())
def sigma_yy():
return np.kron(sigma_y(), sigma_y())
def sigma_yz():
return np.kron(sigma_y(), sigma_z())
def sigma_z0():
return np.kron(sigma_z(), sigma_0())
def sigma_zx():
return np.kron(sigma_z(), sigma_x())
def sigma_zy():
return np.kron(sigma_z(), sigma_y())
def sigma_zz():
return np.kron(sigma_z(), sigma_z())
# Module 2: Fourier_transform
# Fourier transform for discrete lattices
def one_dimensional_fourier_transform(k, unit_cell, hopping):
unit_cell = np.array(unit_cell)
hopping = np.array(hopping)
hamiltonian = unit_cell+hopping*cmath.exp(1j*k)+hopping.transpose().conj()*cmath.exp(-1j*k)
return hamiltonian
def two_dimensional_fourier_transform_for_square_lattice(k1, k2, unit_cell, hopping_1, hopping_2):
unit_cell = np.array(unit_cell)
hopping_1 = np.array(hopping_1)
hopping_2 = np.array(hopping_2)
hamiltonian = unit_cell+hopping_1*cmath.exp(1j*k1)+hopping_1.transpose().conj()*cmath.exp(-1j*k1)+hopping_2*cmath.exp(1j*k2)+hopping_2.transpose().conj()*cmath.exp(-1j*k2)
return hamiltonian
def three_dimensional_fourier_transform_for_cubic_lattice(k1, k2, k3, unit_cell, hopping_1, hopping_2, hopping_3):
unit_cell = np.array(unit_cell)
hopping_1 = np.array(hopping_1)
hopping_2 = np.array(hopping_2)
hopping_3 = np.array(hopping_3)
hamiltonian = unit_cell+hopping_1*cmath.exp(1j*k1)+hopping_1.transpose().conj()*cmath.exp(-1j*k1)+hopping_2*cmath.exp(1j*k2)+hopping_2.transpose().conj()*cmath.exp(-1j*k2)+hopping_3*cmath.exp(1j*k3)+hopping_3.transpose().conj()*cmath.exp(-1j*k3)
return hamiltonian
def one_dimensional_fourier_transform_with_k(unit_cell, hopping):
import functools
hamiltonian_function = functools.partial(guan.one_dimensional_fourier_transform, unit_cell=unit_cell, hopping=hopping)
return hamiltonian_function
def two_dimensional_fourier_transform_for_square_lattice_with_k1_k2(unit_cell, hopping_1, hopping_2):
import functools
hamiltonian_function = functools.partial(guan.two_dimensional_fourier_transform_for_square_lattice, unit_cell=unit_cell, hopping_1=hopping_1, hopping_2=hopping_2)
return hamiltonian_function
def three_dimensional_fourier_transform_for_cubic_lattice_with_k1_k2_k3(unit_cell, hopping_1, hopping_2, hopping_3):
import functools
hamiltonian_function = functools.partial(guan.three_dimensional_fourier_transform_for_cubic_lattice, unit_cell=unit_cell, hopping_1=hopping_1, hopping_2=hopping_2, hopping_3=hopping_3)
return hamiltonian_function
## calculate reciprocal lattice vectors
def calculate_one_dimensional_reciprocal_lattice_vector(a1):
b1 = 2*np.pi/a1
return b1
def calculate_two_dimensional_reciprocal_lattice_vectors(a1, a2):
a1 = np.array(a1)
a2 = np.array(a2)
a1 = np.append(a1, 0)
a2 = np.append(a2, 0)
a3 = np.array([0, 0, 1])
b1 = 2*np.pi*np.cross(a2, a3)/np.dot(a1, np.cross(a2, a3))
b2 = 2*np.pi*np.cross(a3, a1)/np.dot(a1, np.cross(a2, a3))
b1 = np.delete(b1, 2)
b2 = np.delete(b2, 2)
return b1, b2
def calculate_three_dimensional_reciprocal_lattice_vectors(a1, a2, a3):
a1 = np.array(a1)
a2 = np.array(a2)
a3 = np.array(a3)
b1 = 2*np.pi*np.cross(a2, a3)/np.dot(a1, np.cross(a2, a3))
b2 = 2*np.pi*np.cross(a3, a1)/np.dot(a1, np.cross(a2, a3))
b3 = 2*np.pi*np.cross(a1, a2)/np.dot(a1, np.cross(a2, a3))
return b1, b2, b3
def calculate_one_dimensional_reciprocal_lattice_vector_with_sympy(a1):
import sympy
b1 = 2*sympy.pi/a1
return b1
def calculate_two_dimensional_reciprocal_lattice_vectors_with_sympy(a1, a2):
import sympy
a1 = sympy.Matrix(1, 3, [a1[0], a1[1], 0])
a2 = sympy.Matrix(1, 3, [a2[0], a2[1], 0])
a3 = sympy.Matrix(1, 3, [0, 0, 1])
cross_a2_a3 = a2.cross(a3)
cross_a3_a1 = a3.cross(a1)
b1 = 2*sympy.pi*cross_a2_a3/a1.dot(cross_a2_a3)
b2 = 2*sympy.pi*cross_a3_a1/a1.dot(cross_a2_a3)
b1 = sympy.Matrix(1, 2, [b1[0], b1[1]])
b2 = sympy.Matrix(1, 2, [b2[0], b2[1]])
return b1, b2
def calculate_three_dimensional_reciprocal_lattice_vectors_with_sympy(a1, a2, a3):
import sympy
cross_a2_a3 = a2.cross(a3)
cross_a3_a1 = a3.cross(a1)
cross_a1_a2 = a1.cross(a2)
b1 = 2*sympy.pi*cross_a2_a3/a1.dot(cross_a2_a3)
b2 = 2*sympy.pi*cross_a3_a1/a1.dot(cross_a2_a3)
b3 = 2*sympy.pi*cross_a1_a2/a1.dot(cross_a2_a3)
return b1, b2, b3
# Module 3: Hamiltonian of finite size systems
def hamiltonian_of_finite_size_system_along_one_direction(N, on_site=0, hopping=1, period=0):
on_site = np.array(on_site)
hopping = np.array(hopping)
if on_site.shape==():
dim = 1
else:
dim = on_site.shape[0]
hamiltonian = np.zeros((N*dim, N*dim), dtype=complex)
for i0 in range(N):
hamiltonian[i0*dim+0:i0*dim+dim, i0*dim+0:i0*dim+dim] = on_site
for i0 in range(N-1):
hamiltonian[i0*dim+0:i0*dim+dim, (i0+1)*dim+0:(i0+1)*dim+dim] = hopping
hamiltonian[(i0+1)*dim+0:(i0+1)*dim+dim, i0*dim+0:i0*dim+dim] = hopping.transpose().conj()
if period == 1:
hamiltonian[(N-1)*dim+0:(N-1)*dim+dim, 0:dim] = hopping
hamiltonian[0:dim, (N-1)*dim+0:(N-1)*dim+dim] = hopping.transpose().conj()
return hamiltonian
def hamiltonian_of_finite_size_system_along_two_directions_for_square_lattice(N1, N2, on_site=0, hopping_1=1, hopping_2=1, period_1=0, period_2=0):
on_site = np.array(on_site)
hopping_1 = np.array(hopping_1)
hopping_2 = np.array(hopping_2)
if on_site.shape==():
dim = 1
else:
dim = on_site.shape[0]
hamiltonian = np.zeros((N1*N2*dim, N1*N2*dim), dtype=complex)
for i1 in range(N1):
for i2 in range(N2):
hamiltonian[i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim, i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim] = on_site
for i1 in range(N1-1):
for i2 in range(N2):
hamiltonian[i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim, (i1+1)*N2*dim+i2*dim+0:(i1+1)*N2*dim+i2*dim+dim] = hopping_1
hamiltonian[(i1+1)*N2*dim+i2*dim+0:(i1+1)*N2*dim+i2*dim+dim, i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim] = hopping_1.transpose().conj()
for i1 in range(N1):
for i2 in range(N2-1):
hamiltonian[i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim, i1*N2*dim+(i2+1)*dim+0:i1*N2*dim+(i2+1)*dim+dim] = hopping_2
hamiltonian[i1*N2*dim+(i2+1)*dim+0:i1*N2*dim+(i2+1)*dim+dim, i1*N2*dim+i2*dim+0:i1*N2*dim+i2*dim+dim] = hopping_2.transpose().conj()
if period_1 == 1:
for i2 in range(N2):
hamiltonian[(N1-1)*N2*dim+i2*dim+0:(N1-1)*N2*dim+i2*dim+dim, i2*dim+0:i2*dim+dim] = hopping_1
hamiltonian[i2*dim+0:i2*dim+dim, (N1-1)*N2*dim+i2*dim+0:(N1-1)*N2*dim+i2*dim+dim] = hopping_1.transpose().conj()
if period_2 == 1:
for i1 in range(N1):
hamiltonian[i1*N2*dim+(N2-1)*dim+0:i1*N2*dim+(N2-1)*dim+dim, i1*N2*dim+0:i1*N2*dim+dim] = hopping_2
hamiltonian[i1*N2*dim+0:i1*N2*dim+dim, i1*N2*dim+(N2-1)*dim+0:i1*N2*dim+(N2-1)*dim+dim] = hopping_2.transpose().conj()
return hamiltonian
def hamiltonian_of_finite_size_system_along_three_directions_for_cubic_lattice(N1, N2, N3, on_site=0, hopping_1=1, hopping_2=1, hopping_3=1, period_1=0, period_2=0, period_3=0):
on_site = np.array(on_site)
hopping_1 = np.array(hopping_1)
hopping_2 = np.array(hopping_2)
hopping_3 = np.array(hopping_3)
if on_site.shape==():
dim = 1
else:
dim = on_site.shape[0]
hamiltonian = np.zeros((N1*N2*N3*dim, N1*N2*N3*dim), dtype=complex)
for i1 in range(N1):
for i2 in range(N2):
for i3 in range(N3):
hamiltonian[i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim, i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim] = on_site
for i1 in range(N1-1):
for i2 in range(N2):
for i3 in range(N3):
hamiltonian[i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim, (i1+1)*N2*N3*dim+i2*N3*dim+i3*dim+0:(i1+1)*N2*N3*dim+i2*N3*dim+i3*dim+dim] = hopping_1
hamiltonian[(i1+1)*N2*N3*dim+i2*N3*dim+i3*dim+0:(i1+1)*N2*N3*dim+i2*N3*dim+i3*dim+dim, i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim] = hopping_1.transpose().conj()
for i1 in range(N1):
for i2 in range(N2-1):
for i3 in range(N3):
hamiltonian[i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim, i1*N2*N3*dim+(i2+1)*N3*dim+i3*dim+0:i1*N2*N3*dim+(i2+1)*N3*dim+i3*dim+dim] = hopping_2
hamiltonian[i1*N2*N3*dim+(i2+1)*N3*dim+i3*dim+0:i1*N2*N3*dim+(i2+1)*N3*dim+i3*dim+dim, i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim] = hopping_2.transpose().conj()
for i1 in range(N1):
for i2 in range(N2):
for i3 in range(N3-1):
hamiltonian[i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim, i1*N2*N3*dim+i2*N3*dim+(i3+1)*dim+0:i1*N2*N3*dim+i2*N3*dim+(i3+1)*dim+dim] = hopping_3
hamiltonian[i1*N2*N3*dim+i2*N3*dim+(i3+1)*dim+0:i1*N2*N3*dim+i2*N3*dim+(i3+1)*dim+dim, i1*N2*N3*dim+i2*N3*dim+i3*dim+0:i1*N2*N3*dim+i2*N3*dim+i3*dim+dim] = hopping_3.transpose().conj()
if period_1 == 1:
for i2 in range(N2):
for i3 in range(N3):
hamiltonian[(N1-1)*N2*N3*dim+i2*N3*dim+i3*dim+0:(N1-1)*N2*N3*dim+i2*N3*dim+i3*dim+dim, i2*N3*dim+i3*dim+0:i2*N3*dim+i3*dim+dim] = hopping_1
hamiltonian[i2*N3*dim+i3*dim+0:i2*N3*dim+i3*dim+dim, (N1-1)*N2*N3*dim+i2*N3*dim+i3*dim+0:(N1-1)*N2*N3*dim+i2*N3*dim+i3*dim+dim] = hopping_1.transpose().conj()
if period_2 == 1:
for i1 in range(N1):
for i3 in range(N3):
hamiltonian[i1*N2*N3*dim+(N2-1)*N3*dim+i3*dim+0:i1*N2*N3*dim+(N2-1)*N3*dim+i3*dim+dim, i1*N2*N3*dim+i3*dim+0:i1*N2*N3*dim+i3*dim+dim] = hopping_2
hamiltonian[i1*N2*N3*dim+i3*dim+0:i1*N2*N3*dim+i3*dim+dim, i1*N2*N3*dim+(N2-1)*N3*dim+i3*dim+0:i1*N2*N3*dim+(N2-1)*N3*dim+i3*dim+dim] = hopping_2.transpose().conj()
if period_3 == 1:
for i1 in range(N1):
for i2 in range(N2):
hamiltonian[i1*N2*N3*dim+i2*N3*dim+(N3-1)*dim+0:i1*N2*N3*dim+i2*N3*dim+(N3-1)*dim+dim, i1*N2*N3*dim+i2*N3*dim+0:i1*N2*N3*dim+i2*N3*dim+dim] = hopping_3
hamiltonian[i1*N2*N3*dim+i2*N3*dim+0:i1*N2*N3*dim+i2*N3*dim+dim, i1*N2*N3*dim+i2*N3*dim+(N3-1)*dim+0:i1*N2*N3*dim+i2*N3*dim+(N3-1)*dim+dim] = hopping_3.transpose().conj()
return hamiltonian
def hopping_matrix_along_zigzag_direction_for_graphene_ribbon(N):
hopping = np.zeros((4*N, 4*N), dtype=complex)
for i0 in range(N):
hopping[4*i0+1, 4*i0+0] = 1
hopping[4*i0+2, 4*i0+3] = 1
return hopping
def hamiltonian_of_finite_size_system_along_two_directions_for_graphene(N1, N2, period_1=0, period_2=0):
on_site = guan.hamiltonian_of_finite_size_system_along_one_direction(4)
hopping_1 = guan.hopping_matrix_along_zigzag_direction_for_graphene_ribbon(1)
hopping_2 = np.zeros((4, 4), dtype=complex)
hopping_2[3, 0] = 1
hamiltonian = guan.finite_size_along_two_directions_for_square_lattice(N1, N2, on_site, hopping_1, hopping_2, period_1, period_2)
return hamiltonian
# Module 4: Hamiltonian of models in the reciprocal space
def hamiltonian_of_simple_chain(k):
hamiltonian = guan.one_dimensional_fourier_transform(k, unit_cell=0, hopping=1)
return hamiltonian
def hamiltonian_of_square_lattice(k1, k2):
hamiltonian = guan.two_dimensional_fourier_transform_for_square_lattice(k1, k2, unit_cell=0, hopping_1=1, hopping_2=1)
return hamiltonian
def hamiltonian_of_square_lattice_in_quasi_one_dimension(k, N=10):
h00 = np.zeros((N, N), dtype=complex) # hopping in a unit cell
h01 = np.zeros((N, N), dtype=complex) # hopping between unit cells
for i in range(N-1):
h00[i, i+1] = 1
h00[i+1, i] = 1
for i in range(N):
h01[i, i] = 1
hamiltonian = guan.one_dimensional_fourier_transform(k, unit_cell=h00, hopping=h01)
return hamiltonian
def hamiltonian_of_cubic_lattice(k1, k2, k3):
hamiltonian = guan.three_dimensional_fourier_transform_for_cubic_lattice(k1, k2, k3, unit_cell=0, hopping_1=1, hopping_2=1, hopping_3=1)
return hamiltonian
def hamiltonian_of_ssh_model(k, v=0.6, w=1):
hamiltonian = np.zeros((2, 2), dtype=complex)
hamiltonian[0,1] = v+w*cmath.exp(-1j*k)
hamiltonian[1,0] = v+w*cmath.exp(1j*k)
return hamiltonian
def hamiltonian_of_graphene(k1, k2, M=0, t=1, a=1/math.sqrt(3)):
h0 = np.zeros((2, 2), dtype=complex) # mass term
h1 = np.zeros((2, 2), dtype=complex) # nearest hopping
h0[0, 0] = M
h0[1, 1] = -M
h1[1, 0] = t*(cmath.exp(1j*k2*a)+cmath.exp(1j*math.sqrt(3)/2*k1*a-1j/2*k2*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a-1j/2*k2*a))
h1[0, 1] = h1[1, 0].conj()
hamiltonian = h0 + h1
return hamiltonian
def hamiltonian_of_graphene_with_zigzag_in_quasi_one_dimension(k, N=10, M=0, t=1):
h00 = np.zeros((4*N, 4*N), dtype=complex) # hopping in a unit cell
h01 = np.zeros((4*N, 4*N), dtype=complex) # hopping between unit cells
for i in range(N):
h00[i*4+0, i*4+0] = M
h00[i*4+1, i*4+1] = -M
h00[i*4+2, i*4+2] = M
h00[i*4+3, i*4+3] = -M
h00[i*4+0, i*4+1] = t
h00[i*4+1, i*4+0] = t
h00[i*4+1, i*4+2] = t
h00[i*4+2, i*4+1] = t
h00[i*4+2, i*4+3] = t
h00[i*4+3, i*4+2] = t
for i in range(N-1):
h00[i*4+3, (i+1)*4+0] = t
h00[(i+1)*4+0, i*4+3] = t
for i in range(N):
h01[i*4+1, i*4+0] = t
h01[i*4+2, i*4+3] = t
hamiltonian = guan.one_dimensional_fourier_transform(k, unit_cell=h00, hopping=h01)
return hamiltonian
def hamiltonian_of_haldane_model(k1, k2, M=2/3, t1=1, t2=1/3, phi=math.pi/4, a=1/math.sqrt(3)):
h0 = np.zeros((2, 2), dtype=complex) # mass term
h1 = np.zeros((2, 2), dtype=complex) # nearest hopping
h2 = np.zeros((2, 2), dtype=complex) # next nearest hopping
h0[0, 0] = M
h0[1, 1] = -M
h1[1, 0] = t1*(cmath.exp(1j*k2*a)+cmath.exp(1j*math.sqrt(3)/2*k1*a-1j/2*k2*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a-1j/2*k2*a))
h1[0, 1] = h1[1, 0].conj()
h2[0, 0] = t2*cmath.exp(-1j*phi)*(cmath.exp(1j*math.sqrt(3)*k1*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a+1j*3/2*k2*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a-1j*3/2*k2*a))
h2[1, 1] = t2*cmath.exp(1j*phi)*(cmath.exp(1j*math.sqrt(3)*k1*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a+1j*3/2*k2*a)+cmath.exp(-1j*math.sqrt(3)/2*k1*a-1j*3/2*k2*a))
hamiltonian = h0 + h1 + h2 + h2.transpose().conj()
return hamiltonian
def hamiltonian_of_haldane_model_in_quasi_one_dimension(k, N=10, M=2/3, t1=1, t2=1/3, phi=math.pi/4):
h00 = np.zeros((4*N, 4*N), dtype=complex) # hopping in a unit cell
h01 = np.zeros((4*N, 4*N), dtype=complex) # hopping between unit cells
for i in range(N):
h00[i*4+0, i*4+0] = M
h00[i*4+1, i*4+1] = -M
h00[i*4+2, i*4+2] = M
h00[i*4+3, i*4+3] = -M
h00[i*4+0, i*4+1] = t1
h00[i*4+1, i*4+0] = t1
h00[i*4+1, i*4+2] = t1
h00[i*4+2, i*4+1] = t1
h00[i*4+2, i*4+3] = t1
h00[i*4+3, i*4+2] = t1
h00[i*4+0, i*4+2] = t2*cmath.exp(-1j*phi)
h00[i*4+2, i*4+0] = h00[i*4+0, i*4+2].conj()
h00[i*4+1, i*4+3] = t2*cmath.exp(-1j*phi)
h00[i*4+3, i*4+1] = h00[i*4+1, i*4+3].conj()
for i in range(N-1):
h00[i*4+3, (i+1)*4+0] = t1
h00[(i+1)*4+0, i*4+3] = t1
h00[i*4+2, (i+1)*4+0] = t2*cmath.exp(1j*phi)
h00[(i+1)*4+0, i*4+2] = h00[i*4+2, (i+1)*4+0].conj()
h00[i*4+3, (i+1)*4+1] = t2*cmath.exp(1j*phi)
h00[(i+1)*4+1, i*4+3] = h00[i*4+3, (i+1)*4+1].conj()
for i in range(N):
h01[i*4+1, i*4+0] = t1
h01[i*4+2, i*4+3] = t1
h01[i*4+0, i*4+0] = t2*cmath.exp(1j*phi)
h01[i*4+1, i*4+1] = t2*cmath.exp(-1j*phi)
h01[i*4+2, i*4+2] = t2*cmath.exp(1j*phi)
h01[i*4+3, i*4+3] = t2*cmath.exp(-1j*phi)
h01[i*4+1, i*4+3] = t2*cmath.exp(1j*phi)
h01[i*4+2, i*4+0] = t2*cmath.exp(-1j*phi)
if i != 0:
h01[i*4+1, (i-1)*4+3] = t2*cmath.exp(1j*phi)
for i in range(N-1):
h01[i*4+2, (i+1)*4+0] = t2*cmath.exp(-1j*phi)
hamiltonian = h00 + h01*cmath.exp(1j*k) + h01.transpose().conj()*cmath.exp(-1j*k)
return hamiltonian
def hamiltonian_of_one_QAH_model(k1, k2, t1=1, t2=1, t3=0.5, m=-1):
hamiltonian = np.zeros((2, 2), dtype=complex)
hamiltonian[0, 1] = 2*t1*math.cos(k1)-1j*2*t1*math.cos(k2)
hamiltonian[1, 0] = 2*t1*math.cos(k1)+1j*2*t1*math.cos(k2)
hamiltonian[0, 0] = m+2*t3*math.sin(k1)+2*t3*math.sin(k2)+2*t2*math.cos(k1+k2)
hamiltonian[1, 1] = -(m+2*t3*math.sin(k1)+2*t3*math.sin(k2)+2*t2*math.cos(k1+k2))
return hamiltonian
def hamiltonian_of_BBH_model(kx, ky, gamma_x=0.5, gamma_y=0.5, lambda_x=1, lambda_y=1):
# label of atoms in a unit cell
# (2) —— (0)
# | |
# (1) —— (3)
hamiltonian = np.zeros((4, 4), dtype=complex)
hamiltonian[0, 2] = gamma_x+lambda_x*cmath.exp(1j*kx)
hamiltonian[1, 3] = gamma_x+lambda_x*cmath.exp(-1j*kx)
hamiltonian[0, 3] = gamma_y+lambda_y*cmath.exp(1j*ky)
hamiltonian[1, 2] = -gamma_y-lambda_y*cmath.exp(-1j*ky)
hamiltonian[2, 0] = np.conj(hamiltonian[0, 2])
hamiltonian[3, 1] = np.conj(hamiltonian[1, 3])
hamiltonian[3, 0] = np.conj(hamiltonian[0, 3])
hamiltonian[2, 1] = np.conj(hamiltonian[1, 2])
return hamiltonian
# Module 5: band structures and wave functions
## band structures
def calculate_eigenvalue(hamiltonian):
if np.array(hamiltonian).shape==():
eigenvalue = np.real(hamiltonian)
else:
eigenvalue, eigenvector = np.linalg.eigh(hamiltonian)
return eigenvalue
def calculate_eigenvalue_with_one_parameter(x_array, hamiltonian_function, print_show=0):
dim_x = np.array(x_array).shape[0]
i0 = 0
if np.array(hamiltonian_function(0)).shape==():
eigenvalue_array = np.zeros((dim_x, 1))
for x0 in x_array:
hamiltonian = hamiltonian_function(x0)
eigenvalue_array[i0, 0] = np.real(hamiltonian)
i0 += 1
else:
dim = np.array(hamiltonian_function(0)).shape[0]
eigenvalue_array = np.zeros((dim_x, dim))
for x0 in x_array:
if print_show==1:
print(x0)
hamiltonian = hamiltonian_function(x0)
eigenvalue, eigenvector = np.linalg.eigh(hamiltonian)
eigenvalue_array[i0, :] = eigenvalue
i0 += 1
return eigenvalue_array
def calculate_eigenvalue_with_two_parameters(x_array, y_array, hamiltonian_function, print_show=0, print_show_more=0):
dim_x = np.array(x_array).shape[0]
dim_y = np.array(y_array).shape[0]
if np.array(hamiltonian_function(0,0)).shape==():
eigenvalue_array = np.zeros((dim_y, dim_x, 1))
i0 = 0
for y0 in y_array:
j0 = 0
for x0 in x_array:
hamiltonian = hamiltonian_function(x0, y0)
eigenvalue_array[i0, j0, 0] = np.real(hamiltonian)
j0 += 1
i0 += 1
else:
dim = np.array(hamiltonian_function(0, 0)).shape[0]
eigenvalue_array = np.zeros((dim_y, dim_x, dim))
i0 = 0
for y0 in y_array:
j0 = 0
if print_show==1:
print(y0)
for x0 in x_array:
if print_show_more==1:
print(x0)
hamiltonian = hamiltonian_function(x0, y0)
eigenvalue, eigenvector = np.linalg.eigh(hamiltonian)
eigenvalue_array[i0, j0, :] = eigenvalue
j0 += 1
i0 += 1
return eigenvalue_array
## wave functions
def calculate_eigenvector(hamiltonian):
eigenvalue, eigenvector = np.linalg.eigh(hamiltonian)
return eigenvector
## find vector with the same gauge
def find_vector_with_the_same_gauge_with_binary_search(vector_target, vector_ref, show_error=1, show_times=0, show_phase=0, n_test=1000, precision=1e-6):
phase_1_pre = 0
phase_2_pre = np.pi
for i0 in range(n_test):
test_1 = np.sum(np.abs(vector_target*cmath.exp(1j*phase_1_pre) - vector_ref))
test_2 = np.sum(np.abs(vector_target*cmath.exp(1j*phase_2_pre) - vector_ref))
if test_1 < precision:
phase = phase_1_pre
if show_times==1:
print('Binary search times=', i0)
break
if i0 == n_test-1:
phase = phase_1_pre
if show_error==1:
print('Gauge not found with binary search times=', i0)
if test_1 < test_2:
if i0 == 0:
phase_1 = phase_1_pre-(phase_2_pre-phase_1_pre)/2
phase_2 = phase_1_pre+(phase_2_pre-phase_1_pre)/2
else:
phase_1 = phase_1_pre
phase_2 = phase_1_pre+(phase_2_pre-phase_1_pre)/2
else:
if i0 == 0:
phase_1 = phase_2_pre-(phase_2_pre-phase_1_pre)/2
phase_2 = phase_2_pre+(phase_2_pre-phase_1_pre)/2
else:
phase_1 = phase_2_pre-(phase_2_pre-phase_1_pre)/2
phase_2 = phase_2_pre
phase_1_pre = phase_1
phase_2_pre = phase_2
vector_target = vector_target*cmath.exp(1j*phase)
if show_phase==1:
print('Phase=', phase)
return vector_target
def find_vector_with_fixed_gauge_by_making_one_component_real(vector, precision=0.005, index=None):
vector = np.array(vector)
if index == None:
index = np.argmax(np.abs(vector))
sign_pre = np.sign(np.imag(vector[index]))
for phase in np.arange(0, 2*np.pi, precision):
sign = np.sign(np.imag(vector[index]*cmath.exp(1j*phase)))
if np.abs(np.imag(vector[index]*cmath.exp(1j*phase))) < 1e-9 or sign == -sign_pre:
break
sign_pre = sign
vector = vector*cmath.exp(1j*phase)
if np.real(vector[index]) < 0:
vector = -vector
return vector
def find_vector_array_with_fixed_gauge_by_making_one_component_real(vector_array, precision=0.005):
vector_sum = 0
Num_k = np.array(vector_array).shape[0]
for i0 in range(Num_k):
vector_sum += np.abs(vector_array[i0])
index = np.argmax(np.abs(vector_sum))
for i0 in range(Num_k):
vector_array[i0] = guan.find_vector_with_fixed_gauge_by_making_one_component_real(vector_array[i0], precision=precision, index=index)
return vector_array
def rotation_of_degenerate_vectors(vector1, vector2, index1, index2, precision=0.01, criterion=0.01, show_theta=0):
vector1 = np.array(vector1)
vector2 = np.array(vector2)
if np.abs(vector1[index2])>criterion or np.abs(vector2[index1])>criterion:
for theta in np.arange(0, 2*math.pi, precision):
if show_theta==1:
print(theta)
for phi1 in np.arange(0, 2*math.pi, precision):
for phi2 in np.arange(0, 2*math.pi, precision):
vector1_test = cmath.exp(1j*phi1)*vector1*math.cos(theta)+cmath.exp(1j*phi2)*vector2*math.sin(theta)
vector2_test = -cmath.exp(-1j*phi2)*vector1*math.sin(theta)+cmath.exp(-1j*phi1)*vector2*math.cos(theta)
if np.abs(vector1_test[index2]).*?
.*?
', content, re.S)[0][3:-4] if show_translation==1: time.sleep(translation_time) print(translation) time.sleep(rest_time) pygame.mixer.music.stop() print() def play_element_words(random_on=0, show_translation=1, show_link=1, translation_time=2, rest_time=1): from bs4 import BeautifulSoup import re import urllib.request import requests import os import pygame import time import ssl import random ssl._create_default_https_context = ssl._create_unverified_context html = urllib.request.urlopen("https://www.guanjihuan.com/archives/10897").read().decode('utf-8') directory = 'prons/' exist_directory = os.path.exists(directory) html_file = urllib.request.urlopen("https://file.guanjihuan.com/words/periodic_table_of_elements/"+directory).read().decode('utf-8') if exist_directory == 0: os.makedirs(directory) soup = BeautifulSoup(html, features='lxml') contents = re.findall('.*?
', content, re.S)[0][3:-4] if show_translation==1: time.sleep(translation_time) print(translation) time.sleep(rest_time) pygame.mixer.music.stop() print()