def Eindex(self, E): """ Return the closest energy index corresponding to the energy ``E`` Parameters ---------- E : float or int if ``int``, return it-self, else return the energy index which is closests to the energy. """ if isinstance(E, Integral): return E elif isinstance(E, _str): # This will always be converted to a float E = float(E) idxE = np.abs(self.E - E).argmin() ret_E = self.E[idxE] if abs(ret_E - E) > 5e-3: warn(self.__class__.__name__ + " requesting energy " + "{0:.5f} eV, found {1:.5f} eV as the closest energy!".format( E, ret_E)) elif abs(ret_E - E) > 1e-3: info(self.__class__.__name__ + " requesting energy " + "{0:.5f} eV, found {1:.5f} eV as the closest energy!".format( E, ret_E)) return idxE
def undo_settings_group(self, group): """ Takes the desired group of settings one step back, but the rest of the settings remain unchanged At the moment it is a 'fake' undo function, since it actually updates the settings. Parameters ----------- group: str the key of the settings group for which you want to undo its values. """ #Get the actual settings for that group actualSettings = self.get_settings_group(group) #Try to find any different values for the settings for i in range(len(self.settings_history)): previousSettings = self.get_settings_group(group, steps_back=i) if previousSettings != actualSettings: return self.update_settings(previousSettings) else: info(f"group={group} was never changed; cannot undo nothing.") return self
def read_hamiltonian(self, **kwargs): """ Returns a Hamiltonian from the underlying NetCDF file """ H = self._read_class_spin(Hamiltonian, **kwargs) S = H._csr._D[:, H.S_idx] Ef = self._value('Ef')[:] * Ry2eV if Ef.size == 1: Ef = np.tile(Ef, 2) else: dEf = np.diff(Ef)[0] info( repr(self) + '.read_hamiltonian found a calculation with spin-dependent Fermi-levels: ' 'dEf={:.4f} eV. ' 'Both spin configurations are shifted to 0. This may change in future ' 'versions of sisl.'.format(dEf)) sp = self._crt_grp(self, 'SPARSE') for i in range(len(H.spin)): # Create new container h = np.array(sp.variables['H'][i, :], np.float64) * Ry2eV # Correct for the Fermi-level, Ef == 0 if i < 2: h -= Ef[i] * S[:] H._csr._D[:, i] = h[:] return H
def _r_basis_orb_indx(self): f = self._tofile(self.get('SystemLabel', default='siesta')) + '.ORB_INDX' if isfile(f): info( SileInfo( 'Siesta basis information is read from {}, the radial functions are in accessible.' .format(f))) return orbindxSileSiesta(f).read_basis() return None
def undo_setting(self, key): """ Undoes only a particular setting and leaves the others unchanged At the moment it is a 'fake' undo function, since it actually updates the settings. Parameters ----------- key: str the key of the setting that you want to undo. """ i = self.settings_history.last_update_for(key) if i is None: info(f"key={key} was never changed; cannot undo nothing.") self.update_settings(key=self.settings_history[key][i]) return self
def trigger_notification(title, message, sound="Submarine"): """ Triggers a notification. Will not do anything in Windows (oops!) Parameters ----------- title: str message: str sound: str """ if sys.platform == 'linux': os.system(f"""notify-send "{title}" "{message}" """) elif sys.platform == 'darwin': sound_string = f'sound name "{sound}"' if sound else '' os.system(f"""osascript -e 'display notification "{message}" with title "{title}" {sound_string}' """) else: info(f"sisl cannot issue notifications through the operating system ({sys.platform})")
def Eindex(self, E): """ Return the closest energy index corresponding to the energy ``E`` Parameters ---------- E : float or int if ``int``, return it-self, else return the energy index which is closests to the energy. """ if isinstance(E, Integral): return E idxE = np.abs(self.E - E).argmin() ret_E = self.E[idxE] if abs(ret_E - E) > 5e-3: warn(self.__class__.__name__ + " requesting energy " + f"{E:.5f} eV, found {ret_E:.5f} eV as the closest energy!") elif abs(ret_E - E) > 1e-3: info(self.__class__.__name__ + " requesting energy " + f"{E:.5f} eV, found {ret_E:.5f} eV as the closest energy!") return idxE
def spoken_message(message): """ Trigger a spoken message. In linux espeak must be installed (sudo apt-get install espeak) Will not do anything in Windows (oops!) Parameters ----------- title: str message: str sound: str """ if sys.platform == 'linux': os.system(f"""espeak -s 150 "{message}" 2>/dev/null""") elif sys.platform == 'darwin': os.system(f"""osascript -e 'say "{message}"' """) else: info(f"sisl cannot issue notifications through the operating system ({sys.platform})")
def undo_settings(self, steps=1, run_updates=True): """ Brings the settings back a number of steps Parameters ------------ steps: int, optional the number of steps you want to go back. run_updates: bool, optional whether we should run updates after updating the settings. If not, the settings will be updated, but you won't see any change in the object. """ try: diff = self.settings_history.diff_keys(-1, -steps - 1) self.settings_history.undo(steps=steps) if run_updates: self._run_updates(diff) except IndexError: info( f"This instance of {self.__class__.__name__} does not " f"contain earlier settings as requested ({steps} step(s) back)" ) return self
def Eindex(self, E): """ Return the closest energy index corresponding to the energy ``E`` Parameters ---------- E : float or int or str if `int`, return it-self, else return the energy index which is closests to the energy. For a `str` it will be parsed to a float and treated as such. """ if isinstance(E, Integral): return E elif isinstance(E, str): E = float(E) idxE = np.abs(self.E - E).argmin() ret_E = self.E[idxE] if abs(ret_E - E) > 5e-3: warn(f"{self.__class__.__name__} requesting energy " f"{E:.5f} eV, found {ret_E:.5f} eV as the closest energy!") elif abs(ret_E - E) > 1e-3: info(f"{self.__class__.__name__} requesting energy " f"{E:.5f} eV, found {ret_E:.5f} eV as the closest energy!") return idxE
def _r_dynamical_matrix_got(self, geometry, **kwargs): """ In case the dynamical matrix is read from the file """ # Easier for creation of the sparsity pattern from scipy.sparse import lil_matrix # Default cutoff eV / Ang ** 2 cutoff = kwargs.get('cutoff', 1.e-4) dtype = kwargs.get('dtype', np.float64) nxyz = geometry.no dyn = lil_matrix((nxyz, nxyz), dtype=dtype) f, _ = self.step_to(self._keys['dyn']) if not f: info( self.__class__.__name__ + ' tries to lookup the Dynamical matrix ' 'using key "' + self._keys['dyn'] + '". ' 'Use .set_dynamical_matrix_key(...) to search for different name.' 'This could not be found found in file: "{}".'.format( self.file)) return None # skip 1 line self.readline() # default range dat = np.empty([nxyz], dtype=dtype) i, j = 0, 0 nxyzm1 = nxyz - 1 while i < nxyz: l = self.readline().strip() if len(l) == 0: break # convert to float list ls = [float(x) for x in l.split()] k = min(12, nxyz - j) # GULP only prints columns corresponding # to a full row. Hence the remaining # data must be nxyz - j - 1 dat[j:j + k] = ls[:k] j += k if j >= nxyz: dyn[i, :] = dat[:] # step row i += 1 # reset column j = 0 # clean-up for memory del dat # Convert to COO matrix format dyn = dyn.tocoo() # Construct mass ** (-.5), so we can check cutoff correctly mass_sqrt = geometry.atoms.mass.repeat(3)**0.5 dyn.data[:] *= mass_sqrt[dyn.row] * mass_sqrt[dyn.col] dyn.data[np_abs(dyn.data) < cutoff] = 0. dyn.data[:] *= 1 / (mass_sqrt[dyn.row] * mass_sqrt[dyn.col]) dyn.eliminate_zeros() return dyn
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 test_info_category(): with pytest.warns(sm.SislInfo): sm.info('Information', sm.SislInfo)
def test_info_specific(): with pytest.warns(sm.SislInfo): sm.info(sm.SislInfo('Info'))
def test_info_method(): with pytest.warns(sm.SislInfo): sm.info('Information')
def initialize(self): """ Initialize the internal data-arrays used for efficient calculation of the real-space quantities This method should first be called *after* all options has been specified. If the user hasn't specified the ``bz`` value as an option this method will update the internal integration Brillouin zone based on the ``dk`` option. """ # Try and guess the directions unfold = self._unfold nsc = self.parent.nsc.copy() axes = self._options['axes'] if axes is None: if nsc[2] == 1: axes = _a.arrayi([0, 1]) elif nsc[1] == 1: axes = _a.arrayi([0, 2]) elif nsc[0] == 1: axes = _a.arrayi([1, 2]) else: raise ValueError(self.__class__.__name__ + '.initialize currently only supports a 2D real-space self-energy, hence the MonkhorstPack grid should reflect only 2 periodic directions.') self._options['axes'] = axes # Check that we have periodicity along the chosen axes nsc_sum = nsc[axes].sum() if nsc_sum == 1: raise ValueError(self.__class__.__name__ + '.initialize found no periodic directions for the chosen integration axes: {} and {}.'.format(*nsc[axes])) elif nsc_sum < 6: raise ValueError((self.__class__.__name__ + '.initialize found one periodic direction ' 'out of two for the chosen integration axes: {} and {}. ' 'For 1D systems the regular surface self-energy method is appropriate.').format(*nsc[axes])) if self._options['semi_axis'] is None and self._options['k_axis'] is None: # None of the axis has been described if nsc[axes[0]] > 3: k_ax = axes[0] s_ax = axes[1] elif nsc[axes[1]] > 3: k_ax = axes[1] s_ax = axes[0] else: # Choose the direction of k to be the smallest shortest sc = self.parent.sc.tile(unfold[axes[0]], axes[0]).tile(unfold[axes[1]], axes[1]) rcell = fnorm(sc.rcell)[axes] k_ax = axes[np.argmax(rcell)] if k_ax == axes[0]: s_ax = axes[1] else: s_ax = axes[0] self._options['semi_axis'] = s_ax self._options['k_axis'] = k_ax s_ax = 'ABC'[s_ax] k_ax = 'ABC'[k_ax] info(self.__class__.__name__ + '.initialize determined the semi-inf- and k-directions to be: {} and {}'.format(s_ax, k_ax)) elif self._options['k_axis'] is None: s_ax = self._options['semi_axis'] if s_ax == axes[0]: k_ax = axes[1] else: k_ax = axes[0] self._options['k_axis'] = k_ax k_ax = 'ABC'[k_ax] info(self.__class__.__name__ + '.initialize determined the k direction to be: {}'.format(k_ax)) elif self._options['semi_axis'] is None: k_ax = self._options['k_axis'] if k_ax == axes[0]: s_ax = axes[1] else: s_ax = axes[0] self._options['semi_axis'] = s_ax s_ax = 'ABC'[s_ax] info(self.__class__.__name__ + '.initialize determined the semi-infinite direction to be: {}'.format(s_ax)) k_ax = self._options['k_axis'] s_ax = self._options['semi_axis'] if nsc[s_ax] != 3: raise ValueError(self.__class__.__name__ + '.initialize found the self-energy direction to be ' 'incompatible with the parent object. It *must* have 3 supercells along the ' 'semi-infinite direction.') P0 = self.real_space_parent() V_atoms = self.real_space_coupling(True)[1] self._calc = { # The below algorithm requires the direction to be negative # if changed, B, C should be reversed below 'SE': RecursiveSI(self.parent, '-' + 'ABC'[s_ax], eta=self._options['eta']), # Used to calculate the real-space self-energy 'P0': P0.Pk(), # in sparse format 'S0': P0.Sk(), # in sparse format # Orbitals in the coupling atoms 'orbs': P0.a2o(V_atoms, True).reshape(-1, 1), } # Update the BrillouinZone integration grid in case it isn't specified if self._options['bz'] is None: # Update the integration grid # Note this integration grid is based on the big system. sc = self.parent.sc.tile(unfold[axes[0]], axes[0]).tile(unfold[axes[1]], axes[1]) rcell = fnorm(sc.rcell)[k_ax] nk = [1] * 3 nk[k_ax] = int(self._options['dk'] * rcell) self._options['bz'] = MonkhorstPack(sc, nk, trs=self._options['trs']) info(self.__class__.__name__ + '.initialize determined the number of k-points: {}'.format(nk[k_ax]))
def _r_dynamical_matrix_got(self, geom, **kwargs): """ In case the dynamical matrix is read from the file """ # Easier for creation of the sparsity pattern from scipy.sparse import lil_matrix # Default cutoff eV / Ang ** 2 cutoff = kwargs.get('cutoff', 1.e-4) dtype = kwargs.get('dtype', np.float64) no = geom.no dyn = lil_matrix((no, no), dtype=dtype) f, _ = self.step_to(self._keys['dyn']) if not f: info( self.__class__.__name__ + ' tries to lookup the Dynamical matrix ' 'using key "' + self._keys['dyn'] + '". ' 'Use .set_dynamical_matrix_key(...) to search for different name.' 'This could not be found found in file: "{}".'.format( self.file)) return None # skip 1 line self.readline() # default range dat = np.empty([no], dtype=dtype) i, j = 0, 0 while True: l = self.readline().strip() if len(l) == 0: break # convert to float list ls = [float(x) for x in l.split()] if j + 12 <= no: # Here the full line can fit for the same row dat[j:j + 12] = ls[:12] j += 12 if j >= no: dyn[i, :] = dat[:] # step row i += 1 # reset column j = 0 if i >= no: break else: # add the values (12 values == 3*4) # for atoms on each line for k in [0, 1, 2, 3]: dat[j:j + 3] = ls[k * 3:(k + 1) * 3] j += 3 if j >= no: # Clear those below the cutoff dyn[i, :] = dat[:] i += 1 j = 0 if i >= no: break # clean-up for memory del dat # Convert to COO matrix format dyn = dyn.tocoo() # Construct mass ** (-.5), so we can check cutoff correctly mass_sqrt = np.array(geom.atoms.mass, np.float64).repeat(3)**0.5 dyn.data[:] *= mass_sqrt[dyn.row] * mass_sqrt[dyn.col] dyn.data[np_abs(dyn.data) < cutoff] = 0. dyn.data[:] *= 1 / (mass_sqrt[dyn.row] * mass_sqrt[dyn.col]) dyn.eliminate_zeros() return dyn
def replace(self, k, mp): r""" Replace a k-point with a new set of k-points from a Monkhorst-Pack grid This method tries to replace an area corresponding to `mp.size` around the k-point `k` such that the k-points are replaced. This enables one to zoom in on specific points in the Brillouin zone for detailed analysis. Parameters ---------- k : array_like k-point in this object to replace mp : MonkhorstPack object containing the replacement k-points. Examples -------- This example creates a zoomed-in view of the :math:`\Gamma`-point by replacing it with a 3x3x3 Monkhorst-Pack grid. >>> sc = SuperCell(1.) >>> mp = MonkhorstPack(sc, [3, 3, 3]) >>> mp.replace([0, 0, 0], MonkhorstPack(sc, [3, 3, 3], size=1./3)) This example creates a zoomed-in view of the :math:`\Gamma`-point by replacing it with a 4x4x4 Monkhorst-Pack grid. >>> sc = SuperCell(1.) >>> mp = MonkhorstPack(sc, [3, 3, 3]) >>> mp.replace([0, 0, 0], MonkhorstPack(sc, [4, 4, 4], size=1./3)) This example creates a zoomed-in view of the :math:`\Gamma`-point by replacing it with a 4x4x1 Monkhorst-Pack grid. >>> sc = SuperCell(1.) >>> mp = MonkhorstPack(sc, [3, 3, 3]) >>> mp.replace([0, 0, 0], MonkhorstPack(sc, [4, 4, 1], size=1./3)) Raises ------ SislError : if the size of the replacement `MonkhorstPack` grid is not compatible with the k-point spacing in this object. """ # First we find all k-points within k +- mp.size # Those are the points we wish to remove. # Secondly we need to ensure that the k-points we remove are occupying *exactly* # the Brillouin zone we wish to replace. if not isinstance(mp, MonkhorstPack): raise ValueError('Object `mp` is not a MonkhorstPack object') # We can easily figure out the BZ that each k-point is averaging k_vol = self._size / self._diag # Compare against the size of this one # Since we can remove more than one k-point, we require that the # size of the replacement MP is an integer multiple of the # k-point volumes. k_int = mp._size / k_vol if not np.allclose(np.rint(k_int), k_int): raise SislError(self.__class__.__name__ + '.reduce could not replace k-point, BZ ' 'volume replaced is not equivalent to the inherent k-point volume.') k_int = np.rint(k_int) # 1. find all k-points k = self.in_primitive(k).reshape(1, 3) dk = (mp._size / 2).reshape(1, 3) # Find all points within [k - dk; k + dk] # Since the volume of each k-point is non-zero we know that no k-points will be located # on the boundary. # This does remove boundary points because we shift everything into the positive # plane. diff_k = self.in_primitive(self.k % 1. - k % 1.) idx = np.logical_and.reduce(np.abs(diff_k) <= dk, axis=1).nonzero()[0] if len(idx) == 0: raise SislError(self.__class__.__name__ + '.reduce could not find any points to replace.') # Now we have the k-points we need to remove # Figure out if the total weight is consistent total_weight = self.weight[idx].sum() replace_weight = mp.weight.sum() if abs(total_weight - replace_weight) < 1e-8: weight_factor = 1. elif abs(total_weight - replace_weight * 2) < 1e-8: weight_factor = 2. if self._trs < 0: info(self.__class__.__name__ + '.reduce assumes that the replaced k-point has double weights.') else: print('k-point to replace:') print(' ', k.ravel()) print('delta-k:') print(' ', dk.ravel()) print('Found k-indices that will be replaced:') print(' ', idx) print('k-points replaced:') print(self.k[idx, :]) raise SislError(self.__class__.__name__ + '.reduce could not assert the weights are consistent during replacement.') self._k = np.delete(self._k, idx, axis=0) self._w = np.delete(self._w, idx) # Append the new k-points and weights self._k = np.concatenate((self._k, mp._k), axis=0) self._w = np.concatenate((self._w, mp._w * weight_factor))
def read_grid(self, spin=0, name='gridfunc', *args, **kwargs): """ Reads a grid in the current Siesta.grid.nc file Enables the reading and processing of the grids created by Siesta Parameters ---------- spin : int or array_like, optional specify the retrieved values name : str, optional the name for the grid-function (do not supply for standard Siesta output) geometry: Geometry, optional add the Geometry to the Grid """ # Default to *index* variable spin = kwargs.get('index', spin) # Determine the name of this file f = osp.basename(self.file) # File names are made up of # ElectrostaticPotential.grid.nc # So the first one should be ElectrostaticPotential try: # <>.grid.nc base = f.split('.')[-3] except: base = 'None' # Unit-conversion BohrC2AngC = Bohr2Ang**3 unit = { 'Rho': 1. / BohrC2AngC, 'DeltaRho': 1. / BohrC2AngC, 'RhoXC': 1. / BohrC2AngC, 'RhoInit': 1. / BohrC2AngC, 'Chlocal': 1. / BohrC2AngC, 'TotalCharge': 1. / BohrC2AngC, 'BaderCharge': 1. / BohrC2AngC, 'ElectrostaticPotential': Ry2eV, 'TotalPotential': Ry2eV, 'Vna': Ry2eV, }.get(base, None) # Fall-back if unit is None: unit = 1. show_info = True else: show_info = False # Swap as we swap back in the end sc = self.read_supercell().swapaxes(0, 2) # Create the grid nx = len(self._dimension('n1')) ny = len(self._dimension('n2')) nz = len(self._dimension('n3')) if name is None: v = self._variable('gridfunc') else: v = self._variable(name) # Create the grid, Siesta uses periodic, always grid = Grid([nz, ny, nx], bc=Grid.PERIODIC, sc=sc, dtype=v.dtype, geometry=kwargs.get("geometry", None)) if v.ndim == 3: grid.grid[:, :, :] = v[:, :, :] * unit elif isinstance(spin, Integral): grid.grid[:, :, :] = v[spin, :, :, :] * unit else: if len(spin) > v.shape[0]: raise SileError( f"{self.__class__.__name__}.read_grid requires spin to be an integer or " "an array of length equal to the number of spin components." ) grid.grid[:, :, :] = v[0, :, :, :] * spin[0] * unit for i, scale in enumerate(spin[1:]): grid.grid[:, :, :] += v[1 + i, :, :, :] * scale * unit if show_info: info( f"{self.__class__.__name__}.read_grid cannot determine the units of the grid. " "The units may not be in sisl units.") # Read the grid, we want the z-axis to be the fastest # looping direction, hence x,y,z == 0,1,2 return grid.swapaxes(0, 2)
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)