EigenstateElectron

class sisl.physics.electron.EigenstateElectron(state, c, parent=None, **info)[source]

Eigen states of electrons with eigenvectors and eigenvalues.

This holds routines that enable the calculation of (projected) density of states, spin moments (spin texture).

Attributes

c

dkind

The data-type of the state (in str)

dtype

Data-type for the state

eig

Eigenvalues for each state

info

parent

shape

Returns the shape of the state

state

Methods

DOS(self, E[, distribution])

Calculate DOS for provided energies, E.

PDOS(self, E[, distribution])

Calculate PDOS for provided energies, E.

Sk(self[, format, spin])

Retrieve the overlap matrix corresponding to the originating parent structure.

__init__(self, state, c[, parent])

Define a state container with a given set of states and coefficients for the states

align(self, other[, copy])

Align other.state with the angles for this state, a copy of other is returned with rotated elements

asCoefficient(self)

asState(self)

change_gauge(self, gauge)

In-place change of the gauge of the state coefficients

copy(self)

Return a copy (only the coefficients and states are copied), parent and info are passed by reference

degenerate(self, eps)

Find degenerate coefficients with a specified precision

expectation(self, A[, diag])

Calculate the expectation value of matrix A

inner(self[, right, diagonal, align])

Return the inner product by \(\mathbf M_{ij} = \langle\psi_i|\psi'_j\rangle\)

inv_eff_mass_tensor(self[, as_matrix, eps])

Calculate inverse effective mass tensor for the states

iter(self[, asarray])

An iterator looping over the states in this system

norm(self)

Return a vector with the norm of each state \(\sqrt{\langle\psi|\psi\rangle}\)

norm2(self[, sum])

Return a vector with the norm of each state \(\langle\psi|\psi\rangle\)

normalize(self)

Return a normalized state where each state has \(|\psi|^2=1\)

occupation(self[, distribution])

Calculate the occupations for the states according to a distribution function

outer(self[, idx])

Return the outer product for the indices idx (or all if None) by \(\sum_i|\psi_i\rangle c_i\langle\psi_i|\)

phase(self[, method, return_indices])

Calculate the Euler angle (phase) for the elements of the state, in the range \(]-\pi;\pi]\)

rotate(self[, phi, individual])

Rotate all states (in-place) to rotate the largest component to be along the angle phi

sort(self[, ascending])

Sort and return a new StateC by sorting the coefficients (default to ascending)

spin_moment(self)

Calculate spin moment from the states

sub(self, idx)

Return a new state with only the specified states

velocity(self[, eps])

Calculate velocity for the states

velocity_matrix(self[, eps])

Calculate velocity matrix for the states

wavefunction(self, grid[, spinor, eta])

Expand the coefficients as the wavefunction on grid as-is

DOS(self, E, distribution='gaussian')[source]

Calculate DOS for provided energies, E.

This routine calls sisl.physics.electron.DOS with appropriate arguments and returns the DOS.

See DOS for argument details.

PDOS(self, E, distribution='gaussian')[source]

Calculate PDOS for provided energies, E.

This routine calls PDOS with appropriate arguments and returns the PDOS.

See PDOS for argument details.

Sk(self, format='csr', spin=None)

Retrieve the overlap matrix corresponding to the originating parent structure.

When self.parent is a Hamiltonian this will return \(\mathbf S(k)\) for the \(k\)-point these eigenstates originate from

Parameters
formatstr, optional

the returned format of the overlap matrix. This only takes effect for non-orthogonal parents.

spinSpin, optional

for non-collinear spin configurations the fake overlap matrix returned will have halve the size of the input matrix. If you want the full overlap matrix, simply do not specify the spin argument.

align(self, other, copy=False)

Align other.state with the angles for this state, a copy of other is returned with rotated elements

States will be rotated by \(\pi\) provided the phase difference between the states are above \(|\Delta\theta| > \pi/2\).

Parameters
otherState

the other state to align onto this state

copybool, optional

sometimes no states require rotation, if this is the case this flag determines whether other will be copied or not

asCoefficient(self)
asState(self)
c
change_gauge(self, gauge)

In-place change of the gauge of the state coefficients

The two gauges are related through:

\[\tilde C_j = e^{i\mathbf k\mathbf r_j} C_j\]

where \(C_j\) belongs to the gauge R and \(\tilde C_j\) is in the gauge r.

Parameters
gauge{‘R’, ‘r’}

specify the new gauge for the state coefficients

copy(self)

Return a copy (only the coefficients and states are copied), parent and info are passed by reference

degenerate(self, eps)

Find degenerate coefficients with a specified precision

Parameters
epsfloat

the precision above which coefficients are not considered degenerate

Returns
list of numpy.ndarraya list of indices
property dkind

The data-type of the state (in str)

property dtype

Data-type for the state

property eig

Eigenvalues for each state

expectation(self, A, diag=True)

Calculate the expectation value of matrix A

The expectation matrix is calculated as:

\[A_{ij} = \langle \psi_i | \mathbf A | \psi_j \rangle\]

If diag is true, only the diagonal elements are returned.

Parameters
Aarray_like

a vector or matrix that expresses the operator A

diagbool, optional

whether only the diagonal elements are calculated or if the full expectation matrix is calculated

Returns
numpy.ndarray

a vector if diag is true, otherwise a matrix with expectation values

info
inner(self, right=None, diagonal=True, align=True)

