def orthogonalize(self): print 'NEGF.GF.orthogonalize: Orthogonalizing device region quantities' self.OrthogonalDeviceRegion = True self.HNO = self.H.copy() # nonorthogonal device Hamiltonian (needed) # Device part Usi = MM.mysqrt(self.S) # Folded S Us = LA.inv(Usi) # Store transformation matrices self.Usi, self.Us = Usi, Us # Transform S and H self.S, self.H = MM.mm(Us, self.S, Us), MM.mm(Us, self.H, Us) # Sigmas/Gammas in pyTBT GF can be smaller than device region # First give them the shape of the device region nnL, nnR = len(self.SigL), len(self.SigR) S1, S2 = N.zeros(self.H.shape, N.complex), N.zeros(self.H.shape, N.complex) S1[0:nnL, 0:nnL], S2[-nnR:, -nnR:] = self.SigL, self.SigR # Resetting Sigmas to orthogonalized quantities self.SigL, self.SigR = MM.mm(Us, S1, Us), MM.mm(Us, S2, Us) # ... now the same for the Gammas G1, G2 = N.zeros(self.H.shape, N.complex), N.zeros(self.H.shape, N.complex) G1[0:nnL, 0:nnL], G2[-nnR:, -nnR:] = self.GamL, self.GamR # Resetting Gammas to orthogonalized quantities self.GamL, self.GamR = MM.mm(Us, G1, Us), MM.mm(Us, G2, Us) # Orthogonalize Greens functions self.Gr = MM.mm(Usi, self.Gr, Usi) self.Ga = MM.dagger(self.Gr)
def calcCurrent(options, basis, H, Y): """ Calculate current density in atomic bonds Y : complex scattering state or Y : A_l or A_r! (for total current) """ if isinstance(Y, MM.SpectralMatrix): Y = MM.mm(Y.L, Y.R) NN=len(H) NN2=options.DeviceAtoms[1]-options.DeviceAtoms[0]+1 Curr=N.zeros((NN2, NN2), N.float) if len(Y.shape)==2: for ii in range(NN): a1=basis.ii[ii]-options.DeviceAtoms[0] for jj in range(NN): a2=basis.ii[jj]-options.DeviceAtoms[0] tmp=H[jj, ii]*Y[ii, jj]/2/N.pi # Note that taking the imaginary part is only the valid # expression for Gamma point calculations Curr[a1, a2]=Curr[a1, a2]+4*N.pi*tmp.imag else: for ii in range(NN): a1=basis.ii[ii]-options.DeviceAtoms[0] for jj in range(NN): a2=basis.ii[jj]-options.DeviceAtoms[0] tmp=H[ii, jj]*N.conjugate(Y[ii])*Y[jj] Curr[a1, a2]=Curr[a1, a2]+4*N.pi*tmp.imag return Curr
def Broaden(options, VV, II): """ Broadening corresponding to Lock in measurements for the conductance and IETS spectra. Also resample II, Pow, and nPh to match a common voltage list """ II = II.copy() II = II.real # First derivative dI and bias list dV dI = (II[1:len(II)] - II[:-1]) / (VV[1] - VV[0]) dV = (VV[1:len(VV)] + VV[:-1]) / 2 # Second derivative and bias ddV ddI = (dI[1:len(dI)] - dI[:-1]) / (VV[1] - VV[0]) ddV = (dV[1:len(dV)] + dV[:-1]) / 2 # Modulation amplitude VA = N.sqrt(2.0) * options.Vrms # New bias ranges for broadening tmp = int(N.floor(VA / (dV[1] - dV[0])) + 1) BdV = dV[tmp:-tmp] BddV = ddV[tmp:-tmp] # Initiate derivatives BdI = 0 * BdV BddI = 0 * BddV # Calculate first derivative with Vrms broadening for iV, V in enumerate(BdV): SIO.printDone(iV, len(BdV), 'Inelastica.Broaden: First-derivative Vrms broadening') wt = (N.array(range(200)) / 200.0 - 0.5) * N.pi VL = V + VA * N.sin(wt) dIL = MM.interpolate(VL, dV, dI) BdI[iV] = 2 / N.pi * N.sum(dIL * (N.cos(wt)**2)) * (wt[1] - wt[0]) # Calculate second derivative with Vrms broadening for iV, V in enumerate(BddV): SIO.printDone(iV, len(BddV), 'Inelastica.Broaden: Second-derivative Vrms broadening') wt = (N.array(range(200)) / 200.0 - 0.5) * N.pi VL = V + VA * N.sin(wt) ddIL = MM.interpolate(VL, ddV, ddI) BddI[iV] = 8.0 / 3.0 / N.pi * N.sum(ddIL * (N.cos(wt)**4)) * (wt[1] - wt[0]) # Reduce to one voltage grid NN = options.biasPoints V = N.linspace(options.minBias, options.maxBias, NN) NI = MM.interpolate(V, VV, II) NdI = MM.interpolate(V, dV, dI) NddI = MM.interpolate(V, ddV, ddI) NBdI = MM.interpolate(V, BdV, BdI) NBddI = MM.interpolate(V, BddV, BddI) return V, NI, NdI, NddI, NBdI, NBddI
def __calcEigChan(self, A1, G2, Left, channels=10): # Calculate Eigenchannels using recipe from PRB # For right eigenchannels, A1=A2, G2=G1 !!! if isinstance(A1, MM.SpectralMatrix): ev, U = LA.eigh(MM.mm(A1.L, A1.R)) else: ev, U = LA.eigh(A1) # This small trick will remove all zero contribution vectors # and will diagonalize the tt matrix in the subspace where there # are values. idx = (ev > 0).nonzero()[0] ev = N.sqrt(ev[idx] / (2 * N.pi)) ev.shape = (1, -1) Utilde = ev * U[:, idx] nuo, nuoL, nuoR = self.nuo, self.nuoL, self.nuoR if Left: tt = MM.mm(MM.dagger(Utilde[nuo - nuoR:nuo, :]), 2 * N.pi * G2, Utilde[nuo - nuoR:nuo, :]) else: tt = MM.mm(MM.dagger(Utilde[:nuoL, :]), 2 * N.pi * G2, Utilde[:nuoL, :]) # Diagonalize (note that this is on a reduced tt matrix (no 0 contributing columns) evF, UF = LA.eigh(tt) EC = MM.mm(Utilde, UF[:, -channels:]).T return EC[::-1, :], evF[::-1] # reverse eigenvalues
def calcTEIG(self, channels=10): # Transmission matrix (complex array) TT = self.TT Trans = N.trace(TT) VC.Check("trans-imaginary-part", Trans.imag, "Transmission has large imaginary part") # Calculate eigenchannel transmissions too tval, tvec = LA.eig(TT) idx = (tval.real).argsort()[::-1] # sort from largest to smallest tval = tval[idx] tvec = tvec[:, idx] # Compute shot noise Smat = MM.mm(TT, N.identity(len(TT)) - TT) sval = N.diag(MM.mm(MM.dagger(tvec), Smat, tvec)) # set up arrays T = N.zeros(channels + 1) SN = N.zeros(channels + 1) T[0] = Trans.real SN[0] = N.trace(Smat).real for i in range(min(channels, len(TT))): T[i + 1] = tval[i].real SN[i + 1] = sval[i].real return T, SN
def genkmesh(self): "Generate mesh for a given sampling of the axes" self.k = [] self.w = [] errorw = [] for i in range(3): # loop over the three k-components if self.type[i].upper() == 'GK' or self.type[i].upper( ) == 'GAUSSKRONROD': self.type[i] = 'GK' if self.Nk[i] > 1: # GK-method fails with fewer points kpts, wgts, ew = MM.GaussKronrod(self.Nk[i]) self.k.append(kpts) self.w.append(wgts) errorw.append(ew) else: print 'Kmesh.py: GK method requires Nk=%i>1' % (self.Nk[i]) kuk elif self.type[i].upper() == 'LIN' or self.type[i].upper( ) == 'LINEAR': self.type[i] = 'LIN' kpts, wgts = generatelinmesh(self.Nk[i]) self.k.append(kpts) self.w.append(wgts) errorw.append(wgts) else: print 'Kmesh.py: Unknown meshtype:', self.type[i].upper() self.Nk[i] = len(self.k[i]) self.NNk = N.prod(self.Nk) print 'Kmesh.py: Generating mesh:' print ' ... type = ', self.type print ' ... Nk = ', self.Nk # repete out in 3D kpts = N.zeros((self.NNk, 3)) # Array of k-points wgts = N.ones((4, self.NNk)) # (wgts, errorw1, errorw2, errorw3) nn = 0 for i in range(self.Nk[0]): for j in range(self.Nk[1]): for k in range(self.Nk[2]): kpts[nn, :] = [self.k[0][i], self.k[1][j], self.k[2][k]] wgts[0, nn] = self.w[0][i] * self.w[1][j] * self.w[2][k] wgts[1, nn] = errorw[0][i] * self.w[1][j] * self.w[2][k] wgts[2, nn] = self.w[0][i] * errorw[1][j] * self.w[2][k] wgts[3, nn] = self.w[0][i] * self.w[1][j] * errorw[2][k] nn += 1 print ' ... NNk = %i, sum(wgts) = %.8f' % (self.NNk, N.sum(wgts[0])) print ' ... sum(errorw) = (%.8f,%.8f,%.8f)' % tuple( N.sum(wgts[i + 1]) for i in range(3)) self.k = kpts self.w = wgts
def ComputeElectronStates(self, kpoint, verbose=True, TSrun=False): if TSrun: # kpoint has unit of '2*pi/a' kpt = kpoint / (2 * N.pi) kpt2 = MM.mm(N.array([self.Sym.a1, self.Sym.a2, self.Sym.a3]), kpt) self.TSHS0.setkpoint(kpt2, atype=N.complex, verbose=verbose) self.h0_k = self.TSHS0.H[:, :, :] self.s0_k = self.TSHS0.S[:, :] else: if verbose: print 'SupercellPhonons.ComputeElectronStates: k = ', kpoint, '(1/Ang)' # Fold onto primitive cell self.h0_k = self.Fold2PrimitiveCell(self.h0, kpoint) self.s0_k = self.Fold2PrimitiveCell(self.s0, kpoint) ev = N.empty((self.nspin, self.rednao), N.float) evec = N.empty((self.nspin, self.rednao, self.rednao), N.complex) for ispin in range(self.nspin): ev[ispin], evec[ispin] = SLA.eigh(self.h0_k[ispin], self.s0_k) return ev, evec
def calcWF2(options, geom, DeviceAtoms, basis, Y, NN, Fold=True, k=[0, 0, 0], a=None, fn=''): """ Calculate wavefunction on real space mesh with optional folding of the periodic boundary conditions. If Fold==True the vectors 'a' will be choosen to be the periodic boundary condition vectors and the real space wavefunction folded into the cell using the k-vector. If Folded==False, you can choose any 'a' but folding by the PBC will not be done. Note: Folding assumes coupling only between N.N. cells INPUT: geom : MakeGeom structure for full geometry DeviceAtoms : [first, last] numbering from 1 to N basis : Basis struct for ONLY the device region Y : Wavefunction for the device region NN : [N1,N2,N3,minN3,maxN3] number of points along [a1, a2, a3] only calculate from minN3 to maxN3 Fold : fold periodic boundary conditions k : k-vector to use for folding [-0.5,0.5] a : if Fold==True, uses geom.pbc, if false: [a1, a2, a3] along these directions, i.e., da1=a1/N1 ... RETURNS: YY : complex wavefunction as N1, N2, N3 matrix """ xyz=N.array(geom.xyz[DeviceAtoms[0]-1:DeviceAtoms[1]]) N1, N2, N3, minN3, maxN3 = NN[0], NN[1], NN[2], NN[3], NN[4] NN3 = maxN3-minN3+1 if Fold: da1, da2, da3 = geom.pbc[0]/N1, geom.pbc[1]/N2, geom.pbc[2]/N3 else: da1, da2, da3 = a[0]/N1, a[1]/N2, a[2]/N3 try: NNY=Y.shape[1] except: NNY=1 # Def cube YY=[N.zeros((N1, N2, NN3), N.complex) for ii in range(NNY)] # Help indices i1 = MM.outerAdd(N.arange(N1), N.zeros(N2), N.zeros(NN3)) i2 = MM.outerAdd(N.zeros(N1), N.arange(N2), N.zeros(NN3)) i3 = MM.outerAdd(N.zeros(N1), N.zeros(N2), N.arange(NN3)+minN3) # Positions x,y,z for points in matrix form rx = i1*da1[0]+i2*da2[0]+i3*da3[0] ry = i1*da1[1]+i2*da2[1]+i3*da3[1] rz = i1*da1[2]+i2*da2[2]+i3*da3[2] if Fold: pbc = N.array([da1*N1, da2*N2, da3*N3]) orig = da3*minN3 b = N.transpose(LA.inv(pbc)) olist = [orig, orig+pbc[0]+pbc[1]+pbc[2]] pairs = N.array([[b[1], b[2]], [b[0], b[2]], [b[0], b[1]]]) pbcpairs = N.array([[pbc[1], pbc[2]], [pbc[0], pbc[2]], [pbc[0], pbc[1]]]) for iiatom in range(DeviceAtoms[1]-DeviceAtoms[0]+1): if iiatom>0: #Wavefunction percent done... SIO.printDone(iiatom, DeviceAtoms[1]-DeviceAtoms[0]+1, 'Wave function') if not Fold: shiftvec = [0, 0, 0] else: # include shift vectors if sphere cuts planes of PBC basisindx = N.where(basis.ii==iiatom+DeviceAtoms[0])[0] basisradius = N.max(basis.coff[basisindx]) minshift, maxshift = [0, 0, 0], [0, 0, 0] for iithiso, thiso in enumerate(olist): c=xyz[iiatom, :]-thiso for ip in range(3): n = N.cross(pbcpairs[ip, 0], pbcpairs[ip, 1]) n = n / LA.norm(n) dist = N.abs(N.dot(n, N.dot(N.dot(pairs[ip], c), pbcpairs[ip])-c)) if dist < basisradius: if iithiso==0: minshift[ip]=-1 else: maxshift[ip]=1 shiftvec = [] for ix in range(minshift[0], maxshift[0]+1): for iy in range(minshift[1], maxshift[1]+1): for iz in range(minshift[2], maxshift[2]+1): shiftvec+=[[ix, iy, iz]] #print(shiftvec) for shift in shiftvec: # Atom position shifted vec=-N.dot(N.array(shift), geom.pbc) phase = N.exp(-1.0j*(2*N.pi*N.dot(N.array(k), N.array(shift)))) # Difference and polar, cylinder distances dx, dy, dz=rx-xyz[iiatom, 0]-vec[0], ry-xyz[iiatom, 1]-vec[1], rz-xyz[iiatom, 2]-vec[2] dr2=dx*dx+dy*dy+dz*dz drho2=dx*dx+dy*dy basisindx = N.where(basis.ii==iiatom+DeviceAtoms[0])[0] old_coff, old_delta = 0.0, 0.0 for basisorb in basisindx: if not (basis.coff[basisorb]==old_coff and basis.delta[basisorb]==old_delta): old_coff, old_delta = basis.coff[basisorb], basis.delta[basisorb] indx = N.where(dr2<basis.coff[basisorb]**2) # Find points close to atom idr, idrho=N.sqrt(dr2[indx]), N.sqrt(drho2[indx]) iri=(idr/basis.delta[basisorb]).astype(N.int) idx, idy, idz = dx[indx], dy[indx], dz[indx] costh = idz/idr cosfi, sinfi = idx/idrho, idy/idrho # Fix divide by zeros indxRzero, indxRhozero = N.where(idr==0.0), N.where(idrho==0.0) costh[indxRzero] = 1.0 cosfi[indxRhozero], sinfi[indxRhozero] = 1.0, 0.0 # Numpy has changed the choose function to crap! RR=N.take(basis.orb[basisorb], iri) # Calculate spherical harmonics if len(idr)>0: l = basis.L[basisorb] m = basis.M[basisorb] if l==3: print('f-shell : l=%i, m=%i (NOT TESTED!!)'%(l, m)) thisSphHar = MM.sphericalHarmonics(l, m, costh, sinfi, cosfi) for iy in range(NNY): YY[iy][indx]=YY[iy][indx]+RR*thisSphHar*Y[basisorb, iy]*phase N1, N2, N3, minN3, maxN3 = NN[0], NN[1], NN[2], NN[3], NN[4] for ii in range(NNY): writenetcdf2(geom, fn+'%i.nc'%ii, YY[ii], N1, N2, N3, minN3, maxN3, geom.pbc, options.DeviceAtoms) return YY
def writeFGRrates(options, GF, hw, NCfile): print 'Inelastica.writeFGRrates: Computing FGR rates' # Eigenchannels GF.calcEigChan(channels=options.numchan) NCfile = NC4.Dataset(options.PhononNetCDF, 'r') print 'Reading ', options.PhononNetCDF outFile = file('%s/%s.IN.FGR' % (options.DestDir, options.systemlabel), 'w') outFile.write('Total transmission [in units of (1/s/eV)] : %e\n' % (PC.unitConv * GF.TeF, )) for ihw in range(len(hw)): SIO.printDone(ihw, len(hw), 'Golden Rate') M = N.array(NCfile.variables['He_ph'][ihw, options.iSpin, :, :], N.complex) try: M += 1.j * N.array( NCfile.variables['ImHe_ph'][ihw, options.iSpin, :, :], N.complex) except: print 'Warning: Variable ImHe_ph not found' rate = N.zeros((len(GF.ECleft), len(GF.ECright)), N.float) totrate = 0.0 inter, intra = 0.0, 0.0 # splitting total rate in two for iL in range(len(GF.ECleft)): for iR in range(len(GF.ECright)): tmp = N.dot(N.conjugate(GF.ECleft[iL]), MM.mm(M, GF.ECright[iR])) rate[iL, iR] = (2 * N.pi)**2 * abs(tmp)**2 totrate += rate[iL, iR] if iL == iR: intra += rate[iL, iR] else: inter += rate[iL, iR] outFile.write( '\nPhonon mode %i : %f eV [Rates in units of (1/s/eV)]\n' % (ihw, hw[ihw])) outFile.write( 'eh-damp : %e (1/s) , heating %e (1/(sV)))\n' % (GF.P1T[ihw] * PC.unitConv * hw[ihw], GF.P2T[ihw] * PC.unitConv)) outFile.write( 'eh-damp 1, 2 (MALMAL, MARMAR): %e (1/s) , %e (1/(s)))\n' % (GF.ehDampL[ihw] * PC.unitConv * hw[ihw], GF.ehDampR[ihw] * PC.unitConv * hw[ihw])) outFile.write('SymI : %e (1/(sV)) , AsymI %e (?))\n' % (GF.nHT[ihw] * PC.unitConv, GF.HT[ihw] * PC.unitConv)) #outFile.write('Elast : %e (1/(sV)) , Inelast %e (1/(sV)))\n' % (GF.nHTel[ihw]*PC.unitConv,GF.nHTin[ihw]*PC.unitConv)) outFile.write('down=left EC, right=right EC\n') if GF.P2T[ihw] > 0.0: if abs(totrate / (GF.P2T[ihw]) - 1) < 0.05: outFile.write('Sum/Tr[MALMAR] , Tr: %1.3f %e\n' % (totrate / (GF.P2T[ihw]), PC.unitConv * GF.P2T[ihw])) else: outFile.write( 'WARNING: !!!! Sum/Tr[MALMAR] , Tr: %2.2e %e\n' % (totrate / (GF.P2T[ihw]), PC.unitConv * GF.P2T[ihw])) else: outFile.write(' Tr: %e\n' % (PC.unitConv * GF.P2T[ihw])) inter = inter / GF.P2T[ihw] intra = intra / GF.P2T[ihw] outFile.write( 'Interchannel ratio: Sum(inter)/Tr[MALMAR] = %.4f \n' % inter) outFile.write( 'Intrachannel ratio: Sum(intra)/Tr[MALMAR] = %.4f \n' % intra) outFile.write( 'Inter+intra ratio: Sum(inter+intra)/Tr[MALMAR] = %.4f \n' % (inter + intra)) for iL in range(len(GF.ECleft)): for iR in range(len(GF.ECright)): outFile.write('%e ' % (PC.unitConv * rate[iL, iR], )) outFile.write('\n') outFile.close()
def calcIETS(options, GFp, GFm, basis, hw): # Calculate product of electronic traces and voltage functions print 'Inelastica.calcIETS: nHTp =', GFp.nHT * PC.unitConv # OK print 'Inelastica.calcIETS: nHTm =', GFm.nHT * PC.unitConv # OK print 'Inelastica.calcIETS: HTp =', GFp.HT # OK print 'Inelastica.calcIETS: HTm =', GFm.HT # OK # Set up grid and Hilbert term kT = options.Temp / 11604.0 # (eV) # Generate grid for numerical integration of Hilbert term max_hw = max(hw) max_win = max(-options.minBias, max_hw) + 20 * kT + 4 * options.Vrms min_win = min(-options.maxBias, -max_hw) - 20 * kT - 4 * options.Vrms pts = int(N.floor((max_win - min_win) / kT * 3)) Egrid = N.linspace(min_win, max_win, pts) print "Inelastica.calcIETS: Hilbert integration grid : %i pts [%f,%f]" % ( pts, min(Egrid), max(Egrid)) NN = options.biasPoints print 'Inelastica.calcIETS: Biaspoints =', NN # Add some points for the Lock in broadening approxdV = (options.maxBias - options.minBias) / (NN - 1) NN += int(((8 * options.Vrms) / approxdV) + .5) Vl = N.linspace(options.minBias - 4 * options.Vrms, options.maxBias + 4 * options.Vrms, NN) # Vector implementation on Vgrid: wp = (1 + N.sign(Vl)) / 2. # weights for positive V wm = (1 - N.sign(Vl)) / 2. # weights for negative V # Mode occupation and power dissipation Pow = N.zeros((len(hw), NN), N.float) # (modes,Vgrid) nPh = N.zeros((len(hw), NN), N.float) t0 = N.clip(Vl / kT, -700, 700) cosh0 = N.cosh(t0) # Vgrid sinh0 = N.sinh(t0) for i in (hw > options.modeCutoff).nonzero()[0]: P1T = wm * GFm.P1T[i] + wp * GFp.P1T[i] P2T = wm * GFm.P2T[i] + wp * GFp.P2T[i] # Bose distribution nB = 1 / (N.exp(N.clip(hw[i] / kT, -300, 300)) - 1) # number t1 = N.clip(hw[i] / (2 * kT), -700, 700) # number coth1 = N.cosh(t1) / N.sinh(t1) # Emission rate and e-h damping damp = P1T * hw[i] / N.pi # Vgrid emis = P2T * (hw[i] * (cosh0 - 1) * coth1 - Vl * sinh0) / (N.cosh(2 * t1) - cosh0) / N.pi # Determine mode occupation if options.PhHeating: nPh[i, :] = emis / (hw[i] * P1T / N.pi + options.PhExtDamp) + nB else: nPh[i, :] = nB # Mode-resolved power dissipation Pow[i, :] = hw[i] * ((nB - nPh[i]) * damp + emis) # Current: non-Hilbert part (InH) InH = N.zeros((NN, ), N.float) # Vgrid IsymF = N.zeros((NN, ), N.float) for i in (hw > options.modeCutoff).nonzero()[0]: nHT = wm * GFm.nHT[i] + wp * GFp.nHT[i] # Vgrid t1 = hw[i] / (2 * kT) # number t1 = N.clip(t1, -700, 700) coth1 = N.cosh(t1) / N.sinh(t1) t2 = (hw[i] + Vl) / (2 * kT) # Vgrid t2 = N.clip(t2, -700, 700) coth2 = N.cosh(t2) / N.sinh(t2) t3 = (hw[i] - Vl) / (2 * kT) # Vgrid t3 = N.clip(t3, -700, 700) coth3 = N.cosh(t3) / N.sinh(t3) # Vgrid # Isym function Isym = 0.5 * (hw[i] + Vl) * (coth1 - coth2) # Vgrid Isym -= 0.5 * (hw[i] - Vl) * (coth1 - coth3) # non-Hilbert part InH += (Isym + 2 * Vl * nPh[i]) * nHT # Vgrid IsymF += Isym # Current: Add Landauer part, GFm.TeF = GFp.TeF InH += GFp.TeF * Vl # Vgrid # Current: Asymmetric/Hilbert part (IH) try: import scipy.special as SS print "Inelastica: Computing asymmetric term using digamma function," print "... see G. Bevilacqua et al., Eur. Phys. J. B (2016) 89: 3" IH = N.zeros((NN, ), N.float) IasymF = N.zeros((NN, ), N.float) for i in (hw > options.modeCutoff).nonzero()[0]: v0 = hw[i] / (2 * N.pi * kT) vp = (hw[i] + Vl) / (2 * N.pi * kT) vm = (hw[i] - Vl) / (2 * N.pi * kT) Iasym = kT * (2 * v0 * SS.psi(1.j * v0) - vp * SS.psi(1.j * vp) - vm * SS.psi(1.j * vm)).real IasymF += Iasym IH += GFp.HT[i] * N.array(Vl > 0.0, dtype=int) * Iasym IH += GFm.HT[i] * N.array(Vl < 0.0, dtype=int) * Iasym except: print "Computing using explit Hilbert transformation" IH = N.zeros((NN, ), N.float) IasymF = N.zeros((NN, ), N.float) # Prepare box/window function on array Vl2 = N.outer(Vl, N.ones(Egrid.shape)) Egrid2 = N.outer(N.ones(Vl.shape), Egrid) # Box/window function nF(E-Vl2)-nF(E-0): kasse = MM.box(0, -Vl2, Egrid2, kT) # (Vgrid,Egrid) ker = None for i in (hw > options.modeCutoff).nonzero()[0]: # Box/window function nF(E-hw)-nF(E+hw) tmp = MM.box(-hw[i], hw[i], Egrid, kT) hilb, ker = MM.Hilbert(tmp, ker) # Egrid # Calculate Iasym for each bias point for j in range(len(Vl)): Iasym = MM.trapez(Egrid, kasse[j] * hilb, equidistant=True).real / 2 IasymF[j] += Iasym if Vl[j] > 0: IH[j] += GFp.HT[i] * Iasym else: IH[j] += GFm.HT[i] * Iasym # Compute inelastic shot noise terms here: absVl = N.absolute(Vl) Inew = N.zeros(len(Vl), N.float) Snew = N.zeros(len(Vl), N.float) print 'Noise factors:' print GFp.dIel print GFp.dIinel print GFp.dSel print GFp.dSinel for i in (hw > options.modeCutoff).nonzero()[0]: # Elastic part Inew += GFp.dIel[i] * Vl Snew += GFp.dSel[i] * absVl # Inelastic part indx = (absVl - hw[i] < 0).nonzero()[0] fct = absVl - hw[i] fct[indx] = 0.0 # set elements to zero Inew += GFp.dIinel[i] * fct * N.sign(Vl) Snew += GFp.dSinel[i] * fct # Get the right units for gamma_eh, gamma_heat gamma_eh_p = N.zeros((len(hw), ), N.float) gamma_eh_m = N.zeros((len(hw), ), N.float) gamma_heat_p = N.zeros((len(hw), ), N.float) gamma_heat_m = N.zeros((len(hw), ), N.float) for i in (hw > options.modeCutoff).nonzero()[0]: # Units [Phonons per Second per dN where dN is number extra phonons] gamma_eh_p[i] = GFp.P1T[i] * hw[i] * PC.unitConv gamma_eh_m[i] = GFm.P1T[i] * hw[i] * PC.unitConv # Units [Phonons per second per eV [eV-ihw] gamma_heat_p[i] = GFp.P2T[i] * PC.unitConv gamma_heat_m[i] = GFm.P2T[i] * PC.unitConv print 'Inelastica.calcIETS: gamma_eh_p =', gamma_eh_p # OK print 'Inelastica.calcIETS: gamma_eh_m =', gamma_eh_m # OK print 'Inelastica.calcIETS: gamma_heat_p =', gamma_heat_p # OK print 'Inelastica.calcIETS: gamma_heat_m =', gamma_heat_m # OK V, I, dI, ddI, BdI, BddI = Broaden(options, Vl, InH + IH) V, Is, dIs, ddIs, BdIs, BddIs = Broaden(options, Vl, IsymF) V, Ia, dIa, ddIa, BdIa, BddIa = Broaden(options, Vl, IasymF) # Interpolate quantities to new V-grid NPow = N.zeros((len(Pow), len(V)), N.float) NnPh = N.zeros((len(Pow), len(V)), N.float) for ii in range(len(Pow)): NPow[ii] = MM.interpolate(V, Vl, Pow[ii]) NnPh[ii] = MM.interpolate(V, Vl, nPh[ii]) # Interpolate inelastic noise NV, NI, NdI, NddI, NBdI, NBddI = Broaden(options, Vl, GFp.TeF * Vl + Inew) NV, NS, NdS, NddS, NBdS, NBddS = Broaden(options, Vl, Snew) print 'Inelastica.calcIETS: V[:5] =', V[:5] # OK print 'Inelastica.calcIETS: V[-5:][::-1] =', V[-5:][::-1] # OK print 'Inelastica.calcIETS: I[:5] =', I[:5] # OK print 'Inelastica.calcIETS: I[-5:][::-1] =', I[-5:][::-1] # OK print 'Inelastica.calcIETS: BdI[:5] =', BdI[:5] # OK print 'Inelastica.calcIETS: BdI[-5:][::-1] =', BdI[-5:][::-1] # OK print 'Inelastica.calcIETS: BddI[:5] =', BddI[:5] # OK print 'Inelastica.calcIETS: BddI[-5:][::-1] =', BddI[-5:][::-1] # OK datafile = '%s/%s.IN' % (options.DestDir, options.systemlabel) # ascii format writeLOEData2Datafile(datafile + 'p', hw, GFp.TeF, GFp.nHT, GFp.HT) writeLOEData2Datafile(datafile + 'm', hw, GFm.TeF, GFm.nHT, GFm.HT) # netcdf format outNC = initNCfile(datafile, hw, V) write2NCfile(outNC, BddI / BdI, 'IETS', 'Broadened BddI/BdI [1/V]') write2NCfile(outNC, ddI / dI, 'IETS_0', 'Intrinsic ddI/dI [1/V]') write2NCfile(outNC, BdI, 'BdI', 'Broadened BdI, G0') write2NCfile(outNC, BddI, 'BddI', 'Broadened BddI, G0') write2NCfile(outNC, I, 'I', 'Intrinsic I, G0 V') write2NCfile(outNC, dI, 'dI', 'Intrinsic dI, G0') write2NCfile(outNC, ddI, 'ddI', 'Intrinsic ddI, G0/V') if options.LOEscale == 0.0: write2NCfile( outNC, GFp.nHT, 'ISymTr', 'Trace giving Symmetric current contribution (prefactor to universal function)' ) write2NCfile( outNC, GFp.HT, 'IAsymTr', 'Trace giving Asymmetric current contribution (prefactor to universal function)' ) write2NCfile(outNC, gamma_eh_p, 'gamma_eh', 'e-h damping [*deltaN=1/Second]') write2NCfile(outNC, gamma_heat_p, 'gamma_heat', 'Phonon heating [*(bias-hw) (eV) = 1/Second]') # New stuff related to the noise implementation write2NCfile( outNC, NI, 'Inew', 'Intrinsic Inew (new implementation incl. elastic renormalization, T=0)' ) write2NCfile( outNC, NdI, 'dInew', 'Intrinsic dInew (new implementation incl. elastic renormalization, T=0)' ) write2NCfile( outNC, NddI, 'ddInew', 'Intrinsic ddInew (new implementation incl. elastic renormalization, T=0)' ) write2NCfile( outNC, NddI / NdI, 'IETSnew_0', 'Intrinsic ddInew/dInew (new implementation incl. elastic renormalization, T=0) [1/V]' ) write2NCfile( outNC, NBdI, 'BdInew', 'Broadened BdInew (new implementation incl. elastic renormalization, T=0)' ) write2NCfile( outNC, NBddI, 'BddInew', 'Broadened BddInew (new implementation incl. elastic renormalization, T=0)' ) write2NCfile( outNC, NBddI / NBdI, 'IETSnew', 'Broadened BddInew/BdInew (new implementation incl. elastic renormalization, T=0) [1/V]' ) write2NCfile( outNC, NdS, 'dSnew', 'Inelastic first-derivative of the shot noise dSnew (T=0)') else: write2NCfile( outNC, GFp.nHT, 'ISymTr_p', 'Trace giving Symmetric current contribution (prefactor to universal function)' ) write2NCfile( outNC, GFp.HT, 'IAsymTr_p', 'Trace giving Asymmetric current contribution (prefactor to universal function)' ) write2NCfile( outNC, GFm.nHT, 'ISymTr_m', 'Trace giving Symmetric current contribution (prefactor to universal function)' ) write2NCfile( outNC, GFm.HT, 'IAsymTr_m', 'Trace giving Asymmetric current contribution (prefactor to universal function)' ) write2NCfile(outNC, gamma_eh_p, 'gamma_eh_p', 'e-h damping [*deltaN=1/Second]') write2NCfile(outNC, gamma_heat_p, 'gamma_heat_p', 'Phonon heating [*(bias-hw) (eV) = 1/Second]') write2NCfile(outNC, gamma_eh_m, 'gamma_eh_m', 'e-h damping [*deltaN=1/Second]') write2NCfile(outNC, gamma_heat_m, 'gamma_heat_m', 'Phonon heating [*(bias-hw) (eV) = 1/Second]') # Phonon occupations and power balance write2NCfile(outNC, NnPh, 'nPh', 'Number of phonons') write2NCfile(outNC, N.sum(NnPh, axis=0), 'nPh_tot', 'Total number of phonons') write2NCfile(outNC, NPow, 'Pow', 'Mode-resolved power balance') write2NCfile(outNC, N.sum(NPow, axis=0), 'Pow_tot', 'Total power balance') # Write "universal functions" write2NCfile(outNC, dIs, 'dIs', 'dIsym function') write2NCfile(outNC, dIa, 'dIa', 'dIasym function') write2NCfile(outNC, ddIs, 'ddIs', 'ddIasym function') write2NCfile(outNC, ddIa, 'ddIa', 'ddIasym function') # Write energy reference where Greens functions are evaluated outNC.createDimension('number', 1) tmp = outNC.createVariable('EnergyRef', 'd', ('number', )) tmp[:] = N.array(options.energy) # Write LOEscale tmp = outNC.createVariable('LOEscale', 'd', ('number', )) tmp[:] = N.array(options.LOEscale) # Write k-point outNC.createDimension('vector', 3) tmp = outNC.createVariable('kpoint', 'd', ('vector', )) tmp[:] = N.array(options.kpoint) outNC.close() return V, I, dI, ddI, BdI, BddI
def calcTraces(options, GF1, GF2, basis, NCfile, ihw): # Calculate various traces over the electronic structure # Electron-phonon couplings ihw = int(ihw) M = N.array(NCfile.variables['He_ph'][ihw, options.iSpin, :, :], N.complex) try: M += 1.j * N.array( NCfile.variables['ImHe_ph'][ihw, options.iSpin, :, :], N.complex) except: print 'Warning: Variable ImHe_ph not found' # Calculation of intermediate quantity MARGLGM = MM.mm(M, GF1.ARGLG, M) MARGLGM2 = MM.mm(M, GF2.ARGLG, M) # LOE expressions in compact form t1 = MM.mm(MARGLGM, GF2.AR) t2 = MM.mm(MARGLGM2, GF1.AL) # Note that compared with Eq. (10) of PRB89, 081405 (2014) we here use # the definition B_lambda = MM.trace(t1-dagger(t2)), which in turn gives # ReB = MM.trace(t1).real-MM.trace(t2).real # ImB = MM.trace(t1).imag+MM.trace(t2).imag K23 = MM.trace(t1).imag + MM.trace(t2).imag K4 = MM.trace(MM.mm(M, GF1.ALT, M, GF2.AR)) aK23 = 2 * (MM.trace(t1).real - MM.trace(t2).real) # asymmetric part # Non-Hilbert term defined here with a minus sign GF1.nHT[ihw] = NEGF.AssertReal(K23 + K4, 'nHT[%i]' % ihw) GF1.HT[ihw] = NEGF.AssertReal(aK23, 'HT[%i]' % ihw) # Power, damping and current rates GF1.P1T[ihw] = NEGF.AssertReal(MM.trace(MM.mm(M, GF1.A, M, GF2.A)), 'P1T[%i]' % ihw) GF1.P2T[ihw] = NEGF.AssertReal(MM.trace(MM.mm(M, GF1.AL, M, GF2.AR)), 'P2T[%i]' % ihw) GF1.ehDampL[ihw] = NEGF.AssertReal(MM.trace(MM.mm(M, GF1.AL, M, GF2.AL)), 'ehDampL[%i]' % ihw) GF1.ehDampR[ihw] = NEGF.AssertReal(MM.trace(MM.mm(M, GF1.AR, M, GF2.AR)), 'ehDampR[%i]' % ihw) # Remains from older version (see before rev. 219): #GF.dGnout.append(EC.calcCurrent(options,basis,GF.HNO,mm(Us,-0.5j*(tmp1-dagger(tmp1)),Us))) #GF.dGnin.append(EC.calcCurrent(options,basis,GF.HNO,mm(Us,mm(G,MA1M,Gd)-0.5j*(tmp2-dagger(tmp2)),Us))) # NB: TF Should one use GF.HNO (nonorthogonal) or GF.H (orthogonalized) above? if options.LOEscale == 0.0: # Check against original LOE-WBA formulation isym1 = MM.mm(GF1.ALT, M, GF2.AR, M) isym2 = MM.mm(MM.dagger(GF1.ARGLG), M, GF2.A, M) isym3 = MM.mm(GF1.ARGLG, M, GF2.A, M) isym = MM.trace(isym1) + 1j / 2. * (MM.trace(isym2) - MM.trace(isym3)) print 'LOE-WBA check: Isym diff', K23 + K4 - isym iasym1 = MM.mm(MM.dagger(GF1.ARGLG), M, GF2.AR - GF2.AL, M) iasym2 = MM.mm(GF1.ARGLG, M, GF2.AR - GF2.AL, M) iasym = MM.trace(iasym1) + MM.trace(iasym2) print 'LOE-WBA check: Iasym diff', aK23 - iasym # Compute inelastic shot noise terms according to the papers # Haupt, Novotny & Belzig, PRB 82, 165441 (2010) and # Avriller & Frederiksen, PRB 86, 155411 (2012) # Zero-temperature limit TT = MM.mm(GF1.GammaL, GF1.AR) # this matrix has the correct shape for MM ReGr = (GF1.Gr + GF1.Ga) / 2. tmp = MM.mm(GF1.Gr, M, ReGr, M, GF1.AR) tmp = tmp + MM.dagger(tmp) Tlambda0 = MM.mm(GF1.GammaL, tmp) tmp1 = MM.mm(M, GF1.AR, M) tmp2 = MM.mm(M, GF1.A, M, GF1.Gr, GF1.GammaR) tmp = tmp1 + 1j / 2. * (MM.dagger(tmp2) - tmp2) Tlambda1 = MM.mm(GF1.GammaL, GF1.Gr, tmp, GF1.Ga) MARGL = MM.mm(M, GF1.AR, GF1.GammaL) tmp1 = MM.mm(MARGL, GF1.AR, M) tmp2 = MM.mm(MARGL, GF1.Gr, M, GF1.Gr, GF1.GammaR) tmp = tmp1 + tmp2 tmp = tmp + MM.dagger(tmp) Qlambda = MM.mm(-GF1.Ga, GF1.GammaL, GF1.Gr, tmp) tmp = -2 * TT OneMinusTwoT = tmp + N.identity(len(GF1.GammaL)) # Store relevant traces GF1.dIel[ihw] = NEGF.AssertReal(MM.trace(Tlambda0), 'dIel[%i]' % ihw) GF1.dIinel[ihw] = NEGF.AssertReal(MM.trace(Tlambda1), 'dIinel[%i]' % ihw) GF1.dSel[ihw] = NEGF.AssertReal( MM.trace(MM.mm(OneMinusTwoT, Tlambda0)), 'dSel[%i]' % ihw) GF1.dSinel[ihw] = NEGF.AssertReal( MM.trace(Qlambda + MM.mm(OneMinusTwoT, Tlambda1)), 'dSinel[%i]' % ihw)
def main(options): """ Main routine to compute inelastic transport characteristics (dI/dV, d2I/dV2, IETS etc) Parameters ---------- options : an ``options`` instance """ CF.CreatePipeOutput(options.DestDir + '/' + options.Logfile) VC.OptionsCheck(options, 'Inelastica') CF.PrintMainHeader('Inelastica', options) options.XV = '%s/%s.XV' % (options.head, options.systemlabel) options.geom = MG.Geom(options.XV, BufferAtoms=options.buffer) # Voltage fraction over left-center interface VfracL = options.VfracL # default is 0.5 print 'Inelastica: Voltage fraction over left-center interface: VfracL =', VfracL # Set up electrodes and device Greens function elecL = NEGF.ElectrodeSelfEnergy(options.fnL, options.NA1L, options.NA2L, options.voltage * VfracL) elecL.scaling = options.scaleSigL elecL.semiinf = options.semiinfL elecR = NEGF.ElectrodeSelfEnergy(options.fnR, options.NA1R, options.NA2R, options.voltage * (VfracL - 1.)) elecR.scaling = options.scaleSigR elecR.semiinf = options.semiinfR # Read phonons NCfile = NC4.Dataset(options.PhononNetCDF, 'r') print 'Inelastica: Reading ', options.PhononNetCDF hw = NCfile.variables['hw'][:] # Work with GFs etc for positive (V>0: \mu_L>\mu_R) and negative (V<0: \mu_L<\mu_R) bias voltages GFp = NEGF.GF(options.TSHS, elecL, elecR, Bulk=options.UseBulk, DeviceAtoms=options.DeviceAtoms, BufferAtoms=options.buffer) # Prepare lists for various trace factors #GF.dGnout = [] #GF.dGnin = [] GFp.P1T = N.zeros(len(hw), N.float) # M.A.M.A (total e-h damping) GFp.P2T = N.zeros(len(hw), N.float) # M.AL.M.AR (emission) GFp.ehDampL = N.zeros(len(hw), N.float) # M.AL.M.AL (L e-h damping) GFp.ehDampR = N.zeros(len(hw), N.float) # M.AR.M.AR (R e-h damping) GFp.nHT = N.zeros(len(hw), N.float) # non-Hilbert/Isym factor GFp.HT = N.zeros(len(hw), N.float) # Hilbert/Iasym factor GFp.dIel = N.zeros(len(hw), N.float) GFp.dIinel = N.zeros(len(hw), N.float) GFp.dSel = N.zeros(len(hw), N.float) GFp.dSinel = N.zeros(len(hw), N.float) # GFm = NEGF.GF(options.TSHS, elecL, elecR, Bulk=options.UseBulk, DeviceAtoms=options.DeviceAtoms, BufferAtoms=options.buffer) GFm.P1T = N.zeros(len(hw), N.float) # M.A.M.A (total e-h damping) GFm.P2T = N.zeros(len(hw), N.float) # M.AL.M.AR (emission) GFm.ehDampL = N.zeros(len(hw), N.float) # M.AL.M.AL (L e-h damping) GFm.ehDampR = N.zeros(len(hw), N.float) # M.AR.M.AR (R e-h damping) GFm.nHT = N.zeros(len(hw), N.float) # non-Hilbert/Isym factor GFm.HT = N.zeros(len(hw), N.float) # Hilbert/Iasym factor GFm.dIel = N.zeros(len(hw), N.float) GFm.dIinel = N.zeros(len(hw), N.float) GFm.dSel = N.zeros(len(hw), N.float) GFm.dSinel = N.zeros(len(hw), N.float) # Calculate transmission at Fermi level GFp.calcGF(options.energy + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) L = options.bufferL # Pad lasto with zeroes to enable basis generation... lasto = N.zeros((GFp.HS.nua + L + 1, ), N.int) lasto[L:] = GFp.HS.lasto basis = SIO.BuildBasis(options.fn, options.DeviceAtoms[0] + L, options.DeviceAtoms[1] + L, lasto) basis.ii -= L TeF = MM.trace(GFp.TT).real GFp.TeF = TeF GFm.TeF = TeF # Check consistency of PHrun vs TSrun inputs IntegrityCheck(options, GFp, basis, NCfile) # Calculate trace factors one mode at a time print 'Inelastica: LOEscale =', options.LOEscale if options.LOEscale == 0.0: # LOEscale=0.0 => Original LOE-WBA method, PRB 72, 201101(R) (2005) [cond-mat/0505473]. GFp.calcGF(options.energy + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) GFm.calcGF(options.energy + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) for ihw in (hw > options.modeCutoff).nonzero()[0]: calcTraces(options, GFp, GFm, basis, NCfile, ihw) calcTraces(options, GFm, GFp, basis, NCfile, ihw) writeFGRrates(options, GFp, hw, NCfile) else: # LOEscale=1.0 => Generalized LOE, PRB 89, 081405(R) (2014) [arXiv:1312.7625] for ihw in (hw > options.modeCutoff).nonzero()[0]: GFp.calcGF(options.energy + hw[ihw] * options.LOEscale * VfracL + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) GFm.calcGF(options.energy + hw[ihw] * options.LOEscale * (VfracL - 1.) + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) calcTraces(options, GFp, GFm, basis, NCfile, ihw) if VfracL != 0.5: GFp.calcGF(options.energy - hw[ihw] * options.LOEscale * (VfracL - 1.) + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) GFm.calcGF(options.energy - hw[ihw] * options.LOEscale * VfracL + options.eta * 1.0j, options.kpoint[0:2], ispin=options.iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) calcTraces(options, GFm, GFp, basis, NCfile, ihw) # Multiply traces with voltage-dependent functions data = calcIETS(options, GFp, GFm, basis, hw) NCfile.close() NEGF.SavedSig.close() CF.PrintMainFooter('Inelastica') return data
from Inelastica import MiscMath as mm import numpy as N from Inelastica import WriteXMGR as XMGR import scipy.special as ss # We try a simple Gaussian function: x0 = 50 x = N.linspace(-x0, x0, 1e5) gauss = -1j * N.pi**-0.5 * N.exp(-x**2) # Output from our Hilbert function: Hg, ker = mm.Hilbert(gauss) # We can compare with this: # Hilbert transform of a Gaussian function is related to Faddeva/w(z) functions: # https://en.wikipedia.org/wiki/Dawson_function ex = -1j * N.pi**0.5 * ss.wofz(x) # Collect data in a plot data1 = XMGR.XYset(x, gauss.imag, legend='Gauss', Lwidth=2) data2 = XMGR.XYset(x, ex.real, legend='Faddeva', Lwidth=2) data3 = XMGR.XYset(x, N.pi * Hg.imag, legend='Inelastica', Lwidth=1) graph1 = XMGR.Graph(data1, data2, data3) graph1.SetXaxis(autoscale=True) graph1.SetYaxis(autoscale=True) graph1.ShowLegend() # Compute error/difference between the two methods: err = ex.real - N.pi * Hg.imag data4 = XMGR.XYset(x, err, legend='Difference', Lwidth=2) graph2 = XMGR.Graph(data4)
def main(options): """ Main routine to compute elastic transmission probabilities etc. Parameters ---------- options : an ``options`` instance """ CF.CreatePipeOutput(options.DestDir + '/' + options.Logfile) VC.OptionsCheck(options, 'pyTBT') CF.PrintMainHeader('pyTBT', options) # K-points if options.Gk1 > 1: Nk1, t1 = options.Gk1, 'GK' else: Nk1, t1 = options.Nk1, 'LIN' if options.Gk2 > 1: Nk2, t2 = options.Gk2, 'GK' else: Nk2, t2 = options.Nk2, 'LIN' # Generate full k-mesh: mesh = Kmesh.kmesh(Nk1, Nk2, Nk3=1, meshtype=[t1, t2, 'LIN'], invsymmetry=not options.skipsymmetry) mesh.mesh2file( '%s/%s.%ix%i.mesh' % (options.DestDir, options.systemlabel, mesh.Nk[0], mesh.Nk[1])) # Setup self-energies and device GF elecL = NEGF.ElectrodeSelfEnergy(options.fnL, options.NA1L, options.NA2L, options.voltage / 2.) elecL.scaling = options.scaleSigL elecL.semiinf = options.semiinfL elecR = NEGF.ElectrodeSelfEnergy(options.fnR, options.NA1R, options.NA2R, -options.voltage / 2.) elecR.scaling = options.scaleSigR elecR.semiinf = options.semiinfR DevGF = NEGF.GF(options.TSHS, elecL, elecR, Bulk=options.UseBulk, DeviceAtoms=options.DeviceAtoms, BufferAtoms=options.buffer) nspin = DevGF.HS.nspin # k-sample only self-energies? if options.singlejunction: elecL.mesh = mesh mesh = Kmesh.kmesh(3, 3, 1) if options.dos: DOSL = N.zeros((nspin, len(options.Elist), DevGF.nuo), N.float) DOSR = N.zeros((nspin, len(options.Elist), DevGF.nuo), N.float) # MPSH projections? MPSHL = N.zeros((nspin, len(options.Elist), DevGF.nuo), N.float) MPSHR = N.zeros((nspin, len(options.Elist), DevGF.nuo), N.float) # evaluate eigenstates at Gamma import scipy.linalg as SLA DevGF.setkpoint(N.zeros(2)) ev0, es0 = SLA.eigh(DevGF.H, DevGF.S) print 'MPSH eigenvalues:', ev0 #print 'MPSH eigenvector normalizations:',N.diag(MM.mm(MM.dagger(es0),DevGF.S,es0)).real # right # Loop over spin for iSpin in range(nspin): # initialize transmission and shot noise arrays Tkpt = N.zeros((len(options.Elist), mesh.NNk, options.numchan + 1), N.float) SNkpt = N.zeros((len(options.Elist), mesh.NNk, options.numchan + 1), N.float) # prepare output files outFile = options.DestDir + '/%s.%ix%i' % (options.systemlabel, mesh.Nk[0], mesh.Nk[1]) if nspin < 2: thisspinlabel = outFile else: thisspinlabel = outFile + ['.UP', '.DOWN'][iSpin] fo = open(thisspinlabel + '.AVTRANS', 'write') fo.write('# Nk1(%s)=%i Nk2(%s)=%i eta=%.2e etaLead=%.2e\n' % (mesh.type[0], mesh.Nk[0], mesh.type[1], mesh.Nk[1], options.eta, options.etaLead)) fo.write('# E Ttot(E) Ti(E)(i=1-%i) RelErrorTtot(E)\n' % options.numchan) foSN = open(thisspinlabel + '.AVNOISE', 'write') foSN.write('# Nk1(%s)=%i Nk2(%s)=%i eta=%.2e etaLead=%.2e\n' % (mesh.type[0], mesh.Nk[0], mesh.type[1], mesh.Nk[1], options.eta, options.etaLead)) foSN.write('# E SNtot(E) SNi(E)(i=1-%i)\n' % options.numchan) foFF = open(thisspinlabel + '.FANO', 'write') foFF.write('# Nk1(%s)=%i Nk2(%s)=%i eta=%.2e etaLead=%.2e\n' % (mesh.type[0], mesh.Nk[0], mesh.type[1], mesh.Nk[1], options.eta, options.etaLead)) foFF.write('# E Fano factor \n') # Loop over energy for ie, ee in enumerate(options.Elist): Tavg = N.zeros((options.numchan + 1, len(mesh.w)), N.float) SNavg = N.zeros((options.numchan + 1, len(mesh.w)), N.float) AavL = N.zeros((DevGF.nuo, DevGF.nuo), N.complex) AavR = N.zeros((DevGF.nuo, DevGF.nuo), N.complex) # Loops over k-points for ik in range(mesh.NNk): DevGF.calcGF(ee + options.eta * 1.0j, mesh.k[ik, :2], ispin=iSpin, etaLead=options.etaLead, useSigNCfiles=options.signc, SpectralCutoff=options.SpectralCutoff) # Transmission and shot noise T, SN = DevGF.calcTEIG(options.numchan) for iw in range(len(mesh.w)): Tavg[:, iw] += T * mesh.w[iw, ik] SNavg[:, iw] += SN * mesh.w[iw, ik] Tkpt[ie, ik] = T SNkpt[ie, ik] = SN # DOS calculation: if options.dos: if options.SpectralCutoff > 0.0: AavL += mesh.w[0, ik] * MM.mm(DevGF.AL.L, DevGF.AL.R, DevGF.S) AavR += mesh.w[0, ik] * MM.mm(DevGF.AR.L, DevGF.AR.R, DevGF.S) else: AavL += mesh.w[0, ik] * MM.mm(DevGF.AL, DevGF.S) AavR += mesh.w[0, ik] * MM.mm(DevGF.AR, DevGF.S) # Print calculated quantities err = (N.abs(Tavg[0, 0] - Tavg[0, 1]) + N.abs(Tavg[0, 0] - Tavg[0, 2])) / 2 relerr = err / Tavg[0, 0] print 'ispin= %i, e= %.4f, Tavg= %.8f, RelErr= %.1e' % ( iSpin, ee, Tavg[0, 0], relerr) transline = '\n%.10f ' % ee noiseline = '\n%.10f ' % ee for ichan in range(options.numchan + 1): if ichan == 0: transline += '%.8e ' % Tavg[ichan, 0] noiseline += '%.8e ' % SNavg[ichan, 0] else: transline += '%.4e ' % Tavg[ichan, 0] noiseline += '%.4e ' % SNavg[ichan, 0] transline += '%.2e ' % relerr fo.write(transline) foSN.write(noiseline) foFF.write('\n%.10f %.8e' % (ee, SNavg[0, 0] / Tavg[0, 0])) # Partial density of states: if options.dos: DOSL[iSpin, ie, :] += N.diag(AavL).real / (2 * N.pi) DOSR[iSpin, ie, :] += N.diag(AavR).real / (2 * N.pi) MPSHL[iSpin, ie, :] += N.diag(MM.mm(MM.dagger(es0), AavL, es0)).real / (2 * N.pi) MPSHR[iSpin, ie, :] += N.diag(MM.mm(MM.dagger(es0), AavR, es0)).real / (2 * N.pi) print 'ispin= %i, e= %.4f, DOSL= %.4f, DOSR= %.4f' % ( iSpin, ee, N.sum(DOSL[iSpin, ie, :]), N.sum(DOSR[iSpin, ie, :])) fo.write('\n') fo.close() foSN.write('\n') foSN.close() foFF.write('\n') foFF.close() # Write k-point-resolved transmission fo = open(thisspinlabel + '.TRANS', 'write') for ik in range(mesh.NNk): w = mesh.w[:, ik] fo.write('\n\n# k = %f, %f w = %f %f %f %f' % (mesh.k[ik, 0], mesh.k[ik, 1], w[0], w[1], w[2], w[3])) for ie, ee in enumerate(options.Elist): transline = '\n%.10f ' % ee for ichan in range(options.numchan + 1): if ichan == 0: transline += '%.8e ' % Tkpt[ie, ik, ichan] else: transline += '%.4e ' % Tkpt[ie, ik, ichan] fo.write(transline) fo.write('\n') fo.close() # Write k-point-resolved shot noise fo = open(thisspinlabel + '.NOISE', 'write') for ik in range(mesh.NNk): w = mesh.w[:, ik] fo.write('\n\n# k = %f, %f w = %f %f %f %f' % (mesh.k[ik, 0], mesh.k[ik, 1], w[0], w[1], w[2], w[3])) for ie, ee in enumerate(options.Elist): noiseline = '\n%.10f ' % ee for ichan in range(options.numchan + 1): if ichan == 0: noiseline += '%.8e ' % SNkpt[ie, ik, ichan] else: noiseline += '%.4e ' % SNkpt[ie, ik, ichan] fo.write(noiseline) fo.write('\n') fo.close() # End loop over spin NEGF.SavedSig.close() # Make sure saved Sigma is written to file if options.dos: # Read basis L = options.bufferL # Pad lasto with zeroes to enable basis generation... lasto = N.zeros((DevGF.HS.nua + L + 1, ), N.int) lasto[L:] = DevGF.HS.lasto basis = SIO.BuildBasis(options.fn, 1 + L, DevGF.HS.nua + L, lasto) basis.ii -= L WritePDOS(outFile + '.PDOS.gz', options, DevGF, DOSL + DOSR, basis) WritePDOS(outFile + '.PDOSL.gz', options, DevGF, DOSL, basis) WritePDOS(outFile + '.PDOSR.gz', options, DevGF, DOSR, basis) WriteMPSH(outFile + '.MPSH.gz', options, DevGF, MPSHL + MPSHR, ev0) WriteMPSH(outFile + '.MPSHL.gz', options, DevGF, MPSHL, ev0) WriteMPSH(outFile + '.MPSHR.gz', options, DevGF, MPSHR, ev0) CF.PrintMainFooter('pyTBT')
def main(options): CF.CreatePipeOutput(options.DestDir + '/' + options.Logfile) #VC.OptionsCheck(options,'Phonons') CF.PrintMainHeader('Bandstructures', options) try: fdf = glob.glob(options.onlyTSdir + '/RUN.fdf') TSrun = True except: fdf = glob.glob(options.FCwildcard + '/RUN.fdf') # This should be made an input flag TSrun = False SCDM = Supercell_DynamicalMatrix(fdf, TSrun) # Write high-symmetry path WritePath(options.DestDir + '/symmetry-path', SCDM.Sym.path, options.steps) # Write mesh k1, k2, k3 = ast.literal_eval(options.mesh) rvec = 2 * N.pi * N.array([SCDM.Sym.b1, SCDM.Sym.b2, SCDM.Sym.b3]) import Inelastica.physics.mesh as Kmesh # Full mesh kmesh = Kmesh.kmesh(2**k1, 2**k2, 2**k3, meshtype=['LIN', 'LIN', 'LIN'], invsymmetry=False) WriteKpoints(options.DestDir + '/mesh_%ix%ix%i' % tuple(kmesh.Nk), N.dot(kmesh.k, rvec)) # Mesh reduced by inversion symmetry kmesh = Kmesh.kmesh(2**k1, 2**k2, 2**k3, meshtype=['LIN', 'LIN', 'LIN'], invsymmetry=True) WriteKpoints(options.DestDir + '/mesh_%ix%ix%i_invsym' % tuple(kmesh.Nk), N.dot(kmesh.k, rvec)) # Evaluate electron k-points if options.kfile: # Prepare Hamiltonian etc in Gamma for whole supercell natoms = SIO.GetFDFlineWithDefault(fdf[0], 'NumberOfAtoms', int, -1, 'Error') SCDM.PrepareGradients(options.onlySdir, N.array([0., 0., 0.]), 1, natoms, AbsEref=False, atype=N.complex, TSrun=TSrun) SCDM.nao = SCDM.h0.shape[-1] SCDM.FirstOrb = SCDM.OrbIndx[0][0] # First atom = 1 SCDM.LastOrb = SCDM.OrbIndx[SCDM.Sym.basis.NN - 1][1] # Last atom = Sym.NN SCDM.rednao = SCDM.LastOrb + 1 - SCDM.FirstOrb # Read kpoints kpts, dk, klabels, kticks = ReadKpoints(options.kfile) if klabels: # Only write ascii if labels exist WriteKpoints(options.DestDir + '/kpoints', kpts, klabels) # Prepare netcdf ncfn = options.DestDir + '/Electrons.nc' ncf = NC4.Dataset(ncfn, 'w') # Grid ncf.createDimension('gridpts', len(kpts)) ncf.createDimension('vector', 3) grid = ncf.createVariable('grid', 'd', ('gridpts', 'vector')) grid[:] = kpts grid.units = '1/Angstrom' # Geometry ncf.createDimension('atoms', SCDM.Sym.basis.NN) xyz = ncf.createVariable('xyz', 'd', ('atoms', 'vector')) xyz[:] = SCDM.Sym.basis.xyz xyz.units = 'Angstrom' pbc = ncf.createVariable('pbc', 'd', ('vector', 'vector')) pbc.units = 'Angstrom' pbc[:] = [SCDM.Sym.a1, SCDM.Sym.a2, SCDM.Sym.a3] rvec1 = ncf.createVariable('rvec', 'd', ('vector', 'vector')) rvec1.units = '1/Angstrom (incl. factor 2pi)' rvec1[:] = rvec ncf.sync() # Loop over kpoints for i, k in enumerate(kpts): if i < 100: # Print only for the first 100 points ev, evec = SCDM.ComputeElectronStates(k, verbose=True, TSrun=TSrun) else: ev, evec = SCDM.ComputeElectronStates(k, verbose=False, TSrun=TSrun) # otherwise something simple if i % 100 == 0: print '%i out of %i k-points computed' % (i, len(kpts)) if i == 0: ncf.createDimension('nspin', SCDM.nspin) ncf.createDimension('orbs', SCDM.rednao) if options.nbands and options.nbands < SCDM.rednao: nbands = options.nbands else: nbands = SCDM.rednao ncf.createDimension('bands', nbands) evals = ncf.createVariable('eigenvalues', 'd', ('gridpts', 'nspin', 'bands')) evals.units = 'eV' evecsRe = ncf.createVariable( 'eigenvectors.re', 'd', ('gridpts', 'nspin', 'orbs', 'bands')) evecsIm = ncf.createVariable( 'eigenvectors.im', 'd', ('gridpts', 'nspin', 'orbs', 'bands')) # Check eigenvectors print 'SupercellPhonons: Checking eigenvectors at', k for j in range(SCDM.nspin): ev2 = N.diagonal( MM.mm(MM.dagger(evec[j]), SCDM.h0_k[j], evec[j])) print ' ... spin %i: Allclose=' % j, N.allclose(ev[j], ev2, atol=1e-5, rtol=1e-3) ncf.sync() # Write to NetCDF evals[i, :] = ev[:, :nbands] evecsRe[i, :] = evec[:, :, :nbands].real evecsIm[i, :] = evec[:, :, :nbands].imag ncf.sync() # Include basis orbitals in netcdf file if SCDM.Sym.basis.NN == len(SCDM.OrbIndx): lasto = N.zeros(SCDM.Sym.basis.NN + 1, N.float) lasto[:SCDM.Sym.basis.NN] = SCDM.OrbIndx[:SCDM.Sym.basis.NN, 0] lasto[SCDM.Sym.basis.NN] = SCDM.OrbIndx[SCDM.Sym.basis.NN - 1, 1] + 1 else: lasto = SCDM.OrbIndx[:SCDM.Sym.basis.NN + 1, 0] orbbasis = SIO.BuildBasis(fdf[0], 1, SCDM.Sym.basis.NN, lasto) # Note that the above basis is for the geometry with an atom FC-moved in z. #print dir(orbbasis) #print orbbasis.xyz # Hence, this is not the correct geometry of the basis atoms! center = ncf.createVariable('orbcenter', 'i', ('orbs', )) center[:] = N.array(orbbasis.ii - 1, dtype='int32') center.description = 'Atom index (counting from 0) of the orbital center' nn = ncf.createVariable('N', 'i', ('orbs', )) nn[:] = N.array(orbbasis.N, dtype='int32') ll = ncf.createVariable('L', 'i', ('orbs', )) ll[:] = N.array(orbbasis.L, dtype='int32') mm = ncf.createVariable('M', 'i', ('orbs', )) mm[:] = N.array(orbbasis.M, dtype='int32') # Cutoff radius and delta Rc = ncf.createVariable('Rc', 'd', ('orbs', )) Rc[:] = orbbasis.coff Rc.units = 'Angstrom' delta = ncf.createVariable('delta', 'd', ('orbs', )) delta[:] = orbbasis.delta delta.units = 'Angstrom' # Radial components of the orbitals ntb = len(orbbasis.orb[0]) ncf.createDimension('ntb', ntb) rii = ncf.createVariable('rii', 'd', ('orbs', 'ntb')) rii[:] = N.outer(orbbasis.delta, N.arange(ntb)) rii.units = 'Angstrom' radialfct = ncf.createVariable('radialfct', 'd', ('orbs', 'ntb')) radialfct[:] = orbbasis.orb # Sort eigenvalues to connect crossing bands? if options.sorting: for i in range(SCDM.nspin): evals[:, i, :] = SortBands(evals[:, i, :]) # Produce nice plots if labels exist if klabels: if SCDM.nspin == 1: PlotElectronBands(options.DestDir + '/Electrons.agr', dk, evals[:, 0, :], kticks) elif SCDM.nspin == 2: PlotElectronBands(options.DestDir + '/Electrons.UP.agr', dk, evals[:, 0, :], kticks) PlotElectronBands(options.DestDir + '/Electrons.DOWN.agr', dk, evals[:, 1, :], kticks) ncf.close() if TSrun: # only electronic calculation return SCDM.Sym.path # Compute phonon eigenvalues if options.qfile: SCDM.SymmetrizeFC(options.radius) SCDM.SetMasses() qpts, dq, qlabels, qticks = ReadKpoints(options.qfile) if qlabels: # Only write ascii if labels exist WriteKpoints(options.DestDir + '/qpoints', qpts, qlabels) # Prepare netcdf ncfn = options.DestDir + '/Phonons.nc' ncf = NC4.Dataset(ncfn, 'w') # Grid ncf.createDimension('gridpts', len(qpts)) ncf.createDimension('vector', 3) grid = ncf.createVariable('grid', 'd', ('gridpts', 'vector')) grid[:] = qpts grid.units = '1/Angstrom' # Geometry ncf.createDimension('atoms', SCDM.Sym.basis.NN) xyz = ncf.createVariable('xyz', 'd', ('atoms', 'vector')) xyz[:] = SCDM.Sym.basis.xyz xyz.units = 'Angstrom' pbc = ncf.createVariable('pbc', 'd', ('vector', 'vector')) pbc.units = 'Angstrom' pbc[:] = [SCDM.Sym.a1, SCDM.Sym.a2, SCDM.Sym.a3] rvec1 = ncf.createVariable('rvec', 'd', ('vector', 'vector')) rvec1.units = '1/Angstrom (incl. factor 2pi)' rvec1[:] = rvec ncf.sync() # Loop over q for i, q in enumerate(qpts): if i < 100: # Print only for the first 100 points hw, U = SCDM.ComputePhononModes_q(q, verbose=True) else: hw, U = SCDM.ComputePhononModes_q(q, verbose=False) # otherwise something simple if i % 100 == 0: print '%i out of %i q-points computed' % (i, len(qpts)) if i == 0: ncf.createDimension('bands', len(hw)) ncf.createDimension('displ', len(hw)) evals = ncf.createVariable('eigenvalues', 'd', ('gridpts', 'bands')) evals.units = 'eV' evecsRe = ncf.createVariable('eigenvectors.re', 'd', ('gridpts', 'bands', 'displ')) evecsIm = ncf.createVariable('eigenvectors.im', 'd', ('gridpts', 'bands', 'displ')) # Check eigenvectors print 'SupercellPhonons.Checking eigenvectors at', q tmp = MM.mm(N.conjugate(U), SCDM.FCtilde, N.transpose(U)) const = PC.hbar2SI * (1e20 / (PC.eV2Joule * PC.amu2kg))**0.5 hw2 = const * N.diagonal(tmp)**0.5 # Units in eV print ' ... Allclose=', N.allclose(hw, N.absolute(hw2), atol=1e-5, rtol=1e-3) ncf.sync() evals[i] = hw evecsRe[i] = U.real evecsIm[i] = U.imag ncf.sync() # Sort eigenvalues to connect crossing bands? if options.sorting: evals = SortBands(evals) # Produce nice plots if labels exist if qlabels: PlotPhononBands(options.DestDir + '/Phonons.agr', dq, N.array(evals[:]), qticks) ncf.close() # Compute e-ph couplings if options.kfile and options.qfile: SCDM.ReadGradients(AbsEref=False) ncf = NC4.Dataset(options.DestDir + '/EPH.nc', 'w') ncf.createDimension('kpts', len(kpts)) ncf.createDimension('qpts', len(qpts)) ncf.createDimension('modes', len(hw)) ncf.createDimension('nspin', SCDM.nspin) ncf.createDimension('bands', SCDM.rednao) ncf.createDimension('vector', 3) kgrid = ncf.createVariable('kpts', 'd', ('kpts', 'vector')) kgrid[:] = kpts qgrid = ncf.createVariable('qpts', 'd', ('qpts', 'vector')) qgrid[:] = qpts evalfkq = ncf.createVariable('evalfkq', 'd', ('kpts', 'qpts', 'nspin', 'bands')) # First (second) band index n (n') is the initial (final) state, i.e., # Mkq(k,q,mode,spin,n,n') := < n',k+q | dV_q(mode) | n,k > MkqAbs = ncf.createVariable( 'Mkqabs', 'd', ('kpts', 'qpts', 'modes', 'nspin', 'bands', 'bands')) GkqAbs = ncf.createVariable( 'Gkqabs', 'd', ('kpts', 'qpts', 'modes', 'nspin', 'bands', 'bands')) ncf.sync() # Loop over k-points for i, k in enumerate(kpts): kpts[i] = k # Compute initial electronic states evi, eveci = SCDM.ComputeElectronStates(k, verbose=True) # Loop over q-points for j, q in enumerate(qpts): # Compute phonon modes hw, U = SCDM.ComputePhononModes_q(q, verbose=True) # Compute final electronic states evf, evecf = SCDM.ComputeElectronStates(k + q, verbose=True) evalfkq[i, j, :] = evf # Compute electron-phonon couplings m, g = SCDM.ComputeEPHcouplings_kq( k, q) # (modes,nspin,bands,bands) # Data to file # M (modes,spin,i,l) = m(modes,k,j) init(i,j) final(k,l) # 0 1 2 0,1 0 1 # ^-------^ # ^----------------------^ for ispin in range(SCDM.nspin): evecfd = MM.dagger(evecf[ispin]) # (bands,bands) M = N.tensordot(N.tensordot(m[:, ispin], eveci[ispin], axes=[2, 0]), evecfd, axes=[1, 1]) G = N.tensordot(N.tensordot(g[:, ispin], eveci[ispin], axes=[2, 0]), evecfd, axes=[1, 1]) MkqAbs[i, j, :, ispin] = N.absolute(M) GkqAbs[i, j, :, ispin] = N.absolute(G) ncf.sync() ncf.close() return SCDM.Sym.path
def calcWF(options, geom, basis, Y): """ Calculate wavefunction, returns: YY : complex wavefunction on regular grid dstep : stepsize origo : vector nx, ny, nz : number of grid points """ xyz=N.array(geom.xyz[options.DeviceAtoms[0]-1:options.DeviceAtoms[1]]) # Size of cube xmin, xmax = min(xyz[:, 0])-5.0, max(xyz[:, 0])+5.0 ymin, ymax = min(xyz[:, 1])-5.0, max(xyz[:, 1])+5.0 zmin, zmax = min(xyz[:, 2])-5.0, max(xyz[:, 2])+5.0 xl, yl, zl = xmax-xmin, ymax-ymin, zmax-zmin dx, dy, dz = options.res, options.res, options.res nx, ny, nz = int(xl/dx)+1, int(yl/dy)+1, int(zl/dz)+1 origo = N.array([xmin, ymin, zmin], N.float) # Def cube YY=N.zeros((nx, ny, nz), N.complex) rx=N.array(range(nx), N.float)*dx+origo[0] ry=N.array(range(ny), N.float)*dy+origo[1] rz=N.array(range(nz), N.float)*dz+origo[2] for ii, Yval in enumerate(Y): if ii>0:# and ii%(int(len(Y)/10))==0: SIO.printDone(ii, len(Y), 'Wavefunction') rax, ray, raz=basis.xyz[ii, 0], basis.xyz[ii, 1], basis.xyz[ii, 2] # Only calulate in subset ixmin, ixmax = int((rax-origo[0]-basis.coff[ii])/dx), \ int((rax-origo[0]+basis.coff[ii])/dx) iymin, iymax = int((ray-origo[1]-basis.coff[ii])/dy), \ int((ray-origo[1]+basis.coff[ii])/dy) izmin, izmax = int((raz-origo[2]-basis.coff[ii])/dz), \ int((raz-origo[2]+basis.coff[ii])/dz) ddx, ddy, ddz=rx[ixmin:ixmax]-rax, ry[iymin:iymax]-ray, rz[izmin:izmax]-raz dr=N.sqrt(MM.outerAdd(ddx*ddx, ddy*ddy, ddz*ddz)) drho=N.sqrt(MM.outerAdd(ddx*ddx, ddy*ddy, 0*ddz)) imax=(basis.coff[ii]-2*basis.delta[ii])/basis.delta[ii] ri=dr/basis.delta[ii] ri=N.where(ri<imax, ri, imax) ri=ri.astype(N.int) costh = MM.outerAdd(0*ddx, 0*ddy, ddz)/dr cosfi, sinfi = MM.outerAdd(ddx, 0*ddy, 0*ddz)/drho, MM.outerAdd(0*ddx, ddy, 0*ddz)/drho # Numpy has changed the choose function to crap! RR=N.take(basis.orb[ii], ri) # Calculate spherical harmonics l = basis.L[ii] m = basis.M[ii] if l==3: print 'f-shell : l=%i, m=%i (NOT TESTED!!)'%(l, m) thisSphHar = MM.sphericalHarmonics(l, m, costh, sinfi, cosfi) YY[ixmin:ixmax, iymin:iymax, izmin:izmax]=YY[ixmin:ixmax, iymin:iymax, izmin:izmax]+\ RR*thisSphHar*Yval print "Wave function norm on real space grid:", N.sum(YY.conjugate()*YY)*dx*dy*dz return YY, options.res, origo, nx, ny, nz
def calcGF(self, ee, kpoint, ispin=0, etaLead=0.0, useSigNCfiles=False, SpectralCutoff=0.0): "Calculate GF etc at energy ee and 2d k-point" nuo, nuoL, nuoR = self.nuo, self.nuoL, self.nuoR nuo0, nuoL0, nuoR0 = self.nuo0, self.nuoL0, self.nuoR0 FoldedL, FoldedR = self.FoldedL, self.FoldedR devSt, devEnd = self.DeviceOrbs[0], self.DeviceOrbs[1] # Determine whether electrode self-energies should be k-sampled or not try: mesh = self.elecL.mesh # a mesh was attached except: mesh = False # Calculate electrode self-energies if mesh: try: self.SigAvg # Averaged self-energies exist except: self.SigAvg = [False, -1] if self.SigAvg[0] == ee and self.SigAvg[1] == ispin: # We have already the averaged self-energies print 'NEGF: Reusing sampled electrode self-energies', mesh.Nk, mesh.type, 'for ispin= %i e= %f' % ( ispin, ee) else: # k-sampling performed over folded electrode self-energies print 'NEGF: Sampling electrode self-energies', mesh.Nk, mesh.type, 'for ispin= %i e= %f' % ( ispin, ee) self.calcSigLR(ee, mesh.k[0, :2], ispin, etaLead, useSigNCfiles, SpectralCutoff) AvgSigL = mesh.w[0, 0] * self.SigL AvgSigR = mesh.w[0, 0] * self.SigR for i in range(1, len(mesh.k)): self.calcSigLR(ee, mesh.k[i, :2], ispin, etaLead, useSigNCfiles, SpectralCutoff) AvgSigL += mesh.w[0, i] * self.SigL AvgSigR += mesh.w[0, i] * self.SigR # We now simply continue with the averaged self-energies self.SigL = AvgSigL self.SigR = AvgSigR self.SigAvg = [ee, ispin] else: # We sample k-points the usual way self.calcSigLR(ee, kpoint, ispin, etaLead, useSigNCfiles) # Ready to calculate Gr self.setkpoint(kpoint, ispin) eSmH = ee * self.S - self.H if FoldedL: eSmH[0:nuoL, 0:nuoL] = eSmH[0:nuoL, 0:nuoL] - self.SigL else: if self.Bulk: eSmH[0:nuoL, 0:nuoL] = self.SigL # SGF^1 else: eSmH[0:nuoL, 0:nuoL] = eSmH[0:nuoL, 0:nuoL] - self.SigL if FoldedR: eSmH[nuo - nuoR:nuo, nuo - nuoR:nuo] = eSmH[nuo - nuoR:nuo, nuo - nuoR:nuo] - self.SigR else: if self.Bulk: eSmH[nuo - nuoR:nuo, nuo - nuoR:nuo] = self.SigR # SGF^1 else: eSmH[nuo - nuoR:nuo, nuo - nuoR:nuo] = eSmH[nuo - nuoR:nuo, nuo - nuoR:nuo] - self.SigR self.Gr = LA.inv(eSmH) self.Ga = MM.dagger(self.Gr) # Calculate spectral functions if SpectralCutoff > 0.0: self.AL = MM.SpectralMatrix(MM.mm(self.Gr[:, 0:nuoL], self.GamL, self.Ga[0:nuoL, :]), cutoff=SpectralCutoff) tmp = MM.mm(self.GamL, self.Gr[0:nuoL, :]) self.ALT = MM.SpectralMatrix(MM.mm(self.Ga[:, 0:nuoL], tmp), cutoff=SpectralCutoff) self.AR = MM.SpectralMatrix(MM.mm(self.Gr[:, nuo - nuoR:nuo], self.GamR, self.Ga[nuo - nuoR:nuo, :]), cutoff=SpectralCutoff) self.ARGLG = MM.mm(self.AR.L, self.AR.R[:, 0:nuoL], tmp) self.A = self.AL + self.AR # transmission matrix AL.GamR self.TT = MM.mm(self.AL.R[:, nuo - nuoR:nuo], self.GamR, self.AL.L[nuo - nuoR:nuo, :]) else: self.AL = MM.mm(self.Gr[:, 0:nuoL], self.GamL, self.Ga[0:nuoL, :]) tmp = MM.mm(self.GamL, self.Gr[0:nuoL, :]) self.ALT = MM.mm(self.Ga[:, 0:nuoL], tmp) self.AR = MM.mm(self.Gr[:, nuo - nuoR:nuo], self.GamR, self.Ga[nuo - nuoR:nuo, :]) self.ARGLG = MM.mm(self.AR[:, 0:nuoL], tmp) self.A = self.AL + self.AR # transmission matrix AL.GamR self.TT = MM.mm(self.AL[nuo - nuoR:nuo, nuo - nuoR:nuo], self.GamR) print 'NEGF.calcGF: Shape of transmission matrix (TT):', self.TT.shape print 'NEGF.calcGF: Energy and total transmission Tr[TT].real:', ee, N.trace( self.TT).real # Write also the Gammas in the full space of Gr/Ga/A # (needed for the inelastic shot noise) self.GammaL = N.zeros(self.Gr.shape, N.complex) self.GammaL[0:nuoL, 0:nuoL] = self.GamL self.GammaR = N.zeros(self.Gr.shape, N.complex) self.GammaR[nuo - nuoR:nuo, nuo - nuoR:nuo] = self.GamR
def calcSigLR(self, ee, kpoint, ispin=0, etaLead=0.0, useSigNCfiles=False, SpectralCutoff=0.0): """ Calculate (folded) self-energy at energy ee and 2d k-point Uses SpectralMatrix format for the spectralfunction matrices, see MiscMath, if cutoff>0.0 """ nuoL, nuoR = self.nuoL, self.nuoR nuo0, nuoL0, nuoR0 = self.nuo0, self.nuoL0, self.nuoR0 FoldedL, FoldedR = self.FoldedL, self.FoldedR devSt, devEnd = self.DeviceOrbs[0], self.DeviceOrbs[1] # Calculate Sigma without folding self.setkpoint(kpoint, ispin) SigL0 = self.elecL.getSig(ee, kpoint, left=True, Bulk=self.Bulk, ispin=ispin, etaLead=etaLead, useSigNCfiles=useSigNCfiles) SigR0 = self.elecR.getSig(ee, kpoint, left=False, Bulk=self.Bulk, ispin=ispin, etaLead=etaLead, useSigNCfiles=useSigNCfiles) if FoldedL: # Fold down from nuoL0 to the device region # A11 A12 g11 g12 I 0 # A21 A22 * g21 g22 = 0 I -> # g22 = (A22-A21.A11^-1.A12)^-1 -> # Sigma = A21.A11^-1.A12 (tau=A12) devEndL = self.devEndL # Do folding eSmH = ee * self.S0 - self.H0 eSmHmS = eSmH[0:devEndL, 0:devEndL].copy() if self.Bulk: eSmHmS[0:nuoL0, 0:nuoL0] = SigL0 # SGF^1 else: eSmHmS[0:nuoL0, 0:nuoL0] = eSmHmS[0:nuoL0, 0:nuoL0] - SigL0 tau = eSmHmS[0:devSt - 1, devSt - 1:devEndL].copy() taud = eSmHmS[devSt - 1:devEndL, 0:devSt - 1].copy() inv = LA.inv(eSmHmS[0:devSt - 1, 0:devSt - 1]) eSmHmS[devSt-1:devEndL, devSt-1:devEndL]=eSmHmS[devSt-1:devEndL, devSt-1:devEndL]-\ MM.mm(taud, inv, tau) self.SigL = eSmH[devSt - 1:devEndL, devSt - 1:devEndL] - eSmHmS[devSt - 1:devEndL, devSt - 1:devEndL] else: self.SigL = SigL0 self.GamL = 1.0j * (self.SigL - MM.dagger(self.SigL)) if self.Bulk and not FoldedL: # Reverse sign since SigL is really SGF^-1 self.GamL = -1.0 * self.GamL AssertReal(N.diag(self.GamL), 'GamL') if FoldedR: # Fold down from nuoR0 to the device region devStR = self.devStR eSmH = ee * self.S0 - self.H0 eSmHmS = eSmH[devStR - 1:nuo0, devStR - 1:nuo0].copy() tmpnuo = len(eSmHmS) if self.Bulk: eSmHmS[tmpnuo - nuoR0:tmpnuo, tmpnuo - nuoR0:tmpnuo] = SigR0 # SGF^1 else: eSmHmS[tmpnuo - nuoR0:tmpnuo, tmpnuo - nuoR0:tmpnuo] = eSmHmS[tmpnuo - nuoR0:tmpnuo, tmpnuo - nuoR0:tmpnuo] - SigR0 tau = eSmHmS[0:nuoR, nuoR:tmpnuo].copy() taud = eSmHmS[nuoR:tmpnuo, 0:nuoR].copy() inv = LA.inv(eSmHmS[nuoR:tmpnuo, nuoR:tmpnuo]) eSmHmS[0:nuoR, 0:nuoR] = eSmHmS[0:nuoR, 0:nuoR] - MM.mm(tau, inv, taud) self.SigR = eSmH[devStR - 1:devEnd, devStR - 1:devEnd] - eSmHmS[0:nuoR, 0:nuoR] else: self.SigR = SigR0 self.GamR = 1.0j * (self.SigR - MM.dagger(self.SigR)) if self.Bulk and not FoldedR: # Reverse sign since SigR is really SGF^-1 self.GamR = -1.0 * self.GamR AssertReal(N.diag(self.GamR), 'GamR')
def calcg0_old(self, ee, ispin=0, left=True): """ Only used if SciPy is not installed! For the left surface Green's function (1 is surface layer, 0 is all the other atoms): (E S00-H00 E S01-H01) (g00 g01) ( I 0 ) (E S10-H10 E S11-H11) * (g01 g11) = ( 0 I ) -> call E S - H for t ... t00 g01 + t01 g11 = 0 -> g01 = - t00^-1 t01 g11 t10 g01 + t11 g11 = I -> - t10 t00^-1 t01 g11 + t11 g11 = I -> And we get the surface Green's function: g11 = (t11 - t10 t00^-1 t01)^-1 with the right size of unitcell t00^-1 = g11! g11 = (E S11 - H11 - (E S10 - H10) g11 (E S01 - H01))^-1 In the calculations H01^+ and S01^+ are used instead of S10 and H10. (For complex energies (E S01 -H01)^+ is not (E S10 -H10) because the conjugate of the energy!!!!) For the right surface greens function same but different order on the MM.daggers! i.e., (E S - H - (E S01 - H01) gs (E S01^+ -H01^+) Algorith: Lopez Sancho*2 J Phys F:Met Phys 15 (1985) 851 I'm still very suspicios of this algorithm ... but it works and is really quick! The convergence is always checked against gs (E S - H - (E S01^+ - H01^+) gs (E S01 -H01) ) = I! """ H, S, H01, S01 = self.H[ispin, :, :], self.S, self.H01[ ispin, :, :], self.S01 alpha, beta = MM.dagger(H01) - ee * MM.dagger(S01), H01 - ee * S01 eps, epss = H.copy(), H.copy() converged = False iteration = 0 while not converged: iteration += 1 oldeps, oldepss = eps.copy(), epss.copy() oldalpha, oldbeta = alpha.copy(), beta.copy() tmpa = LA.solve(ee * S - oldeps, oldalpha) tmpb = LA.solve(ee * S - oldeps, oldbeta) alpha, beta = MM.mm(oldalpha, tmpa), MM.mm(oldbeta, tmpb) eps = oldeps + MM.mm(oldalpha, tmpb) + MM.mm(oldbeta, tmpa) if left: epss = oldepss + MM.mm(oldalpha, tmpb) else: epss = oldepss + MM.mm(oldbeta, tmpa) LopezConvTest = N.max(abs(alpha) + abs(beta)) if LopezConvTest < 1.0e-40: gs = LA.inv(ee * S - epss) if left: test = ee * S - H - MM.mm( ee * MM.dagger(S01) - MM.dagger(H01), gs, ee * S01 - H01) else: test = ee * S - H - MM.mm( ee * S01 - H01, gs, ee * MM.dagger(S01) - MM.dagger(H01)) myConvTest = N.max( abs( MM.mm(test, gs) - N.identity((self.HS.nuo), N.complex))) if myConvTest < VC.GetCheck("Lopez-Sancho"): converged = True if myConvTest > VC.GetCheck("Lopez-Sancho-warning"): v = "RIGHT" if left: v = "LEFT" print "WARNING: Lopez-scheme not-so-well converged for " + v + " electrode at E = %.4f eV:" % ee, myConvTest else: VC.Check("Lopez-Sancho", myConvTest, "Error: gs iteration {0}".format(iteration)) return gs
def calcg0(self, ee, ispin=0, left=True): # Calculate surface Green's function # Euro Phys J B 62, 381 (2008) # Inverse of : NOTE, setup for "right" lead. # e-h00 -h01 ... # -h10 e-h00 ... h00, s00, h01, s01 = self.H[ispin, :, :], self.S, self.H01[ ispin, :, :], self.S01 NN, ee = len(h00), N.real(ee) + N.max([N.imag(ee), 1e-8]) * 1.0j if left: h01, s01 = MM.dagger(h01), MM.dagger(s01) # Solve generalized eigen-problem # ( e I - h00 , -I) (eps) (h01 , 0) (eps) # ( h10 , 0) (xi ) = lambda (0 , I) (xi ) a, b = N.zeros((2 * NN, 2 * NN), N.complex), N.zeros((2 * NN, 2 * NN), N.complex) a[0:NN, 0:NN] = ee * s00 - h00 a[0:NN, NN:2 * NN] = -N.eye(NN) a[NN:2 * NN, 0:NN] = MM.dagger(h01) - ee * MM.dagger(s01) b[0:NN, 0:NN] = h01 - ee * s01 b[NN:2 * NN, NN:2 * NN] = N.eye(NN) ev, evec = SLA.eig(a, b) # Select lambda <0 and the eps part of the evec ipiv = N.where(N.abs(ev) < 1.0)[0] ev, evec = ev[ipiv], N.transpose(evec[:NN, ipiv]) # Normalize evec norm = N.sqrt(N.diag(MM.mm(evec, MM.dagger(evec)))) evec = MM.mm(N.diag(1.0 / norm), evec) # E^+ Lambda_+ (E^+)^-1 --->>> g00 EP = N.transpose(evec) FP = MM.mm(EP, N.diag(ev), LA.inv(MM.mm(MM.dagger(EP), EP)), MM.dagger(EP)) g00 = LA.inv(ee * s00 - h00 - MM.mm(h01 - ee * s01, FP)) # Check! err=N.max(N.abs(g00-LA.inv(ee*s00-h00-\ MM.mm(h01-ee*s01, g00, MM.dagger(h01)-ee*MM.dagger(s01))))) if err > 1.0e-8 and left: print "WARNING: Lopez-scheme not-so-well converged for LEFT electrode at E = %.4f eV:" % ee, err if err > 1.0e-8 and not left: print "WARNING: Lopez-scheme not-so-well converged for RIGHT electrode at E = %.4f eV:" % ee, err return g00