root / ase / transport / calculators.py @ 7
Historique | Voir | Annoter | Télécharger (15,15 ko)
| 1 |
import numpy as np |
|---|---|
| 2 |
|
| 3 |
from numpy import linalg |
| 4 |
from ase.transport.selfenergy import LeadSelfEnergy, BoxProbe |
| 5 |
from ase.transport.greenfunction import GreenFunction |
| 6 |
from ase.transport.tools import subdiagonalize, cutcoupling, tri2full, dagger,\ |
| 7 |
rotate_matrix
|
| 8 |
|
| 9 |
|
| 10 |
class TransportCalculator: |
| 11 |
"""Determine transport properties of a device sandwiched between
|
| 12 |
two semi-infinite leads using a Green function method.
|
| 13 |
"""
|
| 14 |
|
| 15 |
def __init__(self, **kwargs): |
| 16 |
"""Create the transport calculator.
|
| 17 |
|
| 18 |
Parameters
|
| 19 |
==========
|
| 20 |
h : (N, N) ndarray
|
| 21 |
Hamiltonian matrix for the central region.
|
| 22 |
s : {None, (N, N) ndarray}, optional
|
| 23 |
Overlap matrix for the central region.
|
| 24 |
Use None for an orthonormal basis.
|
| 25 |
h1 : (N1, N1) ndarray
|
| 26 |
Hamiltonian matrix for lead1.
|
| 27 |
h2 : {None, (N2, N2) ndarray}, optional
|
| 28 |
Hamiltonian matrix for lead2. You may use None if lead1 and lead2
|
| 29 |
are identical.
|
| 30 |
s1 : {None, (N1, N1) ndarray}, optional
|
| 31 |
Overlap matrix for lead1. Use None for an orthonomormal basis.
|
| 32 |
hc1 : {None, (N1, N) ndarray}, optional
|
| 33 |
Hamiltonian coupling matrix between the first principal
|
| 34 |
layer in lead1 and the central region.
|
| 35 |
hc2 : {None, (N2, N} ndarray), optional
|
| 36 |
Hamiltonian coupling matrix between the first principal
|
| 37 |
layer in lead2 and the central region.
|
| 38 |
sc1 : {None, (N1, N) ndarray}, optional
|
| 39 |
Overlap coupling matrix between the first principal
|
| 40 |
layer in lead1 and the central region.
|
| 41 |
sc2 : {None, (N2, N) ndarray}, optional
|
| 42 |
Overlap coupling matrix between the first principal
|
| 43 |
layer in lead2 and the central region.
|
| 44 |
energies : {None, array_like}, optional
|
| 45 |
Energy points for which calculated transport properties are
|
| 46 |
evaluated.
|
| 47 |
eta : {1.0e-5, float}, optional
|
| 48 |
Infinitesimal for the central region Green function.
|
| 49 |
eta1/eta2 : {1.0e-5, float}, optional
|
| 50 |
Infinitesimal for lead1/lead2 Green function.
|
| 51 |
align_bf : {None, int}, optional
|
| 52 |
Use align_bf=m to shift the central region
|
| 53 |
by a constant potential such that the m'th onsite element
|
| 54 |
in the central region is aligned to the m'th onsite element
|
| 55 |
in lead1 principal layer.
|
| 56 |
logfile : {None, str}, optional
|
| 57 |
Write a logfile to file with name `logfile`.
|
| 58 |
Use '-' to write to std out.
|
| 59 |
eigenchannels: {0, int}, optional
|
| 60 |
Number of eigenchannel transmission coefficients to
|
| 61 |
calculate.
|
| 62 |
pdos : {None, (N,) array_like}, optional
|
| 63 |
Specify which basis functions to calculate the
|
| 64 |
projected density of states for.
|
| 65 |
dos : {False, bool}, optional
|
| 66 |
The total density of states of the central region.
|
| 67 |
box: XXX
|
| 68 |
YYY
|
| 69 |
|
| 70 |
If hc1/hc2 are None, they are assumed to be identical to
|
| 71 |
the coupling matrix elements between neareste neighbor
|
| 72 |
principal layers in lead1/lead2.
|
| 73 |
|
| 74 |
Examples
|
| 75 |
========
|
| 76 |
>>> import numpy as np
|
| 77 |
>>> h = np.array((0,)).reshape((1,1))
|
| 78 |
>>> h1 = np.array((0, -1, -1, 0)).reshape(2,2)
|
| 79 |
>>> energies = np.arange(-3, 3, 0.1)
|
| 80 |
>>> calc = TransportCalculator(h=h, h1=h1, energies=energies)
|
| 81 |
>>> T = calc.get_transmission()
|
| 82 |
|
| 83 |
"""
|
| 84 |
|
| 85 |
# The default values for all extra keywords
|
| 86 |
self.input_parameters = {'energies': None, |
| 87 |
'h': None, |
| 88 |
'h1': None, |
| 89 |
'h2': None, |
| 90 |
's': None, |
| 91 |
's1': None, |
| 92 |
's2': None, |
| 93 |
'hc1': None, |
| 94 |
'hc2': None, |
| 95 |
'sc1': None, |
| 96 |
'sc2': None, |
| 97 |
'box': None, |
| 98 |
'align_bf': None, |
| 99 |
'eta1': 1e-5, |
| 100 |
'eta2': 1e-5, |
| 101 |
'eta': 1e-5, |
| 102 |
'logfile': None, # '-', |
| 103 |
'eigenchannels': 0, |
| 104 |
'dos': False, |
| 105 |
'pdos': [],
|
| 106 |
} |
| 107 |
self.initialized = False # Changed Hamiltonians? |
| 108 |
self.uptodate = False # Changed energy grid? |
| 109 |
self.set(**kwargs)
|
| 110 |
|
| 111 |
def set(self, **kwargs): |
| 112 |
for key in kwargs: |
| 113 |
if key in ['h', 'h1', 'h2', 'hc1', 'hc2', |
| 114 |
's', 's1', 's2', 'sc1', 'sc2', |
| 115 |
'eta', 'eta1', 'eta2', 'align_bf', 'box']: |
| 116 |
self.initialized = False |
| 117 |
self.uptodate = False |
| 118 |
break
|
| 119 |
elif key in ['energies', 'eigenchannels', 'dos', 'pdos']: |
| 120 |
self.uptodate = False |
| 121 |
elif key not in self.input_parameters: |
| 122 |
raise KeyError, '\'%s\' not a vaild keyword' % key |
| 123 |
|
| 124 |
self.input_parameters.update(kwargs)
|
| 125 |
log = self.input_parameters['logfile'] |
| 126 |
if log is None: |
| 127 |
class Trash: |
| 128 |
def write(self, s): |
| 129 |
pass
|
| 130 |
def flush(self): |
| 131 |
pass
|
| 132 |
self.log = Trash()
|
| 133 |
elif log == '-': |
| 134 |
from sys import stdout |
| 135 |
self.log = stdout
|
| 136 |
elif 'logfile' in kwargs: |
| 137 |
self.log = open(log, 'w') |
| 138 |
|
| 139 |
def initialize(self): |
| 140 |
if self.initialized: |
| 141 |
return
|
| 142 |
|
| 143 |
print >> self.log, '# Initializing calculator...' |
| 144 |
|
| 145 |
p = self.input_parameters
|
| 146 |
if p['s'] == None: |
| 147 |
p['s'] = np.identity(len(p['h'])) |
| 148 |
|
| 149 |
identical_leads = False
|
| 150 |
if p['h2'] == None: |
| 151 |
p['h2'] = p['h1'] # Lead2 is idendical to lead1 |
| 152 |
identical_leads = True
|
| 153 |
|
| 154 |
if p['s1'] == None: |
| 155 |
p['s1'] = np.identity(len(p['h1'])) |
| 156 |
|
| 157 |
if p['s2'] == None and not identical_leads: |
| 158 |
p['s2'] = np.identity(len(p['h2'])) # Orthonormal basis for lead 2 |
| 159 |
else: # Lead2 is idendical to lead1 |
| 160 |
p['s2'] = p['s1'] |
| 161 |
|
| 162 |
|
| 163 |
h_mm = p['h']
|
| 164 |
s_mm = p['s']
|
| 165 |
pl1 = len(p['h1']) / 2 |
| 166 |
pl2 = len(p['h2']) / 2 |
| 167 |
h1_ii = p['h1'][:pl1, :pl1]
|
| 168 |
h1_ij = p['h1'][:pl1, pl1:2 * pl1] |
| 169 |
s1_ii = p['s1'][:pl1, :pl1]
|
| 170 |
s1_ij = p['s1'][:pl1, pl1:2 * pl1] |
| 171 |
h2_ii = p['h2'][:pl2, :pl2]
|
| 172 |
h2_ij = p['h2'][pl2: 2 * pl2, :pl2] |
| 173 |
s2_ii = p['s2'][:pl2, :pl2]
|
| 174 |
s2_ij = p['s2'][pl2: 2 * pl2, :pl2] |
| 175 |
|
| 176 |
if p['hc1'] is None: |
| 177 |
nbf = len(h_mm)
|
| 178 |
h1_im = np.zeros((pl1, nbf), complex)
|
| 179 |
s1_im = np.zeros((pl1, nbf), complex)
|
| 180 |
h1_im[:pl1, :pl1] = h1_ij |
| 181 |
s1_im[:pl1, :pl1] = s1_ij |
| 182 |
p['hc1'] = h1_im
|
| 183 |
p['sc1'] = s1_im
|
| 184 |
else:
|
| 185 |
h1_im = p['hc1']
|
| 186 |
if p['sc1'] is not None: |
| 187 |
s1_im = p['sc1']
|
| 188 |
else:
|
| 189 |
s1_im = np.zeros(h1_im.shape, complex)
|
| 190 |
p['sc1'] = s1_im
|
| 191 |
|
| 192 |
if p['hc2'] is None: |
| 193 |
h2_im = np.zeros((pl2, nbf), complex)
|
| 194 |
s2_im = np.zeros((pl2, nbf), complex)
|
| 195 |
h2_im[-pl2:, -pl2:] = h2_ij |
| 196 |
s2_im[-pl2:, -pl2:] = s2_ij |
| 197 |
p['hc2'] = h2_im
|
| 198 |
p['sc2'] = s2_im
|
| 199 |
else:
|
| 200 |
h2_im = p['hc2']
|
| 201 |
if p['sc2'] is not None: |
| 202 |
s2_im = p['sc2']
|
| 203 |
else:
|
| 204 |
s2_im = np.zeros(h2_im.shape, complex)
|
| 205 |
p['sc2'] = s2_im
|
| 206 |
|
| 207 |
align_bf = p['align_bf']
|
| 208 |
if align_bf != None: |
| 209 |
diff = (h_mm[align_bf, align_bf] - h1_ii[align_bf, align_bf]) \ |
| 210 |
/ s_mm[align_bf, align_bf] |
| 211 |
print >> self.log, '# Aligning scat. H to left lead H. diff=', diff |
| 212 |
h_mm -= diff * s_mm |
| 213 |
|
| 214 |
# Setup lead self-energies
|
| 215 |
# All infinitesimals must be > 0
|
| 216 |
assert np.all(np.array((p['eta'], p['eta1'], p['eta2'])) > 0.0) |
| 217 |
self.selfenergies = [LeadSelfEnergy((h1_ii, s1_ii),
|
| 218 |
(h1_ij, s1_ij), |
| 219 |
(h1_im, s1_im), |
| 220 |
p['eta1']),
|
| 221 |
LeadSelfEnergy((h2_ii, s2_ii), |
| 222 |
(h2_ij, s2_ij), |
| 223 |
(h2_im, s2_im), |
| 224 |
p['eta2'])]
|
| 225 |
box = p['box']
|
| 226 |
if box is not None: |
| 227 |
print 'Using box probe!' |
| 228 |
self.selfenergies.append(
|
| 229 |
BoxProbe(eta=box[0], a=box[1], b=box[2], energies=box[3], |
| 230 |
S=s_mm, T=0.3))
|
| 231 |
|
| 232 |
#setup scattering green function
|
| 233 |
self.greenfunction = GreenFunction(selfenergies=self.selfenergies, |
| 234 |
H=h_mm, |
| 235 |
S=s_mm, |
| 236 |
eta=p['eta'])
|
| 237 |
|
| 238 |
self.initialized = True |
| 239 |
|
| 240 |
def update(self): |
| 241 |
if self.uptodate: |
| 242 |
return
|
| 243 |
|
| 244 |
p = self.input_parameters
|
| 245 |
self.energies = p['energies'] |
| 246 |
nepts = len(self.energies) |
| 247 |
nchan = p['eigenchannels']
|
| 248 |
pdos = p['pdos']
|
| 249 |
self.T_e = np.empty(nepts)
|
| 250 |
if p['dos']: |
| 251 |
self.dos_e = np.empty(nepts)
|
| 252 |
if pdos != []:
|
| 253 |
self.pdos_ne = np.empty((len(pdos), nepts)) |
| 254 |
if nchan > 0: |
| 255 |
self.eigenchannels_ne = np.empty((nchan, nepts))
|
| 256 |
|
| 257 |
for e, energy in enumerate(self.energies): |
| 258 |
Ginv_mm = self.greenfunction.retarded(energy, inverse=True) |
| 259 |
lambda1_mm = self.selfenergies[0].get_lambda(energy) |
| 260 |
lambda2_mm = self.selfenergies[1].get_lambda(energy) |
| 261 |
a_mm = linalg.solve(Ginv_mm, lambda1_mm) |
| 262 |
b_mm = linalg.solve(dagger(Ginv_mm), lambda2_mm) |
| 263 |
T_mm = np.dot(a_mm, b_mm) |
| 264 |
if nchan > 0: |
| 265 |
t_n = linalg.eigvals(T_mm).real |
| 266 |
self.eigenchannels_ne[:, e] = np.sort(t_n)[-nchan:]
|
| 267 |
self.T_e[e] = np.sum(t_n)
|
| 268 |
else:
|
| 269 |
self.T_e[e] = np.trace(T_mm).real
|
| 270 |
|
| 271 |
print >> self.log, energy, self.T_e[e] |
| 272 |
self.log.flush()
|
| 273 |
|
| 274 |
if p['dos']: |
| 275 |
self.dos_e[e] = self.greenfunction.dos(energy) |
| 276 |
|
| 277 |
if pdos != []:
|
| 278 |
self.pdos_ne[:, e] = np.take(self.greenfunction.pdos(energy), |
| 279 |
pdos) |
| 280 |
|
| 281 |
self.uptodate = True |
| 282 |
|
| 283 |
def print_pl_convergence(self): |
| 284 |
self.initialize()
|
| 285 |
pl1 = len(self.input_parameters['h1']) / 2 |
| 286 |
|
| 287 |
h_ii = self.selfenergies[0].h_ii |
| 288 |
s_ii = self.selfenergies[0].s_ii |
| 289 |
ha_ii = self.greenfunction.H[:pl1, :pl1]
|
| 290 |
sa_ii = self.greenfunction.S[:pl1, :pl1]
|
| 291 |
c1 = np.abs(h_ii - ha_ii).max() |
| 292 |
c2 = np.abs(s_ii - sa_ii).max() |
| 293 |
print 'Conv (h,s)=%.2e, %2.e' % (c1, c2) |
| 294 |
|
| 295 |
def plot_pl_convergence(self): |
| 296 |
self.initialize()
|
| 297 |
pl1 = len(self.input_parameters['h1']) / 2 |
| 298 |
hlead = self.selfenergies[0].h_ii.real.diagonal() |
| 299 |
hprincipal = self.greenfunction.H.real.diagonal[:pl1]
|
| 300 |
|
| 301 |
import pylab as pl |
| 302 |
pl.plot(hlead, label='lead')
|
| 303 |
pl.plot(hprincipal, label='principal layer')
|
| 304 |
pl.axis('tight')
|
| 305 |
pl.show() |
| 306 |
|
| 307 |
def get_transmission(self): |
| 308 |
self.initialize()
|
| 309 |
self.update()
|
| 310 |
return self.T_e |
| 311 |
|
| 312 |
def get_dos(self): |
| 313 |
self.initialize()
|
| 314 |
self.update()
|
| 315 |
return self.dos_e |
| 316 |
|
| 317 |
def get_eigenchannels(self, n=None): |
| 318 |
"""Get ``n`` first eigenchannels."""
|
| 319 |
self.initialize()
|
| 320 |
self.update()
|
| 321 |
if n is None: |
| 322 |
n = self.input_parameters['eigenchannels'] |
| 323 |
return self.eigenchannels_ne[:n] |
| 324 |
|
| 325 |
def get_pdos(self): |
| 326 |
self.initialize()
|
| 327 |
self.update()
|
| 328 |
return self.pdos_ne |
| 329 |
|
| 330 |
def subdiagonalize_bfs(self, bfs, apply=False): |
| 331 |
self.initialize()
|
| 332 |
bfs = np.array(bfs) |
| 333 |
p = self.input_parameters
|
| 334 |
h_mm = p['h']
|
| 335 |
s_mm = p['s']
|
| 336 |
ht_mm, st_mm, c_mm, e_m = subdiagonalize(h_mm, s_mm, bfs) |
| 337 |
if apply: |
| 338 |
self.uptodate = False |
| 339 |
h_mm[:] = ht_mm |
| 340 |
s_mm[:] = st_mm |
| 341 |
# Rotate coupling between lead and central region
|
| 342 |
for alpha, sigma in enumerate(self.selfenergies): |
| 343 |
sigma.h_im[:] = np.dot(sigma.h_im, c_mm) |
| 344 |
sigma.s_im[:] = np.dot(sigma.s_im, c_mm) |
| 345 |
|
| 346 |
c_mm = np.take(c_mm, bfs, axis=0)
|
| 347 |
c_mm = np.take(c_mm, bfs, axis=1)
|
| 348 |
return ht_mm, st_mm, e_m, c_mm
|
| 349 |
|
| 350 |
def cutcoupling_bfs(self, bfs, apply=False): |
| 351 |
self.initialize()
|
| 352 |
bfs = np.array(bfs) |
| 353 |
p = self.input_parameters
|
| 354 |
h_pp = p['h'].copy()
|
| 355 |
s_pp = p['s'].copy()
|
| 356 |
cutcoupling(h_pp, s_pp, bfs) |
| 357 |
if apply: |
| 358 |
self.uptodate = False |
| 359 |
p['h'][:] = h_pp
|
| 360 |
p['s'][:] = s_pp
|
| 361 |
for alpha, sigma in enumerate(self.selfenergies): |
| 362 |
for m in bfs: |
| 363 |
sigma.h_im[:, m] = 0.0
|
| 364 |
sigma.s_im[:, m] = 0.0
|
| 365 |
return h_pp, s_pp
|
| 366 |
|
| 367 |
def lowdin_rotation(self, apply=False): |
| 368 |
p = self.input_parameters
|
| 369 |
h_mm = p['h']
|
| 370 |
s_mm = p['s']
|
| 371 |
eig, rot_mm = linalg.eigh(s_mm) |
| 372 |
eig = np.abs(eig) |
| 373 |
rot_mm = np.dot(rot_mm / np.sqrt(eig), dagger(rot_mm)) |
| 374 |
if apply: |
| 375 |
self.uptodate = False |
| 376 |
h_mm[:] = rotate_matrix(h_mm, rot_mm) # rotate C region
|
| 377 |
s_mm[:] = rotate_matrix(s_mm, rot_mm) |
| 378 |
for alpha, sigma in enumerate(self.selfenergies): |
| 379 |
sigma.h_im[:] = np.dot(sigma.h_im, rot_mm) # rotate L-C coupl.
|
| 380 |
sigma.s_im[:] = np.dot(sigma.s_im, rot_mm) |
| 381 |
|
| 382 |
return rot_mm
|
| 383 |
|
| 384 |
def get_left_channels(self, energy, nchan=1): |
| 385 |
self.initialize()
|
| 386 |
g_s_ii = self.greenfunction.retarded(energy)
|
| 387 |
lambda_l_ii = self.selfenergies[0].get_lambda(energy) |
| 388 |
lambda_r_ii = self.selfenergies[1].get_lambda(energy) |
| 389 |
|
| 390 |
if self.greenfunction.S is None: |
| 391 |
s_s_qsrt_ii = s_s_isqrt = np.identity(len(g_s_ii))
|
| 392 |
else:
|
| 393 |
s_mm = self.greenfunction.S
|
| 394 |
s_s_i, s_s_ii = linalg.eig(s_mm) |
| 395 |
s_s_i = np.abs(s_s_i) |
| 396 |
s_s_sqrt_i = np.sqrt(s_s_i) # sqrt of eigenvalues
|
| 397 |
s_s_sqrt_ii = np.dot(s_s_ii * s_s_sqrt_i, dagger(s_s_ii)) |
| 398 |
s_s_isqrt_ii = np.dot(s_s_ii / s_s_sqrt_i, dagger(s_s_ii)) |
| 399 |
|
| 400 |
lambdab_r_ii = np.dot(np.dot(s_s_isqrt_ii, lambda_r_ii),s_s_isqrt_ii) |
| 401 |
a_l_ii = np.dot(np.dot(g_s_ii, lambda_l_ii), dagger(g_s_ii)) |
| 402 |
ab_l_ii = np.dot(np.dot(s_s_sqrt_ii, a_l_ii), s_s_sqrt_ii) |
| 403 |
lambda_i, u_ii = linalg.eig(ab_l_ii) |
| 404 |
ut_ii = np.sqrt(lambda_i / (2.0 * np.pi)) * u_ii
|
| 405 |
m_ii = 2 * np.pi * np.dot(np.dot(dagger(ut_ii), lambdab_r_ii),ut_ii)
|
| 406 |
T_i,c_in = linalg.eig(m_ii) |
| 407 |
T_i = np.abs(T_i) |
| 408 |
|
| 409 |
channels = np.argsort(-T_i)[:nchan] |
| 410 |
c_in = np.take(c_in, channels, axis=1)
|
| 411 |
T_n = np.take(T_i, channels) |
| 412 |
v_in = np.dot(np.dot(s_s_isqrt_ii, ut_ii), c_in) |
| 413 |
|
| 414 |
return T_n, v_in
|