Return the inner product by \(\mathbf M_{ij} = \langle\psi_i|\psi'_j\rangle\)

Parameters
rightState, optional

the right object to calculate the inner product with, if not passed it will do the inner product with itself. This object will always be the left \(\langle\psi_i|\)

diagonalbool, optional

only return the diagonal matrix \(\mathbf M_{ii}\).

alignbool, optional

first align right with the angles for this state (see align)

Returns
numpy.ndarray

a matrix with the sum of outer state products

Notes

This does not take into account a possible overlap matrix when non-orthogonal basis sets are used.

inv_eff_mass_tensor(self, as_matrix=False, eps=0.001)

Calculate inverse effective mass tensor for the states

This routine calls inv_eff_mass with appropriate arguments and returns the state inverse effective mass tensor. I.e. for non-orthogonal basis the overlap matrix and energy values are also passed.

Note that the coefficients associated with the StateCElectron must correspond to the energies of the states.

See inv_eff_mass_tensor for details.

Parameters
as_matrixbool, optional

if true the returned tensor will be a symmetric matrix, otherwise the Voigt tensor is returned.

epsfloat, optional

precision used to find degenerate states.

Notes

The reason for not inverting the mass-tensor is that for systems with limited periodicities some of the diagonal elements of the inverse mass tensor matrix will be 0, in which case the matrix is singular and non-invertible. Therefore it is the users responsibility to remove any of the non-periodic elements from the matrix.

iter(self, asarray=False)

An iterator looping over the states in this system

Parameters
asarraybool, optional

if true the yielded values are the state vectors, i.e. a numpy array. Otherwise an equivalent object is yielded.

Yields
stateState

a state only containing individual elements, if asarray is false

statenumpy.ndarray

a state only containing individual elements, if asarray is true

norm(self)

Return a vector with the norm of each state \(\sqrt{\langle\psi|\psi\rangle}\)

Returns
numpy.ndarray

the normalization for each state

norm2(self, sum=True)

Return a vector with the norm of each state \(\langle\psi|\psi\rangle\)

Parameters
sumbool, optional

if true the summed orbital square is returned (a vector). For false a matrix with normalization squared per orbital is returned.

Returns
numpy.ndarray

the normalization on each orbital for each state

normalize(self)

Return a normalized state where each state has \(|\psi|^2=1\)

This is roughly equivalent to:

>>> state = StateC(np.arange(10), 1)
>>> n = state.norm()
>>> norm_state = StateC(state.state / n.reshape(-1, 1), state.c.copy())
>>> norm_state.c[0] == 1
Returns
statea new state with all states normalized, otherwise equal to this
occupation(self, distribution='fermi_dirac')[source]

Calculate the occupations for the states according to a distribution function

Parameters
distributionstr or func, optional

distribution used to find occupations

Returns
numpy.ndarray

len(self) with occupation values

outer(self, idx=None)

Return the outer product for the indices idx (or all if None) by \(\sum_i|\psi_i\rangle c_i\langle\psi_i|\)

Parameters
idxint or array_like, optional

only perform an outer product of the specified indices, otherwise all states are used

Returns
numpy.ndarraya matrix with the sum of outer state products
parent
phase(self, method='max', return_indices=False)

Calculate the Euler angle (phase) for the elements of the state, in the range \(]-\pi;\pi]\)

Parameters
method{‘max’, ‘all’}

for max, the phase for the element which has the largest absolute magnitude is returned, for all, all phases are calculated

return_indicesbool, optional

return indices for the elements used when method=='max'

rotate(self, phi=0.0, individual=False)

Rotate all states (in-place) to rotate the largest component to be along the angle phi

The states will be rotated according to:

\[S' = S / S^\dagger_{\phi-\mathrm{max}} \exp (i \phi),\]

where \(S^\dagger_{\phi-\mathrm{max}}\) is the phase of the component with the largest amplitude and \(\phi\) is the angle to align on.

Parameters
phifloat, optional

angle to align the state at (in radians), 0 is the positive real axis

individualbool, optional

whether the rotation is per state, or a single maximum component is chosen.

property shape

Returns the shape of the state

sort(self, ascending=True)

Sort and return a new StateC by sorting the coefficients (default to ascending)

Parameters
ascendingbool, optional

sort the contained elements ascending, else they will be sorted descending

spin_moment(self)

Calculate spin moment from the states

This routine calls spin_moment with appropriate arguments and returns the spin moment for the states.

See spin_moment for details.

state
sub(self, idx)

Return a new state with only the specified states

Parameters
idxint or array_like

indices that are retained in the returned object

Returns
StateCa new object with a subset of the states
velocity(self, eps=0.0001)

Calculate velocity for the states

This routine calls velocity with appropriate arguments and returns the velocity for the states. I.e. for non-orthogonal basis the overlap matrix and energy values are also passed.

Note that the coefficients associated with the StateCElectron must correspond to the energies of the states.

See velocity for details.

Parameters
epsfloat, optional

precision used to find degenerate states.

velocity_matrix(self, eps=0.0001)

Calculate velocity matrix for the states

This routine calls velocity_matrix with appropriate arguments and returns the velocity for the states. I.e. for non-orthogonal basis the overlap matrix and energy values are also passed.

Note that the coefficients associated with the StateCElectron must correspond to the energies of the states.

See velocity_matrix for details.

Parameters
epsfloat, optional

precision used to find degenerate states.

wavefunction(self, grid, spinor=0, eta=False)

Expand the coefficients as the wavefunction on grid as-is

See wavefunction for argument details, the arguments not present in this method are automatically passed from this object.