Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1import numpy as np
2from ase.io.jsonio import read_json, write_json
5class STM:
6 def __init__(self, atoms, symmetries=None, use_density=False):
7 """Scanning tunneling microscope.
9 atoms: Atoms object or filename
10 Atoms to scan or name of file to read LDOS from.
11 symmetries: list of int
12 List of integers 0, 1, and/or 2 indicating which surface
13 symmetries have been used to reduce the number of k-points
14 for the DFT calculation. The three integers correspond to
15 the following three symmetry operations::
17 [-1 0] [ 1 0] [ 0 1]
18 [ 0 1] [ 0 -1] [ 1 0]
20 use_density: bool
21 Use the electron density instead of the LDOS.
22 """
24 self.use_density = use_density
26 if isinstance(atoms, str):
27 with open(atoms, 'r') as fd:
28 self.ldos, self.bias, self.cell = read_json(fd,
29 always_array=False)
30 self.atoms = None
31 else:
32 self.atoms = atoms
33 self.cell = atoms.cell
34 self.bias = None
35 self.ldos = None
36 assert not self.cell[2, :2].any() and not self.cell[:2, 2].any()
38 self.symmetries = symmetries or []
40 def calculate_ldos(self, bias):
41 """Calculate local density of states for given bias."""
42 if self.ldos is not None and bias == self.bias:
43 return
45 self.bias = bias
47 calc = self.atoms.calc
49 if self.use_density:
50 self.ldos = calc.get_pseudo_density()
51 return
53 if bias < 0:
54 emin = bias
55 emax = 0.0
56 else:
57 emin = 0
58 emax = bias
60 nbands = calc.get_number_of_bands()
61 weights = calc.get_k_point_weights()
62 nkpts = len(weights)
63 nspins = calc.get_number_of_spins()
64 eigs = np.array([[calc.get_eigenvalues(k, s)
65 for k in range(nkpts)]
66 for s in range(nspins)])
67 eigs -= calc.get_fermi_level()
68 ldos = np.zeros(calc.get_pseudo_wave_function(0, 0, 0).shape)
70 for s in range(nspins):
71 for k in range(nkpts):
72 for n in range(nbands):
73 e = eigs[s, k, n]
74 if emin < e < emax:
75 psi = calc.get_pseudo_wave_function(n, k, s)
76 ldos += weights[k] * (psi * np.conj(psi)).real
78 if 0 in self.symmetries:
79 # (x,y) -> (-x,y)
80 ldos[1:] += ldos[:0:-1].copy()
81 ldos[1:] *= 0.5
83 if 1 in self.symmetries:
84 # (x,y) -> (x,-y)
85 ldos[:, 1:] += ldos[:, :0:-1].copy()
86 ldos[:, 1:] *= 0.5
88 if 2 in self.symmetries:
89 # (x,y) -> (y,x)
90 ldos += ldos.transpose((1, 0, 2)).copy()
91 ldos *= 0.5
93 self.ldos = ldos
95 def write(self, filename):
96 """Write local density of states to JSON file."""
97 write_json(filename, (self.ldos, self.bias, self.cell))
99 def get_averaged_current(self, bias, z):
100 """Calculate avarage current at height z (in Angstrom).
102 Use this to get an idea of what current to use when scanning."""
104 self.calculate_ldos(bias)
105 nz = self.ldos.shape[2]
107 # Find grid point:
108 n = z / self.cell[2, 2] * nz
109 dn = n - np.floor(n)
110 n = int(n) % nz
112 # Average and do linear interpolation:
113 return ((1 - dn) * self.ldos[:, :, n].mean() +
114 dn * self.ldos[:, :, (n + 1) % nz].mean())
116 def scan(self, bias, current, z0=None, repeat=(1, 1)):
117 """Constant current 2-d scan.
119 Returns three 2-d arrays (x, y, z) containing x-coordinates,
120 y-coordinates and heights. These three arrays can be passed to
121 matplotlibs contourf() function like this:
123 >>> import matplotlib.pyplot as plt
124 >>> plt.contourf(x, y, z)
125 >>> plt.show()
127 """
129 self.calculate_ldos(bias)
131 L = self.cell[2, 2]
132 nz = self.ldos.shape[2]
133 h = L / nz
135 ldos = self.ldos.reshape((-1, nz))
137 heights = np.empty(ldos.shape[0])
138 for i, a in enumerate(ldos):
139 heights[i] = find_height(a, current, h, z0)
141 s0 = heights.shape = self.ldos.shape[:2]
142 heights = np.tile(heights, repeat)
143 s = heights.shape
145 ij = np.indices(s, dtype=float).reshape((2, -1)).T
146 x, y = np.dot(ij / s0, self.cell[:2, :2]).T.reshape((2,) + s)
148 return x, y, heights
150 def scan2(self, bias, z, repeat=(1, 1)):
151 """Constant height 2-d scan.
153 Returns three 2-d arrays (x, y, I) containing x-coordinates,
154 y-coordinates and currents. These three arrays can be passed to
155 matplotlibs contourf() function like this:
157 >>> import matplotlib.pyplot as plt
158 >>> plt.contourf(x, y, I)
159 >>> plt.show()
161 """
163 self.calculate_ldos(bias)
165 nz = self.ldos.shape[2]
166 ldos = self.ldos.reshape((-1, nz))
168 I = np.empty(ldos.shape[0])
170 zp = z / self.cell[2, 2] * nz
171 zp = int(zp) % nz
173 for i, a in enumerate(ldos):
174 I[i] = self.find_current(a, zp)
176 s0 = I.shape = self.ldos.shape[:2]
177 I = np.tile(I, repeat)
178 s = I.shape
180 ij = np.indices(s, dtype=float).reshape((2, -1)).T
181 x, y = np.dot(ij / s0, self.cell[:2, :2]).T.reshape((2,) + s)
183 # Returing scan with axes in Angstrom.
184 return x, y, I
186 def linescan(self, bias, current, p1, p2, npoints=50, z0=None):
187 """Constant current line scan.
189 Example::
191 stm = STM(...)
192 z = ... # tip position
193 c = stm.get_averaged_current(-1.0, z)
194 stm.linescan(-1.0, c, (1.2, 0.0), (1.2, 3.0))
195 """
197 heights = self.scan(bias, current, z0)[2]
199 p1 = np.asarray(p1, float)
200 p2 = np.asarray(p2, float)
201 d = p2 - p1
202 s = np.dot(d, d)**0.5
204 cell = self.cell[:2, :2]
205 shape = np.array(heights.shape, float)
206 M = np.linalg.inv(cell)
207 line = np.empty(npoints)
208 for i in range(npoints):
209 p = p1 + i * d / (npoints - 1)
210 q = np.dot(p, M) * shape
211 line[i] = interpolate(q, heights)
212 return np.linspace(0, s, npoints), line
214 def pointcurrent(self, bias, x, y, z):
215 """Current for a single x, y, z position for a given bias."""
217 self.calculate_ldos(bias)
219 nx = self.ldos.shape[0]
220 ny = self.ldos.shape[1]
221 nz = self.ldos.shape[2]
223 # Find grid point:
224 xp = x / np.linalg.norm(self.cell[0]) * nx
225 dx = xp - np.floor(xp)
226 xp = int(xp) % nx
228 yp = y / np.linalg.norm(self.cell[1]) * ny
229 dy = yp - np.floor(yp)
230 yp = int(yp) % ny
232 zp = z / np.linalg.norm(self.cell[2]) * nz
233 dz = zp - np.floor(zp)
234 zp = int(zp) % nz
236 # 3D interpolation of the LDOS at point (x,y,z) at given bias.
237 xyzldos = (((1 - dx) + (1 - dy) + (1 - dz)) * self.ldos[xp, yp, zp] +
238 dx * self.ldos[(xp + 1) % nx, yp, zp] +
239 dy * self.ldos[xp, (yp + 1) % ny, zp] +
240 dz * self.ldos[xp, yp, (zp + 1) % nz])
242 return dos2current(bias, xyzldos)
244 def sts(self, x, y, z, bias0, bias1, biasstep):
245 """Returns the dI/dV curve for position x, y at height z (in Angstrom),
246 for bias from bias0 to bias1 with step biasstep."""
248 biases = np.arange(bias0, bias1+biasstep, biasstep)
249 I = np.zeros(biases.shape)
251 for b in np.arange(len(biases)):
252 print(b, biases[b])
253 I[b] = self.pointcurrent(biases[b], x, y, z)
255 dIdV = np.gradient(I, biasstep)
257 return biases, I, dIdV
259 def find_current(self, ldos, z):
260 """ Finds current for given LDOS at height z."""
261 nz = self.ldos.shape[2]
263 zp = z / self.cell[2, 2] * nz
264 dz = zp - np.floor(zp)
265 zp = int(zp) % nz
267 ldosz = (1 - dz) * ldos[zp] + dz * ldos[(zp + 1) % nz]
269 return dos2current(self.bias, ldosz)
272def dos2current(bias, dos):
273 # Borrowed from gpaw/analyse/simple_stm.py:
274 # The connection between density n and current I
275 # n [e/Angstrom^3] = 0.0002 sqrt(I [nA])
276 # as given in Hofer et al., RevModPhys 75 (2003) 1287
277 return 5000. * dos**2 * (1 if bias > 0 else -1)
280def interpolate(q, heights):
281 qi = q.astype(int)
282 f = q - qi
283 g = 1 - f
284 qi %= heights.shape
285 n0, m0 = qi
286 n1, m1 = (qi + 1) % heights.shape
287 z = (g[0] * g[1] * heights[n0, m0] +
288 f[0] * g[1] * heights[n1, m0] +
289 g[0] * f[1] * heights[n0, m1] +
290 f[0] * f[1] * heights[n1, m1])
291 return z
294def find_height(ldos, current, h, z0=None):
295 if z0 is None:
296 n = len(ldos) - 2
297 else:
298 n = int(z0 / h)
299 while n >= 0:
300 if ldos[n] > current:
301 break
302 n -= 1
303 else:
304 return 0.0
306 c2, c1 = ldos[n:n + 2]
307 return (n + 1 - (current - c1) / (c2 - c1)) * h
310def delta(biases, bias, width):
311 """Return a delta-function centered at 'bias'"""
312 x = -((biases - bias) / width)**2
313 return np.exp(x) / (np.sqrt(np.pi) * width)