Eigenvalues

class Eigenvalues(configuration, kpoints=None, energy_zero_parameter=None, bands_above_fermi_level=None, enable_symmetry=None, processes_per_kpoint=None, method=None, interpolation_method=None, diagonalization_method=None)

Class for calculating the eigenvalues for a configuration.

Parameters:
  • configuration (BulkConfiguration | MoleculeConfiguration) – The configuration with an attached calculator for which to calculate the eigenvalues.

  • kpoints (sequence (size 3) of int | MonkhorstPackGrid | KpointDensity | RegularKpointGrid | sequence of sequence (size 3) of float) – The k-points for which to perform the calculation. This can either be given as a k-point grid or a list of fractional k-points, e.g., [[0.0, 0.0, 0.0], [0.0, 0.0, 0.1], ...].
    Default: The Monkhorst-Pack grid used for the self-consistent calculation.

  • energy_zero_parameter (FermiLevel | AbsoluteEnergy) – Specifies the choice for the energy zero.
    Default: FermiLevel

  • bands_above_fermi_level
    Deprecated: from 2023.09, see diagonalization_method and FullDiagonalizationSolver instead.

  • enable_symmetry (bool) – Enable or disable the use of symmetry
    Default: True

  • processes_per_kpoint (int | None) – The number of processes to use per kpoint. The parameter is only used by the PlaneWaveCalculator, or if an iterative diagonalzation solver is used (see parameter diagonalization_method).
    Default: For PlaneWaveCalculator, the default is the same number of processes per k-point as used in the SCF calculation. Otherwise, the default is 1.

  • method
    Deprecated: from v2023.09, use interpolation_method parameter instead.

  • interpolation_method (None | Full | KDotPExpansion1D | KDotPExpansion3D) –

    The method used for the bandstructure calculation.

    Default is to perform an exact eigenvalue solution for each k-point.

    Alternatively, the k.p expansion method can be use to interpolate the eigenvalues at the requested k-points using the exact eigensolution computed for a small subset of the k-points. There are two modes for the k.p method. 1) Using KDotPExpansion1D, for each bandstructure path segment, a few exact eigensolution are calculated and the eigenvalues of the remaining k-points are interpolated using a 1D correction term. 2) Using KDotPExpansion3D, the eigenvalues for the requested k-points are interpolated from exact solution on a 3D grid. The bands, grid and wave functions of a previous ground state computation can be used in this approach.


    Default: Full meaning that all k-points will be calculated exactly.

  • diagonalization_method

    Method used for diagonalizing the hamiltonian.

    This parameter allows to choose between a full diagonalization solver and an iterative subspace solver. The full diagonalization solver evaluates all bands from the lowest energy one to a given number of bands above fermi level. The iterative subspace solver allows to evaluate a given number of bands around fermi level, or around an energy of choice.

    The full diagonalization solver is more robust, but can be proibitively expensive for very large systems.

    The iterative solver can deal with very large systems (tens of thousands atoms and beyond) and greatly outperforms when calculating a small number of eigenvalues, but it is also inherently less robust.

    Note: the exact method used when selecting FullDiagonalizationSolver is defined by the calculator (see AlgorithmParameters). IterativeDiagonalizationSolver is not supported for PlaneWaveCalculator


    Default: FullDiagonalizationSolver

Type:

FullDiagonalizationSolver | IterativeDiagonalizationSolver

densityOfStates()
Returns:

A DensityOfStates object constructed from the calculated eigenvalues. This is only possible, if the kpoint used in this Eigenvalues object are of the type MonkhorstPackGrid or RegularKpointGrid.

Return type:

DensityOfStates

eigenvalues(spin=None)

Returns the calculated eigenvalues of the Hamiltonian.

Parameters:

spin (Spin.All | Spin.Up | Spin.Down) – The spin of the eigenvalues. For Noncollinear calculations, only Spin.All is valid.
Default: Spin.All

Returns:

The eigenvalues. If spin is Spin.All the shape of the array is (n-spins, n-kpoints, n-bands) for Unpolarized and Polarized calculations, and (n-kpoints, n-bands) for Noncollinear calculations. When spin is Spin.Up or Spin.Down the shape is (n-kpoints, n-bands).

Return type:

PhysicalQuantity of type energy.

energyZero()
Returns:

The energy zero.

Return type:

PhysicalQuantity of type energy.

evaluate(spin=None)

Returns the irreducible kpoints and the corresponding eigenvalues of the Hamiltonian.

Parameters:

spin (Spin.All | Spin.Up | Spin.Down) – The spin of the eigenvalues. For Noncollinear caclulations, only Spin.All is valid.
Default: Spin.All

Returns:

The irriducible kpoints and the corresponding eigenvalues. When spin is Spin.All the shape of the eigenvalues array is (n-spins, n-kpoints, n-bands) for Unpolarized and Polarized calculations, and (n-kpoints, n-bands) for Noncollinear calculations. When spin is Spin.Up or Spin.Down the shape is (n-kpoints, n-bands).

Return type:

numpy.ndarray, PhysicalQuantity of type energy.

fermiLevel(spin=None)
Parameters:

spin (Spin.Up | Spin.Down | Spin.All) – The spin the Fermi level should be returned for. Must be either Spin.Up, Spin.Down, or Spin.All. Only when the band structure is calculated with a fixed spin moment will the Fermi level depend on spin.
Default: Spin.All

Returns:

The Fermi level in absolute energy.

Return type:

PhysicalQuantity of type energy

interpolationMethod()
Returns:

The eigenvalue generation method.

Return type:

None | Full | KDotPExpansion3D

irreducibleKpoints()
Returns:

The irriducible kpoints for which the eigenvalues are calculated.

Return type:

numpy.ndarray

irreducibleWeights()

Get the weights of the irreducible k-points for which the eigenvalues are calculated.

Returns:

The weights.

Return type:

numpy.ndarray

kpoints()
Returns:

The kpoints given as input.

Return type:

class:~.MonkhorstPackGrid | RegularKpointGrid | numpy.ndarray

metatext()
Returns:

The metatext of the object or None if no metatext is present.

Return type:

str | None

method()


Deprecated: Use interpolationMethod() instead.

Returns:

The method used for the bandstructure calculation.

Return type:

None | Full | KDotPExpansion1D | KDotPExpansion3D

nlprint(stream=None)

Print a string containing an ASCII table useful for plotting the AnalysisSpin object.

Parameters:

stream (python stream) – The stream the table should be written to.
Default: NLPrintLogger()

setMetatext(metatext)

Set a given metatext string on the object.

Parameters:

metatext (str | None) – The metatext string that should be set. A value of “None” can be given to remove the current metatext.

spins()
Returns:

The spins associated to this analysis object.

Return type:

list of spin flags.

uniqueString()

Return a unique string representing the state of the object.

Usage Examples

Calculate and save the eigenvalues for a bulk silicon configuration:

# Set up lattice
lattice = FaceCenteredCubic(5.4306 * Angstrom)

# Define elements
elements = [Silicon, Silicon]

# Define coordinates
fractional_coordinates = [[ 0.  ,  0.  ,  0.  ],
                          [ 0.25,  0.25,  0.25]]

# Set up configuration
bulk_configuration = BulkConfiguration(
    bravais_lattice=lattice,
    elements=elements,
    fractional_coordinates=fractional_coordinates,
    )

# k-point sampling
k_point_sampling = MonkhorstPackGrid(
    na=3,
    nb=3,
    nc=3,
    )

numerical_accuracy_parameters = NumericalAccuracyParameters(
    k_point_sampling=k_point_sampling,
    )

# LCAO calculator
calculator = LCAOCalculator(
    numerical_accuracy_parameters=numerical_accuracy_parameters,
    )

bulk_configuration.setCalculator(calculator)
bulk_configuration.update()

# Calculate and save the eigenvalues
eigenvalues = Eigenvalues(bulk_configuration)
nlprint(eigenvalues)
nlsave('eigenvalues.hdf5', eigenvalues)

eigenvalues.py

By default the eigenvalues will be calculated for the same kpoints as used in the Calculator. It is, however, possible to calculate the eigenvalues at an arbitrary list of kpoints or MonkhorstPackGrid or RegularKpointGrid as:

# Calculate the eigenvalues at two selected k-points.
eigenvalues = Eigenvalues(
        bulk_configuration,
        kpoints=[[0.0, 0.0, 0.0], [0.1, 0.2, 0.3]])

# Calculate the eigenvalues using a MonkhorstPackGrid.
eigenvalues = Eigenvalues(
        bulk_configuration,
        kpoints=MonkhorstPackGrid(5, 5, 5))

# Calculate the eigenvalues using a RegularKpointGrid.
eigenvalues = Eigenvalues(
        bulk_configuration,
        kpoints=RegularKpointGrid(5, 5, 5))

If an Eigenvalues object has been calculated using either a MonkhorstPackGrid or RegularKpointGrid it is possible to construct a DensityOfStates object using the already calculated eigenvalues as:

# Calculate the eigenvalues using a MonkhorstPackGrid.
eigenvalues = Eigenvalues(
        bulk_configuration,
        kpoints=MonkhorstPackGrid(5, 5, 5))

# Calculate a DensityOfStates object.
density_of_states = eigenvalues.densityOfStates()