band-unfolding / banduppy

Python version ofthe BandUP code
GNU General Public License v3.0
17 stars 7 forks source link

Wrong supercell to primitive cell k-path transformation in the final stage #28

Closed bmondal94 closed 2 months ago

bmondal94 commented 2 months ago

The matrix multiplication in the _generate_kpoints_line() function should be in the opposite order. Instead of

RecLattice = 2*np.pi*np.dot(transformation_matrix, np.linalg.pinv(sc_lattice)).T
which is equivalent to
RecLattice = 2*np.pi*np.linalg.pinv(sc_lattice).T @ transformation_matrix.T

it should be RecLattice = transformation_matrix.T @ 2*np.pi*np.linalg.pinv(sc_lattice).T

bmondal94 commented 2 months ago

Following the matrix commutativity rule; if the primitive-supercell transformation matrix is diagonal, there will be no problem.

We will fix the problem in the next version. In the mean time, if your transformation matrix is non-diagonal you can use this short python scripts for temporary solution : [This script reads the already unfolding generated bandstructure_unfolded.dat and kpoints_unfolded.dat files, and re-write back the k-points path (only the Angstrom^-1 data which were wrong) after some post-process into those files.]

from ase.io import read
import numpy as np
import banduppy
import os

'''
This script reads the already unfolding generated bandstructure_unfolded.dat and kpoints_unfolded.dat files, 
and re-write the k-points path (only the Angstrom^-1 points) after some post-process.
'''
##################################################################################
# The folder path where you have POSCAR, bandstructure_unfolded.dat, and kpoints_unfolded.dat file
bandstr_data_file_path = <put_your_path_here>

transformation_matrix = <your_transformation_matrix> # e.g.: np.array([[-1, 1, 1], [1, -1, 1], [1, 1, -1]])

breakTHRESH = 0.5 #<== this is the same option that you used during unfolding

##################################################################################
# Reads old unfolded kpoints and bandstructure files.
pc_kpts_coord = np.loadtxt(os.path.join(bandstr_data_file_path, 'kpoints_unfolded.dat'))[:, 2:]
unfolded_bandstructure_ = np.loadtxt(os.path.join(bandstr_data_file_path, 'bandstructure_unfolded.dat'))
# If you have ASE, then perfect. Else just read 3 lattice vextors in a numpy array
# Supercell lattice vectors in real space
real_lattice_sc = read(f'{bandstr_data_file_path}/CONTCAR', format='vasp').cell

##################################################################################
# Reciprocal lattice super cell
reciprocal_lattice_sc = 2*np.pi*np.linalg.pinv(real_lattice_sc).T
# Reciprocal lattice primitive cell
reciprocal_lattice_pc = transformation_matrix.T @ reciprocal_lattice_sc
# Convert k-points to inverse angstrom unit
KPcart = np.dot(pc_kpts_coord, reciprocal_lattice_pc)
# Generate distance array
K = np.zeros(KPcart.shape[0])
k = np.linalg.norm(KPcart[1:, :] - KPcart[:-1, :], axis=1)
k[k > breakTHRESH] = 0.0
K[1:] = np.cumsum(k)

##################################################################################
unfolded_kpts_dat = np.array([[ii]+[K[ii]]+list(iii) for ii, iii in enumerate(pc_kpts_coord)])
unfolded_bandstructure = np.array([[kl[0], K[int(kl[0])], kl[2], kl[3]] for kl in unfolded_bandstructure_])
##################################################################################
# Re-write the files again and save
save_bandstr = banduppy.SaveBandStructuredata()
save_bandstr.save_unfolded_pc_kpts(unfolded_kpts_dat,
                                    save_dir=bandstr_data_file_path,
                                    file_name='kpoints_unfolded',
                                    file_name_suffix='',
                                    print_information='low')
save_bandstr.save_unfolded_bandstucture(unfolded_bandstructure,
                                        save_dir=bandstr_data_file_path,
                                        file_name='bandstructure_unfolded',
                                        file_name_suffix='',
                                        print_information='low',
                                        is_spinor=False)
##################################################################################
bmondal94 commented 2 months ago

Bug fixed in latest release v0.3.3