def read_hamiltonian(self, **kwargs): """ Returns the electronic structure from the siesta.TSHS file """ tshs_g = self.read_geometry() geom = _geometry_align(tshs_g, kwargs.get('geometry', tshs_g), self.__class__, 'read_hamiltonian') # read the sizes used... sizes = _siesta.read_tshs_sizes(self.file) _bin_check(self, 'read_hamiltonian', 'could not read sizes.') isc = _siesta.read_tshs_cell(self.file, sizes[3])[2].T _bin_check(self, 'read_hamiltonian', 'could not read cell.') spin = sizes[0] no = sizes[2] nnz = sizes[4] ncol, col, dH, dS = _siesta.read_tshs_hs(self.file, spin, no, nnz) _bin_check(self, 'read_hamiltonian', 'could not read Hamiltonian and overlap matrix.') # Check whether it is an orthogonal basis set orthogonal = np.abs(dS).sum() == geom.no # Create the Hamiltonian container H = Hamiltonian(geom, spin, nnzpr=1, orthogonal=orthogonal) # Create the new sparse matrix H._csr.ncol = ncol.astype(np.int32, copy=False) H._csr.ptr = np.insert(np.cumsum(ncol, dtype=np.int32), 0, 0) # Correct fortran indices H._csr.col = col.astype(np.int32, copy=False) - 1 H._csr._nnz = len(col) if orthogonal: H._csr._D = _a.emptyd([nnz, spin]) H._csr._D[:, :] = dH[:, :] else: H._csr._D = _a.emptyd([nnz, spin + 1]) H._csr._D[:, :spin] = dH[:, :] H._csr._D[:, spin] = dS[:] _mat_spin_convert(H) # Convert to sisl supercell _csr_from_sc_off(H.geometry, isc, H._csr) # Find all indices where dS == 1 (remember col is in fortran indices) idx = col[np.isclose(dS, 1.).nonzero()[0]] if np.any(idx > no): print('Number of orbitals: {}'.format(no)) print(idx) raise SileError( str(self) + '.read_hamiltonian could not assert ' 'the supercell connections in the primary unit-cell.') return H
def unfold_points(self, k): r""" Return a list of k-points to be evaluated for this objects unfolding The k-point `k` is with respect to the unfolded geometry. The return list of `k` points are the k-points required to be sampled in the folded geometry. Parameters ---------- k : (3,) of float k-point evaluation corresponding to the unfolded unit-cell Returns ------- k_unfold a list of ``np.prod(self.bloch)`` k-points used for the unfolding """ k = _a.arrayd(k) # Create expansion points B = self._bloch unfold = _a.emptyd([B[2], B[1], B[0], 3]) # Use B-casting rules (much simpler) unfold[:, :, :, 0] = (aranged(B[0]).reshape(1, 1, -1) + k[0]) / B[0] unfold[:, :, :, 1] = (aranged(B[1]).reshape(1, -1, 1) + k[1]) / B[1] unfold[:, :, :, 2] = (aranged(B[2]).reshape(-1, 1, 1) + k[2]) / B[2] # Back-transform shape return unfold.reshape(-1, 3)
def read_overlap(self, **kwargs): """ Returns the overlap matrix from the TranSiesta file """ tshs_g = self.read_geometry() geom = _geometry_align(tshs_g, kwargs.get('geometry', tshs_g), self.__class__, 'read_overlap') # read the sizes used... sizes = _siesta.read_tshs_sizes(self.file) _bin_check(self, 'read_overlap', 'could not read sizes.') isc = _siesta.read_tshs_cell(self.file, sizes[3])[2].T _bin_check(self, 'read_overlap', 'could not read cell.') no = sizes[2] nnz = sizes[4] ncol, col, dS = _siesta.read_tshs_s(self.file, no, nnz) _bin_check(self, 'read_overlap', 'could not read overlap matrix.') # Create the Hamiltonian container S = SparseOrbitalBZ(geom, nnzpr=1) # Create the new sparse matrix S._csr.ncol = ncol.astype(np.int32, copy=False) S._csr.ptr = np.insert(np.cumsum(ncol, dtype=np.int32), 0, 0) # Correct fortran indices S._csr.col = col.astype(np.int32, copy=False) - 1 S._csr._nnz = len(col) S._csr._D = _a.emptyd([nnz, 1]) S._csr._D[:, 0] = dS[:] # Convert to sisl supercell _csr_from_sc_off(S.geometry, isc, S._csr) return S
def lineark(self, ticks=False): """ A 1D array which corresponds to the delta-k values of the path This is meant for plotting Examples -------- >>> p = BandStructure(...) # doctest: +SKIP >>> eigs = Hamiltonian.eigh(p) # doctest: +SKIP >>> for i in range(len(Hamiltonian)): # doctest: +SKIP ... plt.plot(p.lineark(), eigs[:, i]) # doctest: +SKIP >>> p = BandStructure(...) # doctest: +SKIP >>> eigs = Hamiltonian.eigh(p) # doctest: +SKIP >>> lk, kt, kl = p.lineark(True) # doctest: +SKIP >>> plt.xticks(kt, kl) # doctest: +SKIP >>> for i in range(len(Hamiltonian)): # doctest: +SKIP ... plt.plot(lk, eigs[:, i]) # doctest: +SKIP Parameters ---------- ticks : bool, optional if `True` the ticks for the points are also returned lk, xticks, label_ticks, lk = BandStructure.lineark(True) Returns ------- linear_k : The positions in reciprocal space determined by the distance between points k_tick : Linear k-positions of the points, only returned if `ticks` is ``True`` k_label : Labels at `k_tick`, only returned if `ticks` is ``True`` """ # Calculate points k = [self.tocartesian(pnt) for pnt in self.point] # Get difference between points dk = np.diff(k, axis=0) # Calculate the cumultative distance between points k_len = np.insert(_a.cumsumd((dk ** 2).sum(1) ** .5), 0, 0.) xtick = [None] * len(k) # Prepare output array dK = _a.emptyd(len(self)) # short-hand ls = np.linspace xtick = np.insert(_a.cumsumi(self.division), 0, 0) for i in range(len(dk)): n = self.division[i] end = i == len(dk) - 1 dK[xtick[i]:xtick[i+1]] = ls(k_len[i], k_len[i+1], n, dtype=np.float64, endpoint=end) xtick[-1] -= 1 # Get label tick, in case self.name is a single string 'ABCD' label_tick = [a for a in self.name] if ticks: return dK, dK[xtick], label_tick return dK
def __init__(self, parent, nkpt, displacement=None, size=None, trs=True): super(MonkhorstPack, self).__init__(parent) if isinstance(nkpt, Integral): nkpt = np.diag([nkpt] * 3) elif isinstance(nkpt[0], Integral): nkpt = np.diag(nkpt) # Now we have a matrix of k-points if np.any(nkpt - np.diag(np.diag(nkpt)) != 0): raise NotImplementedError( self.__class__.__name__ + " with off-diagonal components is not implemented yet") if displacement is None: displacement = np.zeros(3, np.float64) elif isinstance(displacement, Real): displacement = np.zeros(3, np.float64) + displacement if size is None: size = _a.onesd(3) elif isinstance(size, Real): size = _a.zerosd(3) + size # Retrieve the diagonal number of values Dn = np.diag(nkpt).astype(np.int32) if np.any(Dn) == 0: raise ValueError(self.__class__.__name__ + ' *must* be initialized with ' 'diagonal elements different from 0.') i_trs = -1 if trs: # Figure out which direction to TRS nmax = 0 for i in [0, 1, 2]: if displacement[i] == 0. and Dn[i] > nmax: nmax = Dn[i] i_trs = i # Calculate k-points and weights along all directions kw = [ self.grid(Dn[i], displacement[i], size[i], i == i_trs) for i in (0, 1, 2) ] self._k = _a.emptyd((kw[0][0].size, kw[1][0].size, kw[2][0].size, 3)) self._w = _a.onesd(self._k.shape[:-1]) for i in (0, 1, 2): k = kw[i][0].reshape(-1, 1, 1) w = kw[i][1].reshape(-1, 1, 1) self._k[..., i] = np.rollaxis(k, 0, i + 1) self._w[...] *= np.rollaxis(w, 0, i + 1) del kw self._k.shape = (-1, 3) self._k = np.where(self._k > .5, self._k - 1, self._k) self._w.shape = (-1, )
def _index_shape_ellipsoid(self, ellipsoid): """ Internal routine for ellipsoid shape-indices """ # Figure out the points on the ellipsoid rad1 = pi / 180 theta, phi = ogrid[-pi:pi:rad1, 0:pi:rad1] rxyz = _a.emptyd([theta.size, phi.size, 3]) rxyz[..., 2] = cos(phi) sin(phi, out=phi) rxyz[..., 0] = cos(theta) * phi rxyz[..., 1] = sin(theta) * phi rxyz = dot(rxyz, ellipsoid._v) + ellipsoid.center.reshape(1, 3) del theta, phi # Get all indices of the ellipsoid circumference return self.index(rxyz)
def solve_lagrange(self): r""" Calculate the coefficients according to Pulay's method, return everything + Lagrange multiplier """ hist = self.history n_h = len(hist) metric = self._metric if n_h == 0: # Externally the coefficients should reflect the weight per previous iteration. # The mixing weight is an additional parameter return _a.arrayd([1.]), 100. elif n_h == 1: return _a.arrayd([1.]), metric(hist[0][-1], hist[0][-1]) # Initialize the matrix to be solved against B = _a.emptyd([n_h + 1, n_h + 1]) # Fill matrix B for i in range(n_h): ei = hist[i][-1] B[i, i] = metric(ei, ei) for j in range(i + 1, n_h): ej = hist[j][-1] B[i, j] = metric(ei, ej) B[j, i] = B[i, j] B[:, n_h] = 1. B[n_h, :] = 1. B[n_h, n_h] = 0. # Although B contains 1 and a number on the order of # number of elements (self._hist.size), it seems very # numerically stable. # Create RHS RHS = _a.zerosd(n_h + 1) RHS[-1] = 1 try: # Apparently we cannot use assume_a='sym' # Is this because sym also implies positive definitiness? # However, these are matrices of order ~30, so we don't care c = solve(B, RHS) return c[:-1], -c[-1] except np.linalg.LinAlgError as e: # We have a LinalgError return _a.arrayd([1.]), metric(hist[-1][-1], hist[-1][-1])
def spher2cart(r, theta, phi): r""" Convert spherical coordinates to cartesian coordinates Parameters ---------- r : array_like radius theta : array_like azimuthal angle in the :math:`x-y` plane phi : array_like polar angle from the :math:`z` axis """ rx = r * cos(theta) * sin(phi) R = _a.emptyd(rx.shape + (3, )) R[..., 0] = rx del rx R[..., 1] = r * sin(theta) * sin(phi) R[..., 2] = r * cos(phi) return R
def _index_shape_cuboid(self, cuboid): """ Internal routine for cuboid shape-indices """ # Construct all points on the outer rim of the cuboids min_d = fnorm(self.dcell).min() # Retrieve cuboids edge-lengths v = cuboid.edge_length # Create normalized cuboid vectors (because we expan via the lengths below vn = cuboid._v / fnorm(cuboid._v).reshape(-1, 1) LL = (cuboid.center - cuboid._v.sum(0) / 2).reshape(1, 3) UR = (cuboid.center + cuboid._v.sum(0) / 2).reshape(1, 3) # Create coordinates a = vn[0, :].reshape(1, -1) * _a.aranged(0, v[0] + min_d, min_d).reshape(-1, 1) b = vn[1, :].reshape(1, -1) * _a.aranged(0, v[1] + min_d, min_d).reshape(-1, 1) c = vn[2, :].reshape(1, -1) * _a.aranged(0, v[2] + min_d, min_d).reshape(-1, 1) # Now create all sides sa = a.shape[0] sb = b.shape[0] sc = c.shape[0] def plane(v1, v2): return (v1.reshape(-1, 1, 3) + v2.reshape(1, -1, 3)).reshape(1, -1, 3) # Allocate for the 6 faces of the cuboid rxyz = _a.emptyd([2, sa * sb + sa * sc + sb * sc, 3]) # Define the LL and UR rxyz[0, :, :] = LL rxyz[1, :, :] = UR i = 0 rxyz[:, i:i + sa * sb, :] += plane(a, b) i += sa * sb rxyz[:, i:i + sa * sc, :] += plane(a, c) i += sa * sc rxyz[:, i:i + sb * sc, :] += plane(b, c) del a, b, c, sa, sb, sc rxyz.shape = (-1, 3) # Get all indices of the cuboid planes return self.index(rxyz)
def _index_shape(self, shape): """ Internal routine for shape-indices """ # First grab the sphere, subsequent indices will be reduced # by the actual shape cuboid = shape.toCuboid() ellipsoid = shape.toEllipsoid() if ellipsoid.volume() > cuboid.volume(): idx = self._index_shape_cuboid(cuboid) else: idx = self._index_shape_ellipsoid(ellipsoid) # Get min/max imin = idx.min(0) imax = idx.max(0) del idx dc = self.dcell # Now to find the actual points inside the shape # First create all points in the square and then retrieve all indices # within. ix = _a.aranged(imin[0], imax[0] + 0.5) iy = _a.aranged(imin[1], imax[1] + 0.5) iz = _a.aranged(imin[2], imax[2] + 0.5) output_shape = (ix.size, iy.size, iz.size, 3) rxyz = _a.emptyd(output_shape) ao = add.outer ao(ao(ix * dc[0, 0], iy * dc[1, 0]), iz * dc[2, 0], out=rxyz[:, :, :, 0]) ao(ao(ix * dc[0, 1], iy * dc[1, 1]), iz * dc[2, 1], out=rxyz[:, :, :, 1]) ao(ao(ix * dc[0, 2], iy * dc[1, 2]), iz * dc[2, 2], out=rxyz[:, :, :, 2]) idx = shape.within_index(rxyz.reshape(-1, 3)) del rxyz i = _a.emptyi(output_shape) i[:, :, :, 0] = ix.reshape(-1, 1, 1) i[:, :, :, 1] = iy.reshape(1, -1, 1) i[:, :, :, 2] = iz.reshape(1, 1, -1) del ix, iy, iz i.shape = (-1, 3) i = take(i, idx, axis=0) del idx return i
def fermi_level(self, bz=None, q=None, distribution='fermi_dirac', q_tol=1e-10): """ Calculate the Fermi-level using a Brillouinzone sampling and a target charge The Fermi-level will be calculated using an iterative approach by first calculating all eigenvalues and subsequently fitting the Fermi level to the final charge (`q`). Parameters ---------- bz : Brillouinzone, optional sampled k-points and weights, the ``bz.parent`` will be equal to this object upon return default to Gamma-point q : float, list of float, optional seeked charge, if not set will be equal to ``self.geometry.q0``. If a list of two is passed there will be calculated a Fermi-level per spin-channel. If the Hamiltonian is not spin-polarized the sum of the list will be used and only a single fermi-level will be returned. distribution : str, func, optional used distribution, must accept the keyword ``mu`` as parameter for the Fermi-level q_tol : float, optional tolerance of charge for finding the Fermi-level Returns ------- float or array_like the Fermi-level of the system (or two if two different charges are passed) """ if bz is None: # Gamma-point only from .brillouinzone import BrillouinZone bz = BrillouinZone(self) else: # Overwrite the parent in bz bz.set_parent(self) if q is None: if self.spin.is_unpolarized: q = self.geometry.q0 * 0.5 else: q = self.geometry.q0 # Ensure we have an "array" in case of spin-polarized calculations q = _a.asarrayd(q) if np.any(q <= 0.): raise ValueError( f"{self.__class__.__name__}.fermi_level cannot calculate the Fermi level " "for 0 electrons.") if isinstance(distribution, str): distribution = get_distribution(distribution) # B-cast for easier weights w = bz.weight.reshape(-1, 1) # Internal class to calculate the Fermi-level def _Ef(q, eig): # We could reduce it depending on the temperature, # however the distribution does not have the kT # parameter available. min_Ef, max_Ef = eig.min(), eig.max() nextafter = np.nextafter while nextafter(min_Ef, max_Ef) < max_Ef: Ef = (min_Ef + max_Ef) * 0.5 # Calculate guessed charge qt = (distribution(eig, mu=Ef) * w).sum() if abs(qt - q) < q_tol: return Ef if qt >= q: max_Ef = Ef elif qt <= q: min_Ef = Ef return Ef # Retrieve dispatcher for averaging eigh = bz.apply.array.eigh if self.spin.is_polarized and q.size == 2: if np.any(q >= len(self)): raise ValueError( f"{self.__class__.__name__}.fermi_level cannot calculate the Fermi level " "for electrons ({q}) equal to or above number of orbitals ({len(self)})." ) # We need to do Fermi-level separately since the user requests # separate fillings Ef = _a.emptyd(2) Ef[0] = _Ef(q[0], eigh(spin=0)) Ef[1] = _Ef(q[1], eigh(spin=1)) else: # Ensure a single charge q = q.sum() if q >= len(self): raise ValueError( f"{self.__class__.__name__}.fermi_level cannot calculate the Fermi level " "for electrons ({q}) equal to or above number of orbitals ({len(self)})." ) if self.spin.is_polarized: Ef = _Ef(q, np.concatenate( [eigh(spin=0), eigh(spin=1)], axis=1)) else: Ef = _Ef(q, eigh()) return Ef
def read_energy_density_matrix(self, **kwargs): """ Returns the energy density matrix from the siesta.DM file """ # Now read the sizes used... spin, no, nsc, nnz = _siesta.read_tsde_sizes(self.file) _bin_check(self, 'read_energy_density_matrix', 'could not read energy density matrix sizes.') ncol, col, dEDM = _siesta.read_tsde_edm(self.file, spin, no, nsc, nnz) _bin_check(self, 'read_energy_density_matrix', 'could not read energy density matrix.') # Try and immediately attach a geometry geom = kwargs.get('geometry', kwargs.get('geom', None)) if geom is None: # We truly, have no clue, # Just generate a boxed system xyz = [[x, 0, 0] for x in range(no)] sc = SuperCell([no, 1, 1], nsc=nsc) geom = Geometry(xyz, Atom(1), sc=sc) if nsc[0] != 0 and np.any(geom.nsc != nsc): # We have to update the number of supercells! geom.set_nsc(nsc) if geom.no != no: raise SileError( str(self) + '.read_energy_density_matrix could ' 'not use the passed geometry as the number of atoms or orbitals ' 'is inconsistent with DM file.') # Create the energy density matrix container EDM = EnergyDensityMatrix(geom, spin, nnzpr=1, dtype=np.float64, orthogonal=False) # Create the new sparse matrix EDM._csr.ncol = ncol.astype(np.int32, copy=False) EDM._csr.ptr = np.insert(np.cumsum(ncol, dtype=np.int32), 0, 0) # Correct fortran indices EDM._csr.col = col.astype(np.int32, copy=False) - 1 EDM._csr._nnz = len(col) EDM._csr._D = _a.emptyd([nnz, spin + 1]) EDM._csr._D[:, :spin] = dEDM[:, :] # EDM file does not contain overlap matrix... so neglect it for now. EDM._csr._D[:, spin] = 0. _mat_spin_convert(EDM) # Convert the supercells to sisl supercells if nsc[0] != 0 or geom.no_s >= col.max(): _csr_from_siesta(geom, EDM._csr) else: warn( str(self) + '.read_energy_density_matrix may result in a wrong sparse pattern!' ) return EDM
def fermi_level(self, bz=None, q=None, distribution='fermi_dirac', q_tol=1e-12): """ Calculate the Fermi-level using a Brillouinzone sampling and a target charge The Fermi-level will be calculated using an iterative approach by first calculating all eigenvalues and subsequently fitting the Fermi level to the final charge (`q`). Parameters ---------- bz : Brillouinzone, optional sampled k-points and weights, the ``bz.parent`` will be equal to this object upon return default to Gamma-point q : float, list of float, optional seeked charge, if not set will be equal to ``self.geometry.q0``. If a list of two is passed there will be calculated a Fermi-level per spin-channel. If the Hamiltonian is not spin-polarized the sum of the list will be used and only a single fermi-level will be returned. distribution : str, func, optional used distribution, must accept the keyword ``mu`` as parameter for the Fermi-level q_tol : float, optional tolerance of charge for finding the Fermi-level Returns ------- fermi-level : the Fermi-level of the system (or two if two different charges are passed) """ if bz is None: # Gamma-point only from .brillouinzone import BrillouinZone bz = BrillouinZone(self) else: # Overwrite the parent in bz bz.set_parent(self) # Ensure we are using asarray bz.asarray() if q is None: q = self.geometry.q0 # Ensure we have an "array" in case of spin-polarized calculations q = np.asarray(q) if isinstance(distribution, str): distribution = get_distribution(distribution) # B-cast for easier weights w = bz.weight.reshape(-1, 1) # Internal class to calculate the Fermi-level def _Ef(q, eig): # We could reduce it depending on the temperature, # however the distribution does not have the kT # parameter available. Ef = np.average(eig[:, int(q)]) l_Ef = [] l_q = [] def list_append(q, Ef): l_q.append(q) l_Ef.append(Ef) # Calculate guessed charge qt = (distribution(eig, mu=Ef) * w).sum() while abs(qt - q) > q_tol: # Add to cubic-spline list_append(qt, Ef) # Estimate new Fermi-level if len(l_q) > 1: # We can do a spline interpolation lq = np.array(l_q) idx = np.argsort(lq) lEf = np.array(l_Ef) Ef = CubicSpline(lq[idx], lEf[idx], extrapolate=True)(q) else: # Update limits if qt > q: Ef = Ef - 0.5 elif qt < q: Ef = Ef + 0.5 # Calculate new guessed charge qt = (distribution(eig, mu=Ef) * w).sum() return Ef if self.spin.is_polarized and q.size == 2: # We need to do Fermi-level separately since the user requests # separate fillings Ef = _a.emptyd(2) Ef[0] = _Ef(q[0], bz.eigh(spin=0)) Ef[1] = _Ef(q[1], bz.eigh(spin=1)) else: # Ensure a single charge q = q.sum() if self.spin.is_polarized: Ef = _Ef( q, np.concatenate( [bz.eigh(spin=0), bz.eigh(spin=1)], axis=1)) else: Ef = _Ef(q, bz.eigh()) return Ef
def __init__(self, parent, nkpt, displacement=None, size=None, centered=True, trs=True): super(MonkhorstPack, self).__init__(parent) if isinstance(nkpt, Integral): nkpt = np.diag([nkpt] * 3) elif isinstance(nkpt[0], Integral): nkpt = np.diag(nkpt) # Now we have a matrix of k-points if np.any(nkpt - np.diag(np.diag(nkpt)) != 0): raise NotImplementedError(self.__class__.__name__ + " with off-diagonal components is not implemented yet") if displacement is None: displacement = np.zeros(3, np.float64) elif isinstance(displacement, Real): displacement = np.zeros(3, np.float64) + displacement if size is None: size = _a.onesd(3) elif isinstance(size, Real): size = _a.zerosd(3) + size else: size = _a.arrayd(size) # Retrieve the diagonal number of values Dn = np.diag(nkpt).astype(np.int32) if np.any(Dn) == 0: raise ValueError(self.__class__.__name__ + ' *must* be initialized with ' 'diagonal elements different from 0.') i_trs = -1 if trs: # Figure out which direction to TRS nmax = 0 for i in [0, 1, 2]: if displacement[i] in [0., 0.5] and Dn[i] > nmax: nmax = Dn[i] i_trs = i if nmax == 1: i_trs = -1 if i_trs == -1: # If we still haven't decided (say for weird displacements) # simply take the one with the maximum number of k-points. i_trs = np.argmax(Dn) # Calculate k-points and weights along all directions kw = [self.grid(Dn[i], displacement[i], size[i], centered, i == i_trs) for i in (0, 1, 2)] self._k = _a.emptyd((kw[0][0].size, kw[1][0].size, kw[2][0].size, 3)) self._w = _a.onesd(self._k.shape[:-1]) for i in (0, 1, 2): k = kw[i][0].reshape(-1, 1, 1) w = kw[i][1].reshape(-1, 1, 1) self._k[..., i] = np.rollaxis(k, 0, i + 1) self._w[...] *= np.rollaxis(w, 0, i + 1) del kw self._k.shape = (-1, 3) self._k = np.where(self._k > .5, self._k - 1, self._k) self._w.shape = (-1,) # Store information regarding size and diagonal elements # This information is basically only necessary when # we want to replace special k-points self._diag = Dn # vector self._displ = displacement # vector self._size = size # vector self._centered = centered self._trs = i_trs
def param_circle(self, sc, N_or_dk, kR, normal, origo, loop=False): r""" Create a parameterized k-point list where the k-points are generated on a circle around an origo The generated circle is a perfect circle in the reciprocal space (Cartesian coordinates). To generate a perfect circle in units of the reciprocal lattice vectors one can generate the circle for a diagonal supercell with side-length :math:`2\pi`, see example below. Parameters ---------- sc : SuperCell, or SuperCellChild the supercell used to construct the k-points N_or_dk : int number of k-points generated using the parameterization (if an integer), otherwise it specifies the discretization length on the circle (in 1/Ang), If the latter case will use less than 4 points a warning will be raised and the number of points increased to 4. kR : float radius of the k-point. In 1/Ang normal : array_like of float normal vector to determine the circle plane origo : array_like of float origo of the circle used to generate the circular parameterization loop : bool, optional whether the first and last point are equal Examples -------- >>> sc = SuperCell([1, 1, 10, 90, 90, 60]) >>> bz = BrillouinZone.param_circle(sc, 10, 0.05, [0, 0, 1], [1./3, 2./3, 0]) To generate a circular set of k-points in reduced coordinates (reciprocal >>> sc = SuperCell([1, 1, 10, 90, 90, 60]) >>> bz = BrillouinZone.param_circle(sc, 10, 0.05, [0, 0, 1], [1./3, 2./3, 0]) >>> bz_rec = BrillouinZone.param_circle(2*np.pi, 10, 0.05, [0, 0, 1], [1./3, 2./3, 0]) >>> bz.k[:, :] = bz_rec.k[:, :] Returns ------- BrillouinZone : with the parameterized k-points. """ if isinstance(N_or_dk, Integral): N = N_or_dk else: # Calculate the required number of points N = int(kR ** 2 * np.pi / N_or_dk + 0.5) if N < 4: N = 4 info('BrillouinZone.param_circle increased the number of circle points to 4.') # Conversion object bz = BrillouinZone(sc) normal = _a.arrayd(normal) origo = _a.arrayd(origo) k_n = bz.tocartesian(normal) k_o = bz.tocartesian(origo) # Generate a preset list of k-points on the unit-circle if loop: radians = _a.aranged(N) / (N-1) * 2 * np.pi else: radians = _a.aranged(N) / N * 2 * np.pi k = _a.emptyd([N, 3]) k[:, 0] = np.cos(radians) k[:, 1] = np.sin(radians) k[:, 2] = 0. # Now generate the rotation _, theta, phi = cart2spher(k_n) if theta != 0: pv = _a.arrayd([k_n[0], k_n[1], 0]) pv /= fnorm(pv) q = Quaternion(phi, pv, rad=True) * Quaternion(theta, [0, 0, 1], rad=True) else: q = Quaternion(0., [0, 0, k_n[2] / abs(k_n[2])], rad=True) # Calculate k-points k = q.rotate(k) k *= kR / fnorm(k).reshape(-1, 1) k = bz.toreduced(k + k_o) # The sum of weights is equal to the BZ area W = np.pi * kR ** 2 w = np.repeat([W / N], N) return BrillouinZone(sc, k, w)
def read_geometry(self, ret_dynamic=False): r""" Returns Geometry object from the CONTCAR/POSCAR file Possibly also return the dynamics (if present). Parameters ---------- ret_dynamic : bool, optional also return selective dynamics (if present), if not, None will be returned. """ sc = self.read_supercell() # The species labels are not always included in *CAR line1 = self.readline().split() opt = self.readline().split() try: species = line1 species_count = np.array(opt, np.int32) except: species_count = np.array(line1, np.int32) # We have no species... # We default to consecutive elements in the # periodic table. species = [i + 1 for i in range(len(species_count))] err = '\n'.join([ "POSCAR best format:", " <Specie-1> <Specie-2>", " <#Specie-1> <#Specie-2>", "Format not found, the species are defaulted to the first elements of the periodic table." ]) warn(err) # Create list of atoms to be used subsequently atom = [ Atom[spec] for spec, nsp in zip(species, species_count) for i in range(nsp) ] # Number of atoms na = len(atom) # check whether this is Selective Dynamics opt = self.readline() if opt[0] in 'Ss': dynamics = True # pre-create the dynamic list dynamic = np.empty([na, 3], dtype=np.bool_) opt = self.readline() else: dynamics = False dynamic = None # Check whether this is in fractional or direct # coordinates (Direct == fractional) cart = False if opt[0] in 'CcKk': cart = True xyz = _a.emptyd([na, 3]) for ia in range(na): line = self.readline().split() xyz[ia, :] = list(map(float, line[:3])) if dynamics: dynamic[ia] = list(map(lambda x: x.lower() == 't', line[3:6])) if cart: # The unit of the coordinates are cartesian xyz *= self._scale else: xyz = xyz.dot(sc.cell) # The POT/CONT-CAR does not contain information on the atomic species geom = Geometry(xyz=xyz, atom=atom, sc=sc) if ret_dynamic: return geom, dynamic return geom
def read_data(self, as_dataarray=False): """ Returns data associated with the bands file Parameters -------- as_dataarray: boolean, optional if `True`, the information is returned as an `xarray.DataArray` Ticks (if read) are stored as an attribute of the DataArray (under `array.ticks` and `array.ticklabels`) """ band_lines = False # Luckily the data is in eV Ef = float(self.readline()) # Read the total length of the path (not used) _, _ = map(float, self.readline().split()) l = self.readline() try: _, _ = map(float, l.split()) band_lines = True except: # We are dealing with a band-points file pass # orbitals, n-spin, n-k if band_lines: l = self.readline() no, ns, nk = map(int, l.split()) # Create the data to contain all band points b = _a.emptyd([nk, ns, no]) # for band-lines if band_lines: k = _a.emptyd([nk]) for ik in range(nk): l = [float(x) for x in self.readline().split()] k[ik] = l[0] del l[0] # Now populate the eigenvalues while len(l) < ns * no: l.extend([float(x) for x in self.readline().split()]) l = _a.arrayd(l) l.shape = (ns, no) b[ik, :, :] = l[:, :] - Ef # Now we need to read the labels for the points xlabels = [] labels = [] nl = int(self.readline()) for _ in range(nl): l = self.readline().split() xlabels.append(float(l[0])) labels.append((' '.join(l[1:])).replace("'", '')) vals = (xlabels, labels), k, b else: k = _a.emptyd([nk, 3]) for ik in range(nk): l = [float(x) for x in self.readline().split()] k[ik, :] = l[0:3] del l[2] del l[1] del l[0] # Now populate the eigenvalues while len(l) < ns * no: l.extend([float(x) for x in self.readline().split()]) l = _a.arrayd(l) l.shape = (ns, no) b[ik, :, :] = l[:, :] - Ef vals = k, b if as_dataarray: from xarray import DataArray ticks = { "ticks": xlabels, "ticklabels": labels } if band_lines else {} return DataArray(b, name="Energy", attrs=ticks, coords=[("k", k), ("spin", _a.arangei(0, b.shape[1])), ("band", _a.arangei(0, b.shape[2]))]) return vals
def wavefunction(v, grid, geometry=None, k=None, spinor=0, spin=None, eta=False): r""" Add the wave-function (`Orbital.psi`) component of each orbital to the grid This routine calculates the real-space wave-function components in the specified grid. This is an *in-place* operation that *adds* to the current values in the grid. It may be instructive to check that an eigenstate is normalized: >>> grid = Grid(...) # doctest: +SKIP >>> psi(state, grid) # doctest: +SKIP >>> (np.abs(grid.grid) ** 2).sum() * grid.dvolume == 1. # doctest: +SKIP Note: To calculate :math:`\psi(\mathbf r)` in a unit-cell different from the originating geometry, simply pass a grid with a unit-cell smaller than the originating supercell. The wavefunctions are calculated in real-space via: .. math:: \psi(\mathbf r) = \sum_i\phi_i(\mathbf r) |\psi\rangle_i \exp(-i\mathbf k \mathbf R) While for non-collinear/spin-orbit calculations the wavefunctions are determined from the spinor component (`spinor`) .. math:: \psi_{\alpha/\beta}(\mathbf r) = \sum_i\phi_i(\mathbf r) |\psi_{\alpha/\beta}\rangle_i \exp(-i\mathbf k \mathbf R) where ``spinor in [0, 1]`` determines :math:`\alpha` or :math:`\beta`, respectively. Notes ----- Currently this method only works for :math:`\Gamma` states Parameters ---------- v : array_like coefficients for the orbital expansion on the real-space grid. If `v` is a complex array then the `grid` *must* be complex as well. grid : Grid grid on which the wavefunction will be plotted. If multiple eigenstates are in this object, they will be summed. geometry : Geometry, optional geometry where the orbitals are defined. This geometry's orbital count must match the number of elements in `v`. If this is ``None`` the geometry associated with `grid` will be used instead. k : array_like, optional k-point associated with wavefunction, by default the inherent k-point used to calculate the eigenstate will be used (generally shouldn't be used unless the `EigenstateElectron` object has not been created via `Hamiltonian.eigenstate`). spinor : int, optional the spinor for non-collinear/spin-orbit calculations. This is only used if the eigenstate object has been created from a parent object with a `Spin` object contained, *and* if the spin-configuration is non-collinear or spin-orbit coupling. Default to the first spinor component. spin : Spin, optional specification of the spin configuration of the orbital coefficients. This only has influence for non-collinear wavefunctions where `spinor` choice is important. eta : bool, optional Display a console progressbar. """ if geometry is None: geometry = grid.geometry warn( 'wavefunction was not passed a geometry associated, will use the geometry associated with the Grid.' ) if geometry is None: raise SislError( 'wavefunction did not find a usable Geometry through keywords or the Grid!' ) # In case the user has passed several vectors we sum them to plot the summed state if v.ndim == 2: v = v.sum(0) if spin is None: if len(v) // 2 == geometry.no: # We can see from the input that the vector *must* be a non-collinear calculation v = v.reshape(-1, 2)[:, spinor] info( 'wavefunction assumes the input wavefunction coefficients to originate from a non-collinear calculation!' ) elif spin.kind > Spin.POLARIZED: # For non-collinear cases the user selects the spinor component. v = v.reshape(-1, 2)[:, spinor] if len(v) != geometry.no: raise ValueError( "wavefunction require wavefunction coefficients corresponding to number of orbitals in the geometry." ) # Check for k-points k = _a.asarrayd(k) kl = (k**2).sum()**0.5 has_k = kl > 0.000001 if has_k: raise NotImplementedError( 'wavefunction for k != Gamma does not produce correct wavefunctions!' ) # Check that input/grid makes sense. # If the coefficients are complex valued, then the grid *has* to be # complex valued. # Likewise if a k-point has been passed. is_complex = np.iscomplexobj(v) or has_k if is_complex and not np.iscomplexobj(grid.grid): raise SislError( "wavefunction input coefficients are complex, while grid only contains real." ) if is_complex: psi_init = _a.zerosz else: psi_init = _a.zerosd # Extract sub variables used throughout the loop shape = _a.asarrayi(grid.shape) dcell = grid.dcell ic = grid.sc.icell * shape.reshape(1, -1) geom_shape = dot(ic, geometry.cell.T).T # In the following we don't care about division # So 1) save error state, 2) turn off divide by 0, 3) calculate, 4) turn on old error state old_err = np.seterr(divide='ignore', invalid='ignore') addouter = add.outer def idx2spherical(ix, iy, iz, offset, dc, R): """ Calculate the spherical coordinates from indices """ rx = addouter(addouter(ix * dc[0, 0], iy * dc[1, 0]), iz * dc[2, 0] - offset[0]).ravel() ry = addouter(addouter(ix * dc[0, 1], iy * dc[1, 1]), iz * dc[2, 1] - offset[1]).ravel() rz = addouter(addouter(ix * dc[0, 2], iy * dc[1, 2]), iz * dc[2, 2] - offset[2]).ravel() # Total size of the indices n = rx.shape[0] # Reduce our arrays to where the radius is "fine" idx = indices_le(rx**2 + ry**2 + rz**2, R**2) rx = rx[idx] ry = ry[idx] rz = rz[idx] xyz_to_spherical_cos_phi(rx, ry, rz) return n, idx, rx, ry, rz # Figure out the max-min indices with a spacing of 1 radian rad1 = pi / 180 theta, phi = ogrid[-pi:pi:rad1, 0:pi:rad1] cphi, sphi = cos(phi), sin(phi) ctheta_sphi = cos(theta) * sphi stheta_sphi = sin(theta) * sphi del sphi nrxyz = (theta.size, phi.size, 3) del theta, phi, rad1 # First we calculate the min/max indices for all atoms idx_mm = _a.emptyi([geometry.na, 2, 3]) rxyz = _a.emptyd(nrxyz) rxyz[..., 0] = ctheta_sphi rxyz[..., 1] = stheta_sphi rxyz[..., 2] = cphi # Reshape rxyz.shape = (-1, 3) idx = dot(ic, rxyz.T) idxm = idx.min(1).reshape(1, 3) idxM = idx.max(1).reshape(1, 3) del ctheta_sphi, stheta_sphi, cphi, idx, rxyz, nrxyz origo = grid.sc.origo.reshape(1, -1) for atom, ia in geometry.atom.iter(True): if len(ia) == 0: continue R = atom.maxR() # Now do it for all the atoms to get indices of the middle of # the atoms # The coordinates are relative to origo, so we need to shift (when writing a grid # it is with respect to origo) xyz = geometry.xyz[ia, :] - origo idx = dot(ic, xyz.T).T # Get min-max for all atoms idx_mm[ia, 0, :] = idxm * R + idx idx_mm[ia, 1, :] = idxM * R + idx # Now we have min-max for all atoms # When we run the below loop all indices can be retrieved by looking # up in the above table. # Before continuing, we can easily clean up the temporary arrays del origo, idx aranged = _a.aranged # In case this grid does not have a Geometry associated # We can *perhaps* easily attach a geometry with the given # atoms in the unit-cell sc = grid.sc.copy() if grid.geometry is None: # Create the actual geometry that encompass the grid ia, xyz, _ = geometry.within_inf(sc) if len(ia) > 0: grid.set_geometry(Geometry(xyz, geometry.atom[ia], sc=sc)) # Instead of looping all atoms in the supercell we find the exact atoms # and their supercell indices. add_R = _a.zerosd(3) + geometry.maxR() # Calculate the required additional vectors required to increase the fictitious # supercell by add_R in each direction. # For extremely skewed lattices this will be way too much, hence we make # them square. o = sc.toCuboid(True) sc = SuperCell(o._v, origo=o.origo) + np.diag(2 * add_R) sc.origo -= add_R # Retrieve all atoms within the grid supercell # (and the neighbours that connect into the cell) IA, XYZ, ISC = geometry.within_inf(sc) r_k = dot(geometry.rcell, k) r_k_cell = dot(r_k, geometry.cell) phase = 1 # Retrieve progressbar eta = tqdm_eta(len(IA), 'wavefunction', 'atom', eta) # Loop over all atoms in the grid-cell for ia, xyz, isc in zip(IA, XYZ - grid.origo.reshape(1, 3), ISC): # Get current atom atom = geometry.atom[ia] # Extract maximum R R = atom.maxR() if R <= 0.: warn("Atom '{}' does not have a wave-function, skipping atom.". format(atom)) eta.update() continue # Get indices in the supercell grid idx = (isc.reshape(3, 1) * geom_shape).sum(0) idxm = floor(idx_mm[ia, 0, :] + idx).astype(int32) idxM = ceil(idx_mm[ia, 1, :] + idx).astype(int32) + 1 # Fast check whether we can skip this point if idxm[0] >= shape[0] or idxm[1] >= shape[1] or idxm[2] >= shape[2] or \ idxM[0] <= 0 or idxM[1] <= 0 or idxM[2] <= 0: eta.update() continue # Truncate values if idxm[0] < 0: idxm[0] = 0 if idxM[0] > shape[0]: idxM[0] = shape[0] if idxm[1] < 0: idxm[1] = 0 if idxM[1] > shape[1]: idxM[1] = shape[1] if idxm[2] < 0: idxm[2] = 0 if idxM[2] > shape[2]: idxM[2] = shape[2] # Now idxm/M contains min/max indices used # Convert to spherical coordinates n, idx, r, theta, phi = idx2spherical(aranged(idxm[0], idxM[0]), aranged(idxm[1], idxM[1]), aranged(idxm[2], idxM[2]), xyz, dcell, R) # Get initial orbital io = geometry.a2o(ia) if has_k: phase = np.exp(-1j * (dot(r_k_cell, isc))) # TODO # Possibly the phase should be an additional # array for the position in the unit-cell! # + np.exp(-1j * dot(r_k, spher2cart(r, theta, np.arccos(phi)).T) ) # Allocate a temporary array where we add the psi elements psi = psi_init(n) # Loop on orbitals on this atom, grouped by radius for os in atom.iter(True): # Get the radius of orbitals (os) oR = os[0].R if oR <= 0.: warn( "Orbital(s) '{}' does not have a wave-function, skipping orbital!" .format(os)) # Skip these orbitals io += len(os) continue # Downsize to the correct indices if R - oR < 1e-6: idx1 = idx.view() r1 = r.view() theta1 = theta.view() phi1 = phi.view() else: idx1 = indices_le(r, oR) # Reduce arrays r1 = r[idx1] theta1 = theta[idx1] phi1 = phi[idx1] idx1 = idx[idx1] # Loop orbitals with the same radius for o in os: # Evaluate psi component of the wavefunction and add it for this atom psi[idx1] += o.psi_spher(r1, theta1, phi1, cos_phi=True) * (v[io] * phase) io += 1 # Clean-up del idx1, r1, theta1, phi1, idx, r, theta, phi # Convert to correct shape and add the current atom contribution to the wavefunction psi.shape = idxM - idxm grid.grid[idxm[0]:idxM[0], idxm[1]:idxM[1], idxm[2]:idxM[2]] += psi # Clean-up del psi # Step progressbar eta.update() eta.close() # Reset the error code for division np.seterr(**old_err)