root / ase / dft / wannier.py @ 20
Historique | Voir | Annoter | Télécharger (28,24 ko)
1 | 1 | tkerber | """ Maximally localized Wannier Functions
|
---|---|---|---|
2 | 1 | tkerber |
|
3 | 1 | tkerber | Find the set of maximally localized Wannier functions
|
4 | 1 | tkerber | using the spread functional of Marzari and Vanderbilt
|
5 | 1 | tkerber | (PRB 56, 1997 page 12847).
|
6 | 1 | tkerber | """
|
7 | 1 | tkerber | import numpy as np |
8 | 1 | tkerber | from time import time |
9 | 1 | tkerber | from math import sqrt, pi |
10 | 1 | tkerber | from pickle import dump, load |
11 | 1 | tkerber | from ase.parallel import paropen |
12 | 1 | tkerber | from ase.calculators.dacapo import Dacapo |
13 | 1 | tkerber | from ase.dft.kpoints import get_monkhorst_shape |
14 | 1 | tkerber | from ase.transport.tools import dagger, normalize |
15 | 1 | tkerber | |
16 | 1 | tkerber | dag = dagger |
17 | 1 | tkerber | |
18 | 1 | tkerber | |
19 | 1 | tkerber | def gram_schmidt(U): |
20 | 1 | tkerber | """Orthonormalize columns of U according to the Gram-Schmidt procedure."""
|
21 | 1 | tkerber | for i, col in enumerate(U.T): |
22 | 1 | tkerber | for col2 in U.T[:i]: |
23 | 1 | tkerber | col -= col2 * np.dot(col2.conj(), col) |
24 | 1 | tkerber | col /= np.linalg.norm(col) |
25 | 1 | tkerber | |
26 | 1 | tkerber | |
27 | 1 | tkerber | def lowdin(U, S=None): |
28 | 1 | tkerber | """Orthonormalize columns of U according to the Lowdin procedure.
|
29 | 1 | tkerber |
|
30 | 1 | tkerber | If the overlap matrix is know, it can be specified in S.
|
31 | 1 | tkerber | """
|
32 | 1 | tkerber | if S is None: |
33 | 1 | tkerber | S = np.dot(dag(U), U) |
34 | 1 | tkerber | eig, rot = np.linalg.eigh(S) |
35 | 1 | tkerber | rot = np.dot(rot / np.sqrt(eig), dag(rot)) |
36 | 1 | tkerber | U[:] = np.dot(U, rot) |
37 | 1 | tkerber | |
38 | 1 | tkerber | |
39 | 1 | tkerber | def neighbor_k_search(k_c, G_c, kpt_kc, tol=1e-4): |
40 | 1 | tkerber | # search for k1 (in kpt_kc) and k0 (in alldir), such that
|
41 | 1 | tkerber | # k1 - k - G + k0 = 0
|
42 | 1 | tkerber | alldir_dc = np.array([[0,0,0],[1,0,0],[0,1,0],[0,0,1], |
43 | 1 | tkerber | [1,1,0],[1,0,1],[0,1,1]], int) |
44 | 1 | tkerber | for k0_c in alldir_dc: |
45 | 1 | tkerber | for k1, k1_c in enumerate(kpt_kc): |
46 | 1 | tkerber | if np.linalg.norm(k1_c - k_c - G_c + k0_c) < tol:
|
47 | 1 | tkerber | return k1, k0_c
|
48 | 1 | tkerber | |
49 | 1 | tkerber | print 'Wannier: Did not find matching kpoint for kpt=', k_c |
50 | 1 | tkerber | print 'Probably non-uniform k-point grid' |
51 | 1 | tkerber | raise NotImplementedError |
52 | 1 | tkerber | |
53 | 1 | tkerber | |
54 | 1 | tkerber | def calculate_weights(cell_cc): |
55 | 1 | tkerber | """ Weights are used for non-cubic cells, see PRB **61**, 10040"""
|
56 | 1 | tkerber | alldirs_dc = np.array([[1, 0, 0], [0, 1, 0], [0, 0, 1], |
57 | 1 | tkerber | [1, 1, 0], [1, 0, 1], [0, 1, 1]], dtype=int) |
58 | 1 | tkerber | g = np.dot(cell_cc, cell_cc.T) |
59 | 1 | tkerber | # NOTE: Only first 3 of following 6 weights are presently used:
|
60 | 1 | tkerber | w = np.zeros(6)
|
61 | 1 | tkerber | w[0] = g[0, 0] - g[0, 1] - g[0, 2] |
62 | 1 | tkerber | w[1] = g[1, 1] - g[0, 1] - g[1, 2] |
63 | 1 | tkerber | w[2] = g[2, 2] - g[0, 2] - g[1, 2] |
64 | 1 | tkerber | w[3] = g[0, 1] |
65 | 1 | tkerber | w[4] = g[0, 2] |
66 | 1 | tkerber | w[5] = g[1, 2] |
67 | 1 | tkerber | # Make sure that first 3 Gdir vectors are included -
|
68 | 1 | tkerber | # these are used to calculate Wanniercenters.
|
69 | 1 | tkerber | Gdir_dc = alldirs_dc[:3]
|
70 | 1 | tkerber | weight_d = w[:3]
|
71 | 1 | tkerber | for d in range(3, 6): |
72 | 1 | tkerber | if abs(w[d]) > 1e-5: |
73 | 1 | tkerber | Gdir_dc = np.concatenate(Gdir_dc, alldirs_dc[d]) |
74 | 1 | tkerber | weight_d = np.concatenate(weight_d, w[d]) |
75 | 1 | tkerber | weight_d /= max(abs(weight_d)) |
76 | 1 | tkerber | return weight_d, Gdir_dc
|
77 | 1 | tkerber | |
78 | 1 | tkerber | |
79 | 1 | tkerber | def random_orthogonal_matrix(dim, seed=None, real=False): |
80 | 1 | tkerber | """Generate a random orthogonal matrix"""
|
81 | 1 | tkerber | if seed is not None: |
82 | 1 | tkerber | np.random.seed(seed) |
83 | 1 | tkerber | |
84 | 1 | tkerber | H = np.random.rand(dim, dim) |
85 | 1 | tkerber | np.add(dag(H), H, H) |
86 | 1 | tkerber | np.multiply(.5, H, H)
|
87 | 1 | tkerber | |
88 | 1 | tkerber | if real:
|
89 | 1 | tkerber | gram_schmidt(H) |
90 | 1 | tkerber | return H
|
91 | 1 | tkerber | else:
|
92 | 1 | tkerber | val, vec = np.linalg.eig(H) |
93 | 1 | tkerber | return np.dot(vec * np.exp(1.j * val), dag(vec)) |
94 | 1 | tkerber | |
95 | 1 | tkerber | |
96 | 1 | tkerber | def steepest_descent(func, step=.005, tolerance=1e-6, **kwargs): |
97 | 1 | tkerber | fvalueold = 0.
|
98 | 1 | tkerber | fvalue = fvalueold + 10
|
99 | 1 | tkerber | count=0
|
100 | 1 | tkerber | while abs((fvalue - fvalueold) / fvalue) > tolerance: |
101 | 1 | tkerber | fvalueold = fvalue |
102 | 1 | tkerber | dF = func.get_gradients() |
103 | 1 | tkerber | func.step(dF * step, **kwargs) |
104 | 1 | tkerber | fvalue = func.get_functional_value() |
105 | 1 | tkerber | count += 1
|
106 | 1 | tkerber | print 'SteepestDescent: iter=%s, value=%s' % (count, fvalue) |
107 | 1 | tkerber | |
108 | 1 | tkerber | |
109 | 1 | tkerber | def md_min(func, step=.25, tolerance=1e-6, verbose=False, **kwargs): |
110 | 1 | tkerber | if verbose:
|
111 | 1 | tkerber | print 'Localize with step =', step, 'and tolerance =', tolerance |
112 | 1 | tkerber | t = -time() |
113 | 1 | tkerber | fvalueold = 0.
|
114 | 1 | tkerber | fvalue = fvalueold + 10
|
115 | 1 | tkerber | count = 0
|
116 | 1 | tkerber | V = np.zeros(func.get_gradients().shape, dtype=complex)
|
117 | 1 | tkerber | while abs((fvalue - fvalueold) / fvalue) > tolerance: |
118 | 1 | tkerber | fvalueold = fvalue |
119 | 1 | tkerber | dF = func.get_gradients() |
120 | 1 | tkerber | V *= (dF * V.conj()).real > 0
|
121 | 1 | tkerber | V += step * dF |
122 | 1 | tkerber | func.step(V, **kwargs) |
123 | 1 | tkerber | fvalue = func.get_functional_value() |
124 | 1 | tkerber | if fvalue < fvalueold:
|
125 | 1 | tkerber | step *= 0.5
|
126 | 1 | tkerber | count += 1
|
127 | 1 | tkerber | if verbose:
|
128 | 1 | tkerber | print 'MDmin: iter=%s, step=%s, value=%s' % (count, step, fvalue) |
129 | 1 | tkerber | if verbose:
|
130 | 1 | tkerber | t += time() |
131 | 1 | tkerber | print '%d iterations in %0.2f seconds (%0.2f ms/iter), endstep = %s' %( |
132 | 1 | tkerber | count, t, t * 1000. / count, step)
|
133 | 1 | tkerber | |
134 | 1 | tkerber | |
135 | 1 | tkerber | def rotation_from_projection(proj_nw, fixed, ortho=True): |
136 | 1 | tkerber | """Determine rotation and coefficient matrices from projections
|
137 | 1 | tkerber |
|
138 | 1 | tkerber | proj_nw = <psi_n|p_w>
|
139 | 1 | tkerber | psi_n: eigenstates
|
140 | 1 | tkerber | p_w: localized function
|
141 | 1 | tkerber |
|
142 | 1 | tkerber | Nb (n) = Number of bands
|
143 | 1 | tkerber | Nw (w) = Number of wannier functions
|
144 | 1 | tkerber | M (f) = Number of fixed states
|
145 | 1 | tkerber | L (l) = Number of extra degrees of freedom
|
146 | 1 | tkerber | U (u) = Number of non-fixed states
|
147 | 1 | tkerber | """
|
148 | 1 | tkerber | |
149 | 1 | tkerber | Nb, Nw = proj_nw.shape |
150 | 1 | tkerber | M = fixed |
151 | 1 | tkerber | L = Nw - M |
152 | 1 | tkerber | |
153 | 1 | tkerber | U_ww = np.empty((Nw, Nw), dtype=proj_nw.dtype) |
154 | 1 | tkerber | U_ww[:M] = proj_nw[:M] |
155 | 1 | tkerber | |
156 | 1 | tkerber | if L > 0: |
157 | 1 | tkerber | proj_uw = proj_nw[M:] |
158 | 1 | tkerber | eig_w, C_ww = np.linalg.eigh(np.dot(dag(proj_uw), proj_uw)) |
159 | 1 | tkerber | C_ul = np.dot(proj_uw, C_ww[:, np.argsort(-eig_w.real)[:L]]) |
160 | 1 | tkerber | #eig_u, C_uu = np.linalg.eigh(np.dot(proj_uw, dag(proj_uw)))
|
161 | 1 | tkerber | #C_ul = C_uu[:, np.argsort(-eig_u.real)[:L]]
|
162 | 1 | tkerber | |
163 | 1 | tkerber | U_ww[M:] = np.dot(dag(C_ul), proj_uw) |
164 | 1 | tkerber | else:
|
165 | 1 | tkerber | C_ul = np.empty((Nb - M, 0))
|
166 | 1 | tkerber | |
167 | 1 | tkerber | normalize(C_ul) |
168 | 1 | tkerber | if ortho:
|
169 | 1 | tkerber | lowdin(U_ww) |
170 | 1 | tkerber | else:
|
171 | 1 | tkerber | normalize(U_ww) |
172 | 1 | tkerber | |
173 | 1 | tkerber | return U_ww, C_ul
|
174 | 1 | tkerber | |
175 | 1 | tkerber | |
176 | 1 | tkerber | class Wannier: |
177 | 1 | tkerber | """Maximally localized Wannier Functions
|
178 | 1 | tkerber |
|
179 | 1 | tkerber | Find the set of maximally localized Wannier functions using the
|
180 | 1 | tkerber | spread functional of Marzari and Vanderbilt (PRB 56, 1997 page
|
181 | 1 | tkerber | 12847).
|
182 | 1 | tkerber | """
|
183 | 1 | tkerber | |
184 | 1 | tkerber | def __init__(self, nwannier, calc, |
185 | 1 | tkerber | file=None,
|
186 | 1 | tkerber | nbands=None,
|
187 | 1 | tkerber | fixedenergy=None,
|
188 | 1 | tkerber | fixedstates=None,
|
189 | 1 | tkerber | spin=0,
|
190 | 1 | tkerber | initialwannier='random',
|
191 | 1 | tkerber | seed=None,
|
192 | 1 | tkerber | verbose=False):
|
193 | 1 | tkerber | """
|
194 | 1 | tkerber | Required arguments:
|
195 | 1 | tkerber |
|
196 | 1 | tkerber | ``nwannier``: The number of Wannier functions you wish to construct.
|
197 | 1 | tkerber | This must be at least half the number of electrons in the system
|
198 | 1 | tkerber | and at most equal to the number of bands in the calculation.
|
199 | 1 | tkerber |
|
200 | 1 | tkerber | ``calc``: A converged DFT calculator class.
|
201 | 1 | tkerber | If ``file`` arg. is not provided, the calculator *must* provide the
|
202 | 1 | tkerber | method ``get_wannier_localization_matrix``, and contain the
|
203 | 1 | tkerber | wavefunctions (save files with only the density is not enough).
|
204 | 1 | tkerber | If the localization matrix is read from file, this is not needed,
|
205 | 1 | tkerber | unless ``get_function`` or ``write_cube`` is called.
|
206 | 1 | tkerber |
|
207 | 1 | tkerber | Optional arguments:
|
208 | 1 | tkerber |
|
209 | 1 | tkerber | ``nbands``: Bands to include in localization.
|
210 | 1 | tkerber | The number of bands considered by Wannier can be smaller than the
|
211 | 1 | tkerber | number of bands in the calculator. This is useful if the highest
|
212 | 1 | tkerber | bands of the DFT calculation are not well converged.
|
213 | 1 | tkerber |
|
214 | 1 | tkerber | ``spin``: The spin channel to be considered.
|
215 | 1 | tkerber | The Wannier code treats each spin channel independently.
|
216 | 1 | tkerber |
|
217 | 1 | tkerber | ``fixedenergy`` / ``fixedstates``: Fixed part of Heilbert space.
|
218 | 1 | tkerber | Determine the fixed part of Hilbert space by either a maximal
|
219 | 1 | tkerber | energy *or* a number of bands (possibly a list for multiple
|
220 | 1 | tkerber | k-points).
|
221 | 1 | tkerber | Default is None meaning that the number of fixed states is equated
|
222 | 1 | tkerber | to ``nwannier``.
|
223 | 1 | tkerber |
|
224 | 1 | tkerber | ``file``: Read localization and rotation matrices from this file.
|
225 | 1 | tkerber |
|
226 | 1 | tkerber | ``initialwannier``: Initial guess for Wannier rotation matrix.
|
227 | 1 | tkerber | Can be 'bloch' to start from the Bloch states, 'random' to be
|
228 | 1 | tkerber | randomized, or a list passed to calc.get_initial_wannier.
|
229 | 1 | tkerber |
|
230 | 1 | tkerber | ``seed``: Seed for random ``initialwannier``.
|
231 | 1 | tkerber |
|
232 | 1 | tkerber | ``verbose``: True / False level of verbosity.
|
233 | 1 | tkerber | """
|
234 | 1 | tkerber | # Bloch phase sign convention
|
235 | 1 | tkerber | sign = -1
|
236 | 1 | tkerber | classname = calc.__class__.__name__ |
237 | 1 | tkerber | if classname in ['Dacapo', 'Jacapo']: |
238 | 1 | tkerber | print 'Using ' + classname |
239 | 1 | tkerber | sign = +1
|
240 | 1 | tkerber | |
241 | 1 | tkerber | self.nwannier = nwannier
|
242 | 1 | tkerber | self.calc = calc
|
243 | 1 | tkerber | self.spin = spin
|
244 | 1 | tkerber | self.verbose = verbose
|
245 | 1 | tkerber | self.kpt_kc = sign * calc.get_ibz_k_points()
|
246 | 1 | tkerber | assert len(calc.get_bz_k_points()) == len(self.kpt_kc) |
247 | 1 | tkerber | |
248 | 1 | tkerber | self.kptgrid = get_monkhorst_shape(self.kpt_kc) |
249 | 1 | tkerber | self.Nk = len(self.kpt_kc) |
250 | 1 | tkerber | self.unitcell_cc = calc.get_atoms().get_cell()
|
251 | 1 | tkerber | self.largeunitcell_cc = (self.unitcell_cc.T * self.kptgrid).T |
252 | 1 | tkerber | self.weight_d, self.Gdir_dc = calculate_weights(self.largeunitcell_cc) |
253 | 1 | tkerber | self.Ndir = len(self.weight_d) # Number of directions |
254 | 1 | tkerber | |
255 | 1 | tkerber | if nbands is not None: |
256 | 1 | tkerber | self.nbands = nbands
|
257 | 1 | tkerber | else:
|
258 | 1 | tkerber | self.nbands = calc.get_number_of_bands()
|
259 | 1 | tkerber | if fixedenergy is None: |
260 | 1 | tkerber | if fixedstates is None: |
261 | 1 | tkerber | self.fixedstates_k = np.array([nwannier] * self.Nk, int) |
262 | 1 | tkerber | else:
|
263 | 1 | tkerber | if type(fixedstates) is int: |
264 | 1 | tkerber | fixedstates = [fixedstates] * self.Nk
|
265 | 1 | tkerber | self.fixedstates_k = np.array(fixedstates, int) |
266 | 1 | tkerber | else:
|
267 | 1 | tkerber | # Setting number of fixed states and EDF from specified energy.
|
268 | 1 | tkerber | # All states below this energy (relative to Fermi level) are fixed.
|
269 | 1 | tkerber | fixedenergy += calc.get_fermi_level() |
270 | 1 | tkerber | print fixedenergy
|
271 | 1 | tkerber | self.fixedstates_k = np.array(
|
272 | 1 | tkerber | [calc.get_eigenvalues(k, spin).searchsorted(fixedenergy) |
273 | 1 | tkerber | for k in range(self.Nk)], int) |
274 | 1 | tkerber | self.edf_k = self.nwannier - self.fixedstates_k |
275 | 1 | tkerber | if verbose:
|
276 | 1 | tkerber | print 'Wannier: Fixed states : %s' % self.fixedstates_k |
277 | 1 | tkerber | print 'Wannier: Extra degrees of freedom: %s' % self.edf_k |
278 | 1 | tkerber | |
279 | 1 | tkerber | # Set the list of neighboring k-points k1, and the "wrapping" k0,
|
280 | 1 | tkerber | # such that k1 - k - G + k0 = 0
|
281 | 1 | tkerber | #
|
282 | 1 | tkerber | # Example: kpoints = (-0.375,-0.125,0.125,0.375), dir=0
|
283 | 1 | tkerber | # G = [0.25,0,0]
|
284 | 1 | tkerber | # k=0.375, k1= -0.375 : -0.375-0.375-0.25 => k0=[1,0,0]
|
285 | 1 | tkerber | #
|
286 | 1 | tkerber | # For a gamma point calculation k1 = k = 0, k0 = [1,0,0] for dir=0
|
287 | 1 | tkerber | if self.Nk == 1: |
288 | 1 | tkerber | self.kklst_dk = np.zeros((self.Ndir, 1), int) |
289 | 1 | tkerber | k0_dkc = self.Gdir_dc.reshape(-1, 1, 3) |
290 | 1 | tkerber | else:
|
291 | 1 | tkerber | self.kklst_dk = np.empty((self.Ndir, self.Nk), int) |
292 | 1 | tkerber | k0_dkc = np.empty((self.Ndir, self.Nk, 3), int) |
293 | 1 | tkerber | |
294 | 1 | tkerber | # Distance between kpoints
|
295 | 1 | tkerber | kdist_c = np.empty(3)
|
296 | 1 | tkerber | for c in range(3): |
297 | 1 | tkerber | # make a sorted list of the kpoint values in this direction
|
298 | 1 | tkerber | slist = np.argsort(self.kpt_kc[:, c], kind='mergesort') |
299 | 1 | tkerber | skpoints_kc = np.take(self.kpt_kc, slist, axis=0) |
300 | 1 | tkerber | kdist_c[c] = max([skpoints_kc[n + 1, c] - skpoints_kc[n, c] |
301 | 1 | tkerber | for n in range(self.Nk - 1)]) |
302 | 1 | tkerber | |
303 | 1 | tkerber | for d, Gdir_c in enumerate(self.Gdir_dc): |
304 | 1 | tkerber | for k, k_c in enumerate(self.kpt_kc): |
305 | 1 | tkerber | # setup dist vector to next kpoint
|
306 | 1 | tkerber | G_c = np.where(Gdir_c > 0, kdist_c, 0) |
307 | 1 | tkerber | if max(G_c) < 1e-4: |
308 | 1 | tkerber | self.kklst_dk[d, k] = k
|
309 | 1 | tkerber | k0_dkc[d, k] = Gdir_c |
310 | 1 | tkerber | else:
|
311 | 1 | tkerber | self.kklst_dk[d, k], k0_dkc[d, k] = \
|
312 | 1 | tkerber | neighbor_k_search(k_c, G_c, self.kpt_kc)
|
313 | 1 | tkerber | |
314 | 1 | tkerber | # Set the inverse list of neighboring k-points
|
315 | 1 | tkerber | self.invkklst_dk = np.empty((self.Ndir, self.Nk), int) |
316 | 1 | tkerber | for d in range(self.Ndir): |
317 | 1 | tkerber | for k1 in range(self.Nk): |
318 | 1 | tkerber | self.invkklst_dk[d, k1] = self.kklst_dk[d].tolist().index(k1) |
319 | 1 | tkerber | |
320 | 1 | tkerber | Nw = self.nwannier
|
321 | 1 | tkerber | Nb = self.nbands
|
322 | 1 | tkerber | self.Z_dkww = np.empty((self.Ndir, self.Nk, Nw, Nw), complex) |
323 | 1 | tkerber | self.V_knw = np.zeros((self.Nk, Nb, Nw), complex) |
324 | 1 | tkerber | if file is None: |
325 | 1 | tkerber | self.Z_dknn = np.empty((self.Ndir, self.Nk, Nb, Nb), complex) |
326 | 1 | tkerber | for d, dirG in enumerate(self.Gdir_dc): |
327 | 1 | tkerber | for k in range(self.Nk): |
328 | 1 | tkerber | k1 = self.kklst_dk[d, k]
|
329 | 1 | tkerber | k0_c = k0_dkc[d, k] |
330 | 1 | tkerber | self.Z_dknn[d, k] = calc.get_wannier_localization_matrix(
|
331 | 1 | tkerber | nbands=Nb, dirG=dirG, kpoint=k, nextkpoint=k1, |
332 | 1 | tkerber | G_I=k0_c, spin=self.spin)
|
333 | 1 | tkerber | self.initialize(file=file, initialwannier=initialwannier, seed=seed) |
334 | 1 | tkerber | |
335 | 1 | tkerber | def initialize(self, file=None, initialwannier='random', seed=None): |
336 | 1 | tkerber | """Re-initialize current rotation matrix.
|
337 | 1 | tkerber |
|
338 | 1 | tkerber | Keywords are identical to those of the constructor.
|
339 | 1 | tkerber | """
|
340 | 1 | tkerber | Nw = self.nwannier
|
341 | 1 | tkerber | Nb = self.nbands
|
342 | 1 | tkerber | |
343 | 1 | tkerber | if file is not None: |
344 | 1 | tkerber | self.Z_dknn, self.U_kww, self.C_kul = load(paropen(file)) |
345 | 1 | tkerber | elif initialwannier == 'bloch': |
346 | 1 | tkerber | # Set U and C to pick the lowest Bloch states
|
347 | 1 | tkerber | self.U_kww = np.zeros((self.Nk, Nw, Nw), complex) |
348 | 1 | tkerber | self.C_kul = []
|
349 | 1 | tkerber | for U, M, L in zip(self.U_kww, self.fixedstates_k, self.edf_k): |
350 | 1 | tkerber | U[:] = np.identity(Nw, complex)
|
351 | 1 | tkerber | if L > 0: |
352 | 1 | tkerber | self.C_kul.append(
|
353 | 1 | tkerber | np.identity(Nb - M, complex)[:, :L])
|
354 | 1 | tkerber | else:
|
355 | 1 | tkerber | self.C_kul.append([])
|
356 | 1 | tkerber | elif initialwannier == 'random': |
357 | 1 | tkerber | # Set U and C to random (orthogonal) matrices
|
358 | 1 | tkerber | self.U_kww = np.zeros((self.Nk, Nw, Nw), complex) |
359 | 1 | tkerber | self.C_kul = []
|
360 | 1 | tkerber | for U, M, L in zip(self.U_kww, self.fixedstates_k, self.edf_k): |
361 | 1 | tkerber | U[:] = random_orthogonal_matrix(Nw, seed, real=False)
|
362 | 1 | tkerber | if L > 0: |
363 | 1 | tkerber | self.C_kul.append(random_orthogonal_matrix(
|
364 | 1 | tkerber | Nb - M, seed=seed, real=False)[:, :L])
|
365 | 1 | tkerber | else:
|
366 | 1 | tkerber | self.C_kul.append(np.array([]))
|
367 | 1 | tkerber | else:
|
368 | 1 | tkerber | # Use initial guess to determine U and C
|
369 | 1 | tkerber | self.C_kul, self.U_kww = self.calc.initial_wannier( |
370 | 1 | tkerber | initialwannier, self.kptgrid, self.fixedstates_k, |
371 | 1 | tkerber | self.edf_k, self.spin) |
372 | 1 | tkerber | self.update()
|
373 | 1 | tkerber | |
374 | 1 | tkerber | def save(self, file): |
375 | 1 | tkerber | """Save information on localization and rotation matrices to file."""
|
376 | 1 | tkerber | dump((self.Z_dknn, self.U_kww, self.C_kul), paropen(file, 'w')) |
377 | 1 | tkerber | |
378 | 1 | tkerber | def update(self): |
379 | 1 | tkerber | # Update large rotation matrix V (from rotation U and coeff C)
|
380 | 1 | tkerber | for k, M in enumerate(self.fixedstates_k): |
381 | 1 | tkerber | self.V_knw[k, :M] = self.U_kww[k, :M] |
382 | 1 | tkerber | if M < self.nwannier: |
383 | 1 | tkerber | self.V_knw[k, M:] = np.dot(self.C_kul[k], self.U_kww[k, M:]) |
384 | 1 | tkerber | # else: self.V_knw[k, M:] = 0.0
|
385 | 1 | tkerber | |
386 | 1 | tkerber | # Calculate the Zk matrix from the large rotation matrix:
|
387 | 1 | tkerber | # Zk = V^d[k] Zbloch V[k1]
|
388 | 1 | tkerber | for d in range(self.Ndir): |
389 | 1 | tkerber | for k in range(self.Nk): |
390 | 1 | tkerber | k1 = self.kklst_dk[d, k]
|
391 | 1 | tkerber | self.Z_dkww[d, k] = np.dot(dag(self.V_knw[k]), np.dot( |
392 | 1 | tkerber | self.Z_dknn[d, k], self.V_knw[k1])) |
393 | 1 | tkerber | |
394 | 1 | tkerber | # Update the new Z matrix
|
395 | 1 | tkerber | self.Z_dww = self.Z_dkww.sum(axis=1) / self.Nk |
396 | 1 | tkerber | |
397 | 1 | tkerber | def get_centers(self, scaled=False): |
398 | 1 | tkerber | """Calculate the Wannier centers
|
399 | 1 | tkerber |
|
400 | 1 | tkerber | ::
|
401 | 1 | tkerber |
|
402 | 1 | tkerber | pos = L / 2pi * phase(diag(Z))
|
403 | 1 | tkerber | """
|
404 | 1 | tkerber | coord_wc = np.angle(self.Z_dww[:3].diagonal(0, 1, 2)).T / (2 * pi) % 1 |
405 | 1 | tkerber | if not scaled: |
406 | 1 | tkerber | coord_wc = np.dot(coord_wc, self.largeunitcell_cc)
|
407 | 1 | tkerber | return coord_wc
|
408 | 1 | tkerber | |
409 | 1 | tkerber | def get_radii(self): |
410 | 1 | tkerber | """Calculate the spread of the Wannier functions.
|
411 | 1 | tkerber |
|
412 | 1 | tkerber | ::
|
413 | 1 | tkerber |
|
414 | 1 | tkerber | -- / L \ 2 2
|
415 | 1 | tkerber | radius**2 = - > | --- | ln |Z|
|
416 | 1 | tkerber | --d \ 2pi /
|
417 | 1 | tkerber | """
|
418 | 1 | tkerber | r2 = -np.dot(self.largeunitcell_cc.diagonal()**2 / (2 * pi)**2, |
419 | 1 | tkerber | np.log(abs(self.Z_dww[:3].diagonal(0, 1, 2))**2)) |
420 | 1 | tkerber | return np.sqrt(r2)
|
421 | 1 | tkerber | |
422 | 1 | tkerber | def get_spectral_weight(self, w): |
423 | 1 | tkerber | return abs(self.V_knw[:, :, w])**2 / self.Nk |
424 | 1 | tkerber | |
425 | 1 | tkerber | def get_pdos(self, w, energies, width): |
426 | 1 | tkerber | """Projected density of states (PDOS).
|
427 | 1 | tkerber |
|
428 | 1 | tkerber | Returns the (PDOS) for Wannier function ``w``. The calculation
|
429 | 1 | tkerber | is performed over the energy grid specified in energies. The
|
430 | 1 | tkerber | PDOS is produced as a sum of Gaussians centered at the points
|
431 | 1 | tkerber | of the energy grid and with the specified width.
|
432 | 1 | tkerber | """
|
433 | 1 | tkerber | spec_kn = self.get_spectral_weight(w)
|
434 | 1 | tkerber | dos = np.zeros(len(energies))
|
435 | 1 | tkerber | for k, spec_n in enumerate(spec_kn): |
436 | 1 | tkerber | eig_n = self.calc.get_eigenvalues(k=kpt, s=self.spin) |
437 | 1 | tkerber | for weight, eig in zip(spec_n, eig): |
438 | 1 | tkerber | # Add gaussian centered at the eigenvalue
|
439 | 1 | tkerber | x = ((energies - center) / width)**2
|
440 | 1 | tkerber | dos += weight * np.exp(-x.clip(0., 40.)) / (sqrt(pi) * width) |
441 | 1 | tkerber | return dos
|
442 | 1 | tkerber | |
443 | 1 | tkerber | def max_spread(self, directions=[0, 1, 2]): |
444 | 1 | tkerber | """Returns the index of the most delocalized Wannier function
|
445 | 1 | tkerber | together with the value of the spread functional"""
|
446 | 1 | tkerber | d = np.zeros(self.nwannier)
|
447 | 1 | tkerber | for dir in directions: |
448 | 1 | tkerber | d[dir] = np.abs(self.Z_dww[dir].diagonal())**2 *self.weight_d[dir] |
449 | 1 | tkerber | index = np.argsort(d)[0]
|
450 | 1 | tkerber | print 'Index:', index |
451 | 1 | tkerber | print 'Spread:', d[index] |
452 | 1 | tkerber | |
453 | 1 | tkerber | def translate(self, w, R): |
454 | 1 | tkerber | """Translate the w'th Wannier function
|
455 | 1 | tkerber |
|
456 | 1 | tkerber | The distance vector R = [n1, n2, n3], is in units of the basis
|
457 | 1 | tkerber | vectors of the small cell.
|
458 | 1 | tkerber | """
|
459 | 1 | tkerber | for kpt_c, U_ww in zip(self.kpt_kc, self.U_kww): |
460 | 1 | tkerber | U_ww[:, w] *= np.exp(2.j * pi * np.dot(np.array(R), kpt_c))
|
461 | 1 | tkerber | self.update()
|
462 | 1 | tkerber | |
463 | 1 | tkerber | def translate_to_cell(self, w, cell): |
464 | 1 | tkerber | """Translate the w'th Wannier function to specified cell"""
|
465 | 1 | tkerber | scaled_c = np.angle(self.Z_dww[:3, w, w]) * self.kptgrid / (2 * pi) |
466 | 1 | tkerber | trans = np.array(cell) - np.floor(scaled_c) |
467 | 1 | tkerber | self.translate(w, trans)
|
468 | 1 | tkerber | |
469 | 1 | tkerber | def translate_all_to_cell(self, cell=[0, 0, 0]): |
470 | 1 | tkerber | """Translate all Wannier functions to specified cell.
|
471 | 1 | tkerber |
|
472 | 1 | tkerber | Move all Wannier orbitals to a specific unit cell. There
|
473 | 1 | tkerber | exists an arbitrariness in the positions of the Wannier
|
474 | 1 | tkerber | orbitals relative to the unit cell. This method can move all
|
475 | 1 | tkerber | orbitals to the unit cell specified by ``cell``. For a
|
476 | 1 | tkerber | `\Gamma`-point calculation, this has no effect. For a
|
477 | 1 | tkerber | **k**-point calculation the periodicity of the orbitals are
|
478 | 1 | tkerber | given by the large unit cell defined by repeating the original
|
479 | 1 | tkerber | unitcell by the number of **k**-points in each direction. In
|
480 | 1 | tkerber | this case it is usefull to move the orbitals away from the
|
481 | 1 | tkerber | boundaries of the large cell before plotting them. For a bulk
|
482 | 1 | tkerber | calculation with, say 10x10x10 **k** points, one could move
|
483 | 1 | tkerber | the orbitals to the cell [2,2,2]. In this way the pbc
|
484 | 1 | tkerber | boundary conditions will not be noticed.
|
485 | 1 | tkerber | """
|
486 | 1 | tkerber | scaled_wc = np.angle(self.Z_dww[:3].diagonal(0, 1, 2)).T * \ |
487 | 1 | tkerber | self.kptgrid / (2 * pi) |
488 | 1 | tkerber | trans_wc = np.array(cell)[None] - np.floor(scaled_wc)
|
489 | 1 | tkerber | for kpt_c, U_ww in zip(self.kpt_kc, self.U_kww): |
490 | 1 | tkerber | U_ww *= np.exp(2.j * pi * np.dot(trans_wc, kpt_c))
|
491 | 1 | tkerber | self.update()
|
492 | 1 | tkerber | |
493 | 1 | tkerber | def distances(self, R): |
494 | 1 | tkerber | Nw = self.nwannier
|
495 | 1 | tkerber | cen = self.get_centers()
|
496 | 1 | tkerber | r1 = cen.repeat(Nw, axis=0).reshape(Nw, Nw, 3) |
497 | 1 | tkerber | r2 = cen.copy() |
498 | 1 | tkerber | for i in range(3): |
499 | 1 | tkerber | r2 += self.unitcell_cc[i] * R[i]
|
500 | 1 | tkerber | |
501 | 1 | tkerber | r2 = np.swapaxes(r2.repeat(Nw, axis=0).reshape(Nw, Nw, 3), 0, 1) |
502 | 1 | tkerber | return np.sqrt(np.sum((r1 - r2)**2, axis=-1)) |
503 | 1 | tkerber | |
504 | 1 | tkerber | def get_hopping(self, R): |
505 | 1 | tkerber | """Returns the matrix H(R)_nm=<0,n|H|R,m>.
|
506 | 1 | tkerber |
|
507 | 1 | tkerber | ::
|
508 | 1 | tkerber |
|
509 | 1 | tkerber | 1 _ -ik.R
|
510 | 1 | tkerber | H(R) = <0,n|H|R,m> = --- >_ e H(k)
|
511 | 1 | tkerber | Nk k
|
512 | 1 | tkerber |
|
513 | 1 | tkerber | where R is the cell-distance (in units of the basis vectors of
|
514 | 1 | tkerber | the small cell) and n,m are indices of the Wannier functions.
|
515 | 1 | tkerber | """
|
516 | 1 | tkerber | H_ww = np.zeros([self.nwannier, self.nwannier], complex) |
517 | 1 | tkerber | for k, kpt_c in enumerate(self.kpt_kc): |
518 | 1 | tkerber | phase = np.exp(-2.j * pi * np.dot(np.array(R), kpt_c))
|
519 | 1 | tkerber | H_ww += self.get_hamiltonian(k) * phase
|
520 | 1 | tkerber | return H_ww / self.Nk |
521 | 1 | tkerber | |
522 | 1 | tkerber | def get_hamiltonian(self, k=0): |
523 | 1 | tkerber | """Get Hamiltonian at existing k-vector of index k
|
524 | 1 | tkerber |
|
525 | 1 | tkerber | ::
|
526 | 1 | tkerber |
|
527 | 1 | tkerber | dag
|
528 | 1 | tkerber | H(k) = V diag(eps ) V
|
529 | 1 | tkerber | k k k
|
530 | 1 | tkerber | """
|
531 | 1 | tkerber | eps_n = self.calc.get_eigenvalues(kpt=k, spin=self.spin) |
532 | 1 | tkerber | return np.dot(dag(self.V_knw[k]) * eps_n, self.V_knw[k]) |
533 | 1 | tkerber | |
534 | 1 | tkerber | def get_hamiltonian_kpoint(self, kpt_c): |
535 | 1 | tkerber | """Get Hamiltonian at some new arbitrary k-vector
|
536 | 1 | tkerber |
|
537 | 1 | tkerber | ::
|
538 | 1 | tkerber |
|
539 | 1 | tkerber | _ ik.R
|
540 | 1 | tkerber | H(k) = >_ e H(R)
|
541 | 1 | tkerber | R
|
542 | 1 | tkerber |
|
543 | 1 | tkerber | Warning: This method moves all Wannier functions to cell (0, 0, 0)
|
544 | 1 | tkerber | """
|
545 | 1 | tkerber | if self.verbose: |
546 | 1 | tkerber | print 'Translating all Wannier functions to cell (0, 0, 0)' |
547 | 1 | tkerber | self.translate_all_to_cell()
|
548 | 1 | tkerber | max = (self.kptgrid - 1) / 2 |
549 | 1 | tkerber | max += max > 0 |
550 | 1 | tkerber | N1, N2, N3 = max
|
551 | 1 | tkerber | Hk = np.zeros([self.nwannier, self.nwannier], complex) |
552 | 1 | tkerber | for n1 in xrange(-N1, N1 + 1): |
553 | 1 | tkerber | for n2 in xrange(-N2, N2 + 1): |
554 | 1 | tkerber | for n3 in xrange(-N3, N3 + 1): |
555 | 1 | tkerber | R = np.array([n1, n2, n3], float)
|
556 | 1 | tkerber | hop_ww = self.get_hopping(R)
|
557 | 1 | tkerber | phase = np.exp(+2.j * pi * np.dot(R, kpt_c))
|
558 | 1 | tkerber | Hk += hop_ww * phase |
559 | 1 | tkerber | return Hk
|
560 | 1 | tkerber | |
561 | 1 | tkerber | def get_function(self, index, repeat=None): |
562 | 1 | tkerber | """Get Wannier function on grid.
|
563 | 1 | tkerber |
|
564 | 1 | tkerber | Returns an array with the funcion values of the indicated Wannier
|
565 | 1 | tkerber | function on a grid with the size of the *repeated* unit cell.
|
566 | 1 | tkerber |
|
567 | 1 | tkerber | For a calculation using **k**-points the relevant unit cell for
|
568 | 1 | tkerber | eg. visualization of the Wannier orbitals is not the original unit
|
569 | 1 | tkerber | cell, but rather a larger unit cell defined by repeating the
|
570 | 1 | tkerber | original unit cell by the number of **k**-points in each direction.
|
571 | 1 | tkerber | Note that for a `\Gamma`-point calculation the large unit cell
|
572 | 1 | tkerber | coinsides with the original unit cell.
|
573 | 1 | tkerber | The large unitcell also defines the periodicity of the Wannier
|
574 | 1 | tkerber | orbitals.
|
575 | 1 | tkerber |
|
576 | 1 | tkerber | ``index`` can be either a single WF or a coordinate vector in terms
|
577 | 1 | tkerber | of the WFs.
|
578 | 1 | tkerber | """
|
579 | 1 | tkerber | |
580 | 1 | tkerber | # Default size of plotting cell is the one corresponding to k-points.
|
581 | 1 | tkerber | if repeat is None: |
582 | 1 | tkerber | repeat = self.kptgrid
|
583 | 1 | tkerber | N1, N2, N3 = repeat |
584 | 1 | tkerber | |
585 | 1 | tkerber | dim = self.calc.get_number_of_grid_points()
|
586 | 1 | tkerber | largedim = dim * [N1, N2, N3] |
587 | 1 | tkerber | |
588 | 1 | tkerber | wanniergrid = np.zeros(largedim, dtype=complex)
|
589 | 1 | tkerber | for k, kpt_c in enumerate(self.kpt_kc): |
590 | 1 | tkerber | # The coordinate vector of wannier functions
|
591 | 1 | tkerber | if type(index) == int: |
592 | 1 | tkerber | vec_n = self.V_knw[k, :, index]
|
593 | 1 | tkerber | else:
|
594 | 1 | tkerber | vec_n = np.dot(self.V_knw[k], index)
|
595 | 1 | tkerber | |
596 | 1 | tkerber | wan_G = np.zeros(dim, complex)
|
597 | 1 | tkerber | for n, coeff in enumerate(vec_n): |
598 | 1 | tkerber | wan_G += coeff * self.calc.get_pseudo_wave_function(
|
599 | 1 | tkerber | n, k, self.spin, pad=True) |
600 | 1 | tkerber | |
601 | 1 | tkerber | # Distribute the small wavefunction over large cell:
|
602 | 1 | tkerber | for n1 in xrange(N1): |
603 | 1 | tkerber | for n2 in xrange(N2): |
604 | 1 | tkerber | for n3 in xrange(N3): # sign? |
605 | 1 | tkerber | e = np.exp(-2.j * pi * np.dot([n1, n2, n3], kpt_c))
|
606 | 1 | tkerber | wanniergrid[n1 * dim[0]:(n1 + 1) * dim[0], |
607 | 1 | tkerber | n2 * dim[1]:(n2 + 1) * dim[1], |
608 | 1 | tkerber | n3 * dim[2]:(n3 + 1) * dim[2]] += e * wan_G |
609 | 1 | tkerber | |
610 | 1 | tkerber | # Normalization
|
611 | 1 | tkerber | wanniergrid /= np.sqrt(self.Nk)
|
612 | 1 | tkerber | return wanniergrid
|
613 | 1 | tkerber | |
614 | 1 | tkerber | def write_cube(self, index, fname, repeat=None, real=True): |
615 | 1 | tkerber | """Dump specified Wannier function to a cube file"""
|
616 | 1 | tkerber | from ase.io.cube import write_cube |
617 | 1 | tkerber | |
618 | 1 | tkerber | # Default size of plotting cell is the one corresponding to k-points.
|
619 | 1 | tkerber | if repeat is None: |
620 | 1 | tkerber | repeat = self.kptgrid
|
621 | 1 | tkerber | atoms = self.calc.get_atoms() * repeat
|
622 | 1 | tkerber | func = self.get_function(index, repeat)
|
623 | 1 | tkerber | |
624 | 1 | tkerber | # Handle separation of complex wave into real parts
|
625 | 1 | tkerber | if real:
|
626 | 1 | tkerber | if self.Nk == 1: |
627 | 1 | tkerber | func *= np.exp(-1.j * np.angle(func.max()))
|
628 | 1 | tkerber | if 0: assert max(abs(func.imag).flat) < 1e-4 |
629 | 1 | tkerber | func = func.real |
630 | 1 | tkerber | else:
|
631 | 1 | tkerber | func = abs(func)
|
632 | 1 | tkerber | else:
|
633 | 1 | tkerber | phase_fname = fname.split('.')
|
634 | 1 | tkerber | phase_fname.insert(1, 'phase') |
635 | 1 | tkerber | phase_fname = '.'.join(phase_fname)
|
636 | 1 | tkerber | write_cube(phase_fname, atoms, data=np.angle(func)) |
637 | 1 | tkerber | func = abs(func)
|
638 | 1 | tkerber | |
639 | 1 | tkerber | write_cube(fname, atoms, data=func) |
640 | 1 | tkerber | |
641 | 1 | tkerber | def localize(self, step=0.25, tolerance=1e-08, |
642 | 1 | tkerber | updaterot=True, updatecoeff=True): |
643 | 1 | tkerber | """Optimize rotation to give maximal localization"""
|
644 | 1 | tkerber | md_min(self, step, tolerance, verbose=self.verbose, |
645 | 1 | tkerber | updaterot=updaterot, updatecoeff=updatecoeff) |
646 | 1 | tkerber | |
647 | 1 | tkerber | def get_functional_value(self): |
648 | 1 | tkerber | """Calculate the value of the spread functional.
|
649 | 1 | tkerber |
|
650 | 1 | tkerber | ::
|
651 | 1 | tkerber |
|
652 | 1 | tkerber | Tr[|ZI|^2]=sum(I)sum(n) w_i|Z_(i)_nn|^2,
|
653 | 1 | tkerber |
|
654 | 1 | tkerber | where w_i are weights."""
|
655 | 1 | tkerber | a_d = np.sum(np.abs(self.Z_dww.diagonal(0, 1, 2))**2, axis=1) |
656 | 1 | tkerber | return np.dot(a_d, self.weight_d).real |
657 | 1 | tkerber | |
658 | 1 | tkerber | def get_gradients(self): |
659 | 1 | tkerber | # Determine gradient of the spread functional.
|
660 | 1 | tkerber | #
|
661 | 1 | tkerber | # The gradient for a rotation A_kij is::
|
662 | 1 | tkerber | #
|
663 | 1 | tkerber | # dU = dRho/dA_{k,i,j} = sum(I) sum(k')
|
664 | 1 | tkerber | # + Z_jj Z_kk',ij^* - Z_ii Z_k'k,ij^*
|
665 | 1 | tkerber | # - Z_ii^* Z_kk',ji + Z_jj^* Z_k'k,ji
|
666 | 1 | tkerber | #
|
667 | 1 | tkerber | # The gradient for a change of coefficients is::
|
668 | 1 | tkerber | #
|
669 | 1 | tkerber | # dRho/da^*_{k,i,j} = sum(I) [[(Z_0)_{k} V_{k'} diag(Z^*) +
|
670 | 1 | tkerber | # (Z_0_{k''})^d V_{k''} diag(Z)] *
|
671 | 1 | tkerber | # U_k^d]_{N+i,N+j}
|
672 | 1 | tkerber | #
|
673 | 1 | tkerber | # where diag(Z) is a square,diagonal matrix with Z_nn in the diagonal,
|
674 | 1 | tkerber | # k' = k + dk and k = k'' + dk.
|
675 | 1 | tkerber | #
|
676 | 1 | tkerber | # The extra degrees of freedom chould be kept orthonormal to the fixed
|
677 | 1 | tkerber | # space, thus we introduce lagrange multipliers, and minimize instead::
|
678 | 1 | tkerber | #
|
679 | 1 | tkerber | # Rho_L=Rho- sum_{k,n,m} lambda_{k,nm} <c_{kn}|c_{km}>
|
680 | 1 | tkerber | #
|
681 | 1 | tkerber | # for this reason the coefficient gradients should be multiplied
|
682 | 1 | tkerber | # by (1 - c c^d).
|
683 | 1 | tkerber | |
684 | 1 | tkerber | Nb = self.nbands
|
685 | 1 | tkerber | Nw = self.nwannier
|
686 | 1 | tkerber | |
687 | 1 | tkerber | dU = [] |
688 | 1 | tkerber | dC = [] |
689 | 1 | tkerber | for k in xrange(self.Nk): |
690 | 1 | tkerber | M = self.fixedstates_k[k]
|
691 | 1 | tkerber | L = self.edf_k[k]
|
692 | 1 | tkerber | U_ww = self.U_kww[k]
|
693 | 1 | tkerber | C_ul = self.C_kul[k]
|
694 | 1 | tkerber | Utemp_ww = np.zeros((Nw, Nw), complex)
|
695 | 1 | tkerber | Ctemp_nw = np.zeros((Nb, Nw), complex)
|
696 | 1 | tkerber | |
697 | 1 | tkerber | for d, weight in enumerate(self.weight_d): |
698 | 1 | tkerber | if abs(weight) < 1.0e-6: |
699 | 1 | tkerber | continue
|
700 | 1 | tkerber | |
701 | 1 | tkerber | Z_knn = self.Z_dknn[d]
|
702 | 1 | tkerber | diagZ_w = self.Z_dww[d].diagonal()
|
703 | 1 | tkerber | Zii_ww = np.repeat(diagZ_w, Nw).reshape(Nw, Nw) |
704 | 1 | tkerber | k1 = self.kklst_dk[d, k]
|
705 | 1 | tkerber | k2 = self.invkklst_dk[d, k]
|
706 | 1 | tkerber | V_knw = self.V_knw
|
707 | 1 | tkerber | Z_kww = self.Z_dkww[d]
|
708 | 1 | tkerber | |
709 | 1 | tkerber | if L > 0: |
710 | 1 | tkerber | Ctemp_nw += weight * np.dot( |
711 | 1 | tkerber | np.dot(Z_knn[k], V_knw[k1]) * diagZ_w.conj() + |
712 | 1 | tkerber | np.dot(dag(Z_knn[k2]), V_knw[k2]) * diagZ_w, |
713 | 1 | tkerber | dag(U_ww)) |
714 | 1 | tkerber | |
715 | 1 | tkerber | temp = Zii_ww.T * Z_kww[k].conj() - Zii_ww * Z_kww[k2].conj() |
716 | 1 | tkerber | Utemp_ww += weight * (temp - dag(temp)) |
717 | 1 | tkerber | dU.append(Utemp_ww.ravel()) |
718 | 1 | tkerber | if L > 0: |
719 | 1 | tkerber | # Ctemp now has same dimension as V, the gradient is in the
|
720 | 1 | tkerber | # lower-right (Nb-M) x L block
|
721 | 1 | tkerber | Ctemp_ul = Ctemp_nw[M:, M:] |
722 | 1 | tkerber | G_ul = Ctemp_ul - np.dot(np.dot(C_ul, dag(C_ul)), Ctemp_ul) |
723 | 1 | tkerber | dC.append(G_ul.ravel()) |
724 | 1 | tkerber | |
725 | 1 | tkerber | return np.concatenate(dU + dC)
|
726 | 1 | tkerber | |
727 | 1 | tkerber | def step(self, dX, updaterot=True, updatecoeff=True): |
728 | 1 | tkerber | # dX is (A, dC) where U->Uexp(-A) and C->C+dC
|
729 | 1 | tkerber | Nw = self.nwannier
|
730 | 1 | tkerber | Nk = self.Nk
|
731 | 1 | tkerber | M_k = self.fixedstates_k
|
732 | 1 | tkerber | L_k = self.edf_k
|
733 | 1 | tkerber | if updaterot:
|
734 | 1 | tkerber | A_kww = dX[:Nk * Nw**2].reshape(Nk, Nw, Nw)
|
735 | 1 | tkerber | for U, A in zip(self.U_kww, A_kww): |
736 | 1 | tkerber | H = -1.j * A.conj()
|
737 | 1 | tkerber | epsilon, Z = np.linalg.eigh(H) |
738 | 1 | tkerber | # Z contains the eigenvectors as COLUMNS.
|
739 | 1 | tkerber | # Since H = iA, dU = exp(-A) = exp(iH) = ZDZ^d
|
740 | 1 | tkerber | dU = np.dot(Z * np.exp(1.j * epsilon), dag(Z))
|
741 | 1 | tkerber | U[:] = np.dot(U, dU) |
742 | 1 | tkerber | |
743 | 1 | tkerber | if updatecoeff:
|
744 | 1 | tkerber | start = 0
|
745 | 1 | tkerber | for C, unocc, L in zip(self.C_kul, self.nbands - M_k, L_k): |
746 | 1 | tkerber | if L == 0 or unocc == 0: |
747 | 1 | tkerber | continue
|
748 | 1 | tkerber | Ncoeff = L * unocc |
749 | 1 | tkerber | deltaC = dX[Nk * Nw**2 + start: Nk * Nw**2 + start + Ncoeff] |
750 | 1 | tkerber | C += deltaC.reshape(unocc, L) |
751 | 1 | tkerber | gram_schmidt(C) |
752 | 1 | tkerber | start += Ncoeff |
753 | 1 | tkerber | |
754 | 1 | tkerber | self.update() |