3
0
mirror of https://github.com/triqs/dft_tools synced 2024-09-12 05:38:31 +02:00

plovasp: vasp_converter can now deal with H(k) and uncorrelated shells. Removed correleated shells and irrep from H(k) header.

This commit is contained in:
Malte Schüler 2019-07-01 10:51:33 +02:00
parent 3666518bdf
commit e9fd33dffb
2 changed files with 51 additions and 19 deletions

View File

@ -409,16 +409,14 @@ def hk_output(conf_pars, el_struct, pgroups):
Filename is defined by <basename> that is passed from config-file. Filename is defined by <basename> that is passed from config-file.
The Hk for each groups is stored in a '<basename>.hk<Ng>' file. The format is The Hk for each groups is stored in a '<basename>.hk<Ng>' file. The format is
as defined in the Hk dft_tools format similar as defined in the Hk dft_tools format, but does not store info
about correlated shells and irreps
nk # number of k-points nk # number of k-points
n_el # electron density n_el # electron density
n_sh # number of total atomic shells n_sh # number of total atomic shells
at sort l dim # atom, sort, l, dim at sort l dim # atom, sort, l, dim
at sort l dim # atom, sort, l, dim at sort l dim # atom, sort, l, dim
n_cr_sh # number of correlated shells
at sort l dim SO irep # atom, sort, l, dim, SO, irep
n_irrep dim_irrep # number of ireps, dim of irep
After these header lines, the file has to contain the Hamiltonian matrix After these header lines, the file has to contain the Hamiltonian matrix
in orbital space. The standard convention is that you give for each k-point in orbital space. The standard convention is that you give for each k-point
@ -433,7 +431,6 @@ def hk_output(conf_pars, el_struct, pgroups):
hk_fname = conf_pars.general['basename'] + '.hk%i'%(ig + 1) hk_fname = conf_pars.general['basename'] + '.hk%i'%(ig + 1)
print " Storing HK-group file '%s'..."%(hk_fname) print " Storing HK-group file '%s'..."%(hk_fname)
head_corr_shells = []
head_shells = [] head_shells = []
for ish in pgroup.ishells: for ish in pgroup.ishells:
@ -449,11 +446,9 @@ def hk_output(conf_pars, el_struct, pgroups):
# Derive sorts from equivalence classes # Derive sorts from equivalence classes
sh_dict['ion_list'] = ion_output sh_dict['ion_list'] = ion_output
sh_dict['ion_sort'] = shell.ion_sort sh_dict['ion_sort'] = shell.ion_sort
head_shells.append(sh_dict) head_shells.append(sh_dict)
if shell.corr:
head_corr_shells.append(sh_dict)
with open(hk_fname, 'wt') as f: with open(hk_fname, 'wt') as f:
@ -464,10 +459,7 @@ def hk_output(conf_pars, el_struct, pgroups):
f.write('%i # number of shells\n'%len(head_shells)) f.write('%i # number of shells\n'%len(head_shells))
for head in head_shells: for head in head_shells:
f.write('%i %i %i %i # atom sort l dim\n'%(head['ion_list'][0],head['ion_sort'][0],head['lorb'],head['ndim'])) f.write('%i %i %i %i # atom sort l dim\n'%(head['ion_list'][0],head['ion_sort'][0],head['lorb'],head['ndim']))
f.write('%i # number of correlated shells\n'%len(head_corr_shells))
for head in head_corr_shells:
f.write('%i %i %i %i 0 0 # atom sort l dim SO irrep\n'%(head['ion_list'][0],head['ion_sort'][0],head['lorb'],head['ndim']))
f.write('1 5 # number of ireps, dim of irep\n')
norbs = pgroup.hk.shape[2] norbs = pgroup.hk.shape[2]
for isp in xrange(ns_band): for isp in xrange(ns_band):
for ik in xrange(nk): for ik in xrange(nk):

View File

@ -43,7 +43,8 @@ class VaspConverter(ConverterTools):
dft_subgrp = 'dft_input', symmcorr_subgrp = 'dft_symmcorr_input', dft_subgrp = 'dft_input', symmcorr_subgrp = 'dft_symmcorr_input',
parproj_subgrp='dft_parproj_input', symmpar_subgrp='dft_symmpar_input', parproj_subgrp='dft_parproj_input', symmpar_subgrp='dft_symmpar_input',
bands_subgrp = 'dft_bands_input', misc_subgrp = 'dft_misc_input', bands_subgrp = 'dft_bands_input', misc_subgrp = 'dft_misc_input',
transp_subgrp = 'dft_transp_input', repacking = False): transp_subgrp = 'dft_transp_input', repacking = False,
proj_or_hk='proj'):
""" """
Init of the class. Variable filename gives the root of all filenames, e.g. case.ctqmcout, case.h5, and so on. Init of the class. Variable filename gives the root of all filenames, e.g. case.ctqmcout, case.h5, and so on.
""" """
@ -61,7 +62,9 @@ class VaspConverter(ConverterTools):
self.bands_subgrp = bands_subgrp self.bands_subgrp = bands_subgrp
self.misc_subgrp = misc_subgrp self.misc_subgrp = misc_subgrp
self.transp_subgrp = transp_subgrp self.transp_subgrp = transp_subgrp
assert (proj_or_hk == 'proj') or (proj_or_hk == 'hk'), "proj_or_hk has to be 'proj' of 'hk'"
self.proj_or_hk = proj_or_hk
# Checks if h5 file is there and repacks it if wanted: # Checks if h5 file is there and repacks it if wanted:
if (os.path.exists(self.hdf_file) and repacking): if (os.path.exists(self.hdf_file) and repacking):
ConverterTools.repack(self) ConverterTools.repack(self)
@ -147,6 +150,7 @@ class VaspConverter(ConverterTools):
gr_file = self.basename + '.pg%i'%(ig + 1) gr_file = self.basename + '.pg%i'%(ig + 1)
jheader, rf = self.read_header_and_data(gr_file) jheader, rf = self.read_header_and_data(gr_file)
gr_head = json.loads(jheader) gr_head = json.loads(jheader)
e_win = gr_head['ewindow'] e_win = gr_head['ewindow']
nb_max = gr_head['nb_max'] nb_max = gr_head['nb_max']
@ -192,7 +196,7 @@ class VaspConverter(ConverterTools):
n_corr_shells = len(corr_shells) n_corr_shells = len(corr_shells)
#n_shells = n_corr_shells # No non-correlated shells at the moment #n_shells = n_corr_shells # No non-correlated shells at the moment
#shells = corr_shells #shells = corr_shells
n_orbs = sum([sh['dim'] for sh in shells])
# FIXME: atomic sorts in Wien2K are not the same as in VASP. # FIXME: atomic sorts in Wien2K are not the same as in VASP.
# A symmetry analysis from OUTCAR or symmetry file should be used # A symmetry analysis from OUTCAR or symmetry file should be used
# to define equivalence classes of sites. # to define equivalence classes of sites.
@ -231,6 +235,7 @@ class VaspConverter(ConverterTools):
band_window = [numpy.zeros((n_k, 2), dtype=int) for isp in xrange(n_spin_blocs)] band_window = [numpy.zeros((n_k, 2), dtype=int) for isp in xrange(n_spin_blocs)]
n_orbitals = numpy.zeros([n_k, n_spin_blocs], numpy.int) n_orbitals = numpy.zeros([n_k, n_spin_blocs], numpy.int)
for isp in xrange(n_spin_blocs): for isp in xrange(n_spin_blocs):
for ik in xrange(n_k): for ik in xrange(n_k):
ib1, ib2 = int(rf.next()), int(rf.next()) ib1, ib2 = int(rf.next()), int(rf.next())
@ -240,6 +245,27 @@ class VaspConverter(ConverterTools):
for ib in xrange(nb): for ib in xrange(nb):
hopping[ik, isp, ib, ib] = rf.next() hopping[ik, isp, ib, ib] = rf.next()
f_weights[ik, isp, ib] = rf.next() f_weights[ik, isp, ib] = rf.next()
if self.proj_or_hk == 'hk':
hopping = numpy.zeros([n_k, n_spin_blocs, n_orbs, n_orbs], numpy.complex_)
# skip header lines
hk_file = self.basename + '.hk%i'%(ig + 1)
f_hk = open(hk_file, 'rt')
# skip the header (1 line for n_kpoints, n_electrons, n_shells)
# and one line per shell
count = 0
while count < 3 + n_shells:
f_hk.readline()
count += 1
rf_hk = self.read_data(f_hk)
for isp in xrange(n_spin_blocs):
for ik in xrange(n_k):
for ib in xrange(n_orbs):
for jb in xrange(n_orbs):
hopping[ik, isp, ib, jb] = rf_hk.next()
for ib in xrange(n_orbs):
for jb in xrange(n_orbs):
hopping[ik, isp, ib, jb] += 1j*rf_hk.next()
# Projectors # Projectors
# print n_orbitals # print n_orbitals
@ -273,9 +299,9 @@ class VaspConverter(ConverterTools):
# now save only projectors with flag 'corr' to proj_mat # now save only projectors with flag 'corr' to proj_mat
proj_mat = numpy.zeros([n_k, n_spin_blocs, n_corr_shells, sum([crsh['dim'] for crsh in corr_shells]), numpy.max(n_orbitals)], numpy.complex_) proj_mat = numpy.zeros([n_k, n_spin_blocs, n_corr_shells, sum([crsh['dim'] for crsh in corr_shells]), numpy.max(n_orbitals)], numpy.complex_)
addIndex = 0 if self.proj_or_hk == 'proj':
for ish, corr_shell in enumerate(corr_shells): addIndex = 0
if True: for ish, corr_shell in enumerate(corr_shells):
for isp in xrange(n_spin_blocs): for isp in xrange(n_spin_blocs):
for ik in xrange(n_k): for ik in xrange(n_k):
for ion in xrange(len(corr_shell['ion_list'])): for ion in xrange(len(corr_shell['ion_list'])):
@ -285,6 +311,18 @@ class VaspConverter(ConverterTools):
for ib in xrange(n_orbitals[ik, isp]): for ib in xrange(n_orbitals[ik, isp]):
proj_mat[ik,isp,icsh,iclm+addIndex,ib] = proj_mat_all[ik,isp,icsh,ilm,ib] proj_mat[ik,isp,icsh,iclm+addIndex,ib] = proj_mat_all[ik,isp,icsh,ilm,ib]
addIndex += corr_shell['dim'] addIndex += corr_shell['dim']
elif self.proj_or_hk == 'hk':
addIndex = 0
for ish, corr_shell in enumerate(corr_shells):
for isp in xrange(n_spin_blocs):
for ik in xrange(n_k):
for ion in xrange(len(corr_shell['ion_list'])):
icsh = shion_to_shell[ish][ion]
for iclm,ilm in enumerate(xrange(crshorbs_to_globalorbs[ish][0],crshorbs_to_globalorbs[ish][1])):
proj_mat[ik,isp,icsh,iclm+addIndex,ilm] = 1.0
addIndex += corr_shell['dim']
things_to_set = ['n_shells','shells','n_corr_shells','corr_shells','n_spin_blocs','n_orbitals','n_k','SO','SP','energy_unit'] things_to_set = ['n_shells','shells','n_corr_shells','corr_shells','n_spin_blocs','n_orbitals','n_k','SO','SP','energy_unit']
for it in things_to_set: for it in things_to_set:
# print "%s:"%(it), locals()[it] # print "%s:"%(it), locals()[it]
@ -302,7 +340,9 @@ class VaspConverter(ConverterTools):
# The subgroup containing the data. If it does not exist, it is created. If it exists, the data is overwritten! # The subgroup containing the data. If it does not exist, it is created. If it exists, the data is overwritten!
things_to_save = ['energy_unit','n_k','k_dep_projection','SP','SO','charge_below','density_required', things_to_save = ['energy_unit','n_k','k_dep_projection','SP','SO','charge_below','density_required',
'symm_op','n_shells','shells','n_corr_shells','corr_shells','use_rotations','rot_mat', 'symm_op','n_shells','shells','n_corr_shells','corr_shells','use_rotations','rot_mat',
'rot_mat_time_inv','n_reps','dim_reps','T','n_orbitals','proj_mat','proj_mat_all','bz_weights','hopping', 'rot_mat_time_inv','n_reps','dim_reps','T','n_orbitals','proj_mat',
# 'proj_mat_all',
'bz_weights','hopping',
'n_inequiv_shells', 'corr_to_inequiv', 'inequiv_to_corr'] 'n_inequiv_shells', 'corr_to_inequiv', 'inequiv_to_corr']
for it in things_to_save: ar[self.dft_subgrp][it] = locals()[it] for it in things_to_save: ar[self.dft_subgrp][it] = locals()[it]