# This Source Code Form is subject to the terms of the Mozilla Public# License, v. 2.0. If a copy of the MPL was not distributed with this# file, You can obtain one at https://mozilla.org/MPL/2.0/.from__future__importannotationsfromcollectionsimportdefaultdictfromtypingimportOptionalimportnumpyasnpimportscipy.sparseasspsfromscipy.sparseimportlil_matriximportsisl._arrayas_afromsisl._coreimportAtom,AtomicOrbital,Atoms,Geometry,Latticefromsisl.messagesimportwarnfromsisl.physicsimportHamiltonian,Overlapfrom..sileimportadd_sile,sile_fh_openfrom.sileimportSileDFTBdef_lt2full(M):csr=M._csrdiag=csr.diags(csr.diagonal())returnM+M.transpose()-diagclass_realSileDFTB(SileDFTB):def_setup(self,*args,**kwargs):super()._setup(*args,**kwargs)self._comment=["#"]def_r_geometry_info(self):"""Parses the top content of the file for basic info"""# 1. Line will contain the number of atomsline=self.readline()na=int(line)# Create the list of atomic informationnos=[]atom_neighbors=[]foriainrange(na):# 2. contains:line=self.readline()# - the atomic index# - number of neighbors# - number of orbitals on atom_,nneighbor,no=map(int,line.split())nos.append(no)atom_neighbors.append(nneighbor)nos=_a.arrayi(nos)atom_neighbors=_a.arrayi(atom_neighbors)returnna,nos,atom_neighborsdef_r_matrix(self,na,nos,nneighs):# Create some minor variables that aids in the# Hamiltonian construction.no=nos.sum()firsto=np.insert(np.cumsum(nos),0,0)# Create a defaultdict to accommodate arbitrary number of supercells# Typically it doesn't get very big, but this is just to accommodate# different resulting supercells.# The final construction happens in the end.# Remember that DFTB+ only returns the lower triangle of the matrix.Hsc=defaultdict(lambda:lil_matrix((no,no),dtype=np.float64))isc=[0,0,0]foriainrange(na):forinneighinrange(nneighs[ia]):line=self.readline()ia1,ineigh,ja,isc[0],isc[1],isc[2]=map(int,line.split())assertia1-1==iaja-=1# Get the current coupling matrix# The defaultdict will ensure it gets created when neededH=Hsc[tuple(isc)]io=firsto[ia]jo=firsto[ja]joe=jo+nos[ja]foriinrange(nos[ia]):line=self.readline()H[io+i,jo:joe]=list(map(float,line.split()))returnHscdef_get_atoms(self,na,nos):"""Create a ficticious `Atoms` object containing orbitals using the default ordering"""defget_orbital(io):# This is taken directly from the documentation# The documentation lists only these shells to be functionalname=["s","py","pz","px","dxy","dyz","dz2","dxz","dx2-y2","fy(3x2-y2)","fxyz","fz2y","fz3","fz2x","fz(x2-y2)","fx(x2-3y2)",][io]returnAtomicOrbital(name)nos_uniqs=np.unique(nos)defget_atom(ia):nonlocalnos,nos_uniqsno=nos[ia]Z=(nos_uniqs==no).nonzero()[0][0]+1returnAtom(Z,orbitals=[get_orbital(io)forioinrange(no)])returnAtoms([get_atom(ia)foriainrange(na)])@sile_fh_open(True)def_r_file(self,geometry:Optional[Geometry]=None):"""Read content in the current file and return the actual matrices"""na,nos,atom_neighbors=self._r_geometry_info()Msc=self._r_matrix(na,nos,atom_neighbors)# Get supercell sizensc=_a.zerosi(3)foriscinMsc.keys():foriin(0,1,2):nsc[i]=max(abs(isc[i]),nsc[i])# Create the full supercellnsc=nsc*2+1# Create the lattice vectorifgeometryisNone:ifnp.any(nsc>1):warn(f"{self.__class__.__name__}.read_* found a supercell matrix. The returned matrix cannot be used in k-point format since the true atomic coordinates are incorrect, please pass a 'geometry' argument.")lattice=Lattice(na*3+1,nsc=nsc)geometry=Geometry(np.arange(na*3).reshape(na,3),self._get_atoms(na,nos),lattice)else:geometry=geometry.copy()geometry.set_nsc(nsc)# Create the big matrixM=[]foriscingeometry.lattice.sc_off:M.append(Msc[tuple(isc)].tocsr())returngeometry,sps.hstack(M)
[docs]defread_overlap(self,geometry:Optional[Geometry]=None)->Overlap:r"""Parse the output overlap matrix created by DFTB+"""geometry,S=self._r_file(geometry)S.eliminate_zeros()# Convert to classS=Overlap.fromsp(geometry,S)return_lt2full(S)
[docs]defread_overlap(self,geometry:Optional[Geometry]=None)->Overlap:"""Parse the overlap matrix from the ``overreal.dat`` file Parameters ---------- geometry: define the geometry of the Hamiltonian. The data files does *not* contain the geometry information. Hence it can be very useful to retrieve the geometry from somewhere else. """orig_file=self.fileself._file=self.dir_file("overreal.dat")# Read in the overlap matrixgeometry,S=self._r_file(geometry)S.eliminate_zeros()# Return the file-handleself._file=orig_fileS=Overlap.fromsp(geometry,S)return_lt2full(S)
[docs]defread_hamiltonian(self,geometry:Optional[Geometry]=None)->Hamiltonian:r"""Parse the output Hamiltonian created by DFTB+ This will automatically try to discover the ``hamreal[1-4].dat`` and ``overreal.dat`` files in the current directory. As such the single file read is not really done. Parameters ---------- geometry: define the geometry of the Hamiltonian. The data files does *not* contain the geometry information. Hence it can be very useful to retrieve the geometry from somewhere else. """orig_file=self.fileHs=[]foriinrange(1,5):self._file=self.dir_file(f"hamreal{i}.dat")ifnotself._file.exists():continuegeometry,H=self._r_file(geometry)Hs.append(H)# Read in the overlap as wellself._file=self.dir_file("overreal.dat")geometry,S=self._r_file(geometry)# Reset the fileself._file=orig_fileH=Hamiltonian.fromsp(geometry,Hs,S)# Transform back from the charge, x, y, z# to the intrinsic way that sisl handles thingsifH.spin.is_noncolinear:mat=np.empty([4,4])mat[0]=[0.5,0,0,0.5]mat[1]=[0.5,0,0,-0.5]mat[2]=[0,0.5,0,0]mat[3]=[0,0,-0.5,0]H=H.transform(mat)return_lt2full(H)