SurfaceBandstructure

class SurfaceBandstructure(configuration, route=None, points_per_segment=None, kpoints=None, energies=None, projection_list=None, contributions=None, self_energy_calculator=None, energy_zero_parameter=None, infinitesimal=None)

Constructor for the SurfaceBandstructure object.

Parameters:
  • configuration (DeviceConfiguration | SurfaceConfiguration) – The one- or two-probe configuration with attached calculator for which the surface band structure should be calculated.

  • route (list of str) – The route to take through the Brillouin-zone as a list of symmetry points of the unit cell, e.g. ['G', 'X', 'G']. This option is mutually exclusive to kpoints.
    Default: Unit-cell dependent route.

  • points_per_segment (Positive int) – The number of points per segment of the route. This option is mutually exclusive to kpoints.
    Default: 20.

  • kpoints (list of lists of floats) – A list of 3-dimensional fractional k-points at which to calculate the band structure e.g. [[0.0, 0.0, 0.0], [0.0, 0.0, 0.1], ...]. The shape is (K, 3) where K is the number of k-points. This option is mutually exclusive to route, and points_per_segment.
    Default: Unit-cell dependent route.

  • energies (PhysicalQuantity of type energy.) – The energies for which the spectrum of the band structure should be calculated.
    Default: numpy.linspace(-2.0, 2.0, 200)*eV

  • projection_list (ProjectionList) – A projection object list defining a projection.
    Default: If no projection list is specified, all orbitals will be used.

  • contributions (Left | Right | All) – The density contributions to include in the spectral function.
    Default: All

  • self_energy_calculator (DirectSelfEnergy | RecursionSelfEnergy | SparseRecursionSelfEnergy | KrylovSelfEnergy) – The self energy calculator to be used.
    Default: RecursionSelfEnergy(storage_strategy=NoStorage())

  • energy_zero_parameter (AverageFermiLevel | AbsoluteEnergy.) – Specifies the energy zero.
    Default: AverageFermiLevel

  • infinitesimal (PhysicalQuantity of type energy.) – Small energy, used to move the density of states calculation away from the real axis. This is only relevant for recursion-style self-energy calculators.
    Default: 1.0e-6*eV

contributions()
Returns:

The choice of electrodes from which the contributions from the transmission states are used.

Return type:

Left | Right | All

electrodeFermiLevels()

Returns the Fermi level(s) of the electrode(s). For a DeviceConfiguration the left and right electrode Fermi levels are returned, while for a SurfaceConfiguration only the Fermi level of the surface electrode is returned.

Returns:

The Fermi level(s) of the electrode(s).

Return type:

list of PhysicalQuantity of type energy

energies()
Returns:

The energies of the energy spectrum of the surface band structure.

Return type:

PhysicalQuantity of type energy.

energyZero()

The energy zero. It is zero if the energy zero parameter is equal to AbsoluteEnergy, otherwise it is determined by the electrode Fermi level(s). For a DeviceConfiguration it is the average of the left and right electrode Fermi levels. For a SurfaceConfiguration it is equal to the surface Fermi level.

Returns:

The energy zero.

Return type:

PhysicalQuantity of type energy

energyZeroParameter()
Returns:

The specified choice of the energy zero.

Return type:

AverageFermiLevel | AbsoluteEnergy

evaluate(spin=None)

Evaluate the spectral function.

Parameters:

spin (Spin) – The spin component.
Default: Spin.Sum

Returns:

The spectral function. The shape is (K, E) where E is the number of energies and K is the number of k-points. If spin is Spin.All the shape will be (S, K, E) with S being the number of spin components which is 2 for a polarized calculation and 1 otherwise.

Return type:

PhysicalQuantity of type inverse energy.

infinitesimal()
Returns:

The infinitesimal used for calculating the surface band structure.

Return type:

PhysicalQuantity of type energy.

kpoints()
Returns:

The array of 3-dimensional fractional k-points at which the spectrum is calculated. The shape is (K, 3) where K is the number of k-points.

Return type:

numpy.ndarray

metatext()
Returns:

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

Return type:

str | None

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()

projectionList()

The projections to include in the contribution to the surface band structure.

Returns:

The ProjectionList object.

Return type:

ProjectionList

route()
Returns:

The route through the Brillouin-zone as a list of symmetry points of the unit cell.

Return type:

list of str

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.

uniqueString()

Return a unique string representing the state of the object.

Notes

The SurfaceBandstructure object can be used to perform a detailed study of the projected density of states along specified k-point routes. The SurfaceBandstructure can be calculated for both device- and surface configurations as performed in [1] and [2], respectively.

Technically speaking the SurfaceBandstructure doesn’t calculate the bandstructure, but rather the device density of states, \(D(E, k)\) along a route of k-point in the plane perpendicular to the C-axis. As illustrated in [1] and [2] the resulting figures resembles band structure plots, hence the name.

Usage Example

This example script shows how to perform a SurfaceBandstructure calculation on a silver (111) surface. This structure has a characteristic surface state around the \(\Gamma\)-point as also shown in [2] (Fig. 11).

# -*- coding: utf-8 -*-
setVerbosity(MinimalLog)

# -------------------------------------------------------------
# Left Electrode
# -------------------------------------------------------------

# Set up lattice
vector_a = [2.88902617589, 0.0, 0.0]*Angstrom
vector_b = [-1.44451308795, 2.50197006052, 0.0]*Angstrom
vector_c = [0.0, 0.0, 7.07663998448]*Angstrom
left_electrode_lattice = UnitCell(vector_a, vector_b, vector_c)

# Define elements
left_electrode_elements = [Silver, Silver, Silver]

# Define coordinates
left_electrode_coordinates = [[-0.            ,  1.667980040348,  1.179540023274],
                              [-0.            ,  0.            ,  3.538420018102],
                              [ 1.444513087947,  0.833990020174,  5.89730001293 ]]*Angstrom

# Set up configuration
left_electrode = BulkConfiguration(
    bravais_lattice=left_electrode_lattice,
    elements=left_electrode_elements,
    cartesian_coordinates=left_electrode_coordinates
    )

# -------------------------------------------------------------
# Central Region
# -------------------------------------------------------------

# Set up lattice
vector_a = [2.88902617589, 0.0, 0.0]*Angstrom
vector_b = [-1.44451308795, 2.50197006052, 0.0]*Angstrom
vector_c = [0.0, 0.0, 44.7683399716]*Angstrom
central_region_lattice = UnitCell(vector_a, vector_b, vector_c)

# Define elements
central_region_elements = [Silver, Silver, Silver, Silver, Silver, Silver, Silver, Silver,
                           Silver, Silver, Silver]

# Define coordinates
central_region_coordinates = [[ -0.            ,   1.667980040348,   1.179540023274],
                              [ -0.            ,   0.            ,   3.538420018102],
                              [  1.444513087947,   0.833990020174,   5.89730001293 ],
                              [ -0.            ,   1.667980040348,   8.256180007758],
                              [ -0.            ,   0.            ,  10.615060002586],
                              [  1.444513087947,   0.833990020174,  12.973939997414],
                              [ -0.            ,   1.667980040348,  15.332819992242],
                              [ -0.            ,   0.            ,  17.69169998707 ],
                              [  1.444513087947,   0.833990020174,  20.050579981898],
                              [ -0.            ,   1.667980040348,  22.409459976726],
                              [ -0.            ,   0.            ,  24.768339971554]]*Angstrom

# Set up configuration
central_region = BulkConfiguration(
    bravais_lattice=central_region_lattice,
    elements=central_region_elements,
    cartesian_coordinates=central_region_coordinates
    )

surface_configuration = SurfaceConfiguration(
    central_region,
    left_electrode,
    equivalent_electrode_length=7.07663998448*Angstrom,
    )

# -------------------------------------------------------------
# Calculator
# -------------------------------------------------------------
#----------------------------------------
# Basis Set
#----------------------------------------
basis_set = [
    LDABasis.Silver_SingleZetaPolarized,
    ]

#----------------------------------------
# Exchange-Correlation
#----------------------------------------
exchange_correlation = LDA.PZ

#----------------------------------------
# Numerical Accuracy Settings
#----------------------------------------
left_electrode_k_point_sampling = MonkhorstPackGrid(
    na=9,
    nb=9,
    nc=231,
    )
left_electrode_numerical_accuracy_parameters = NumericalAccuracyParameters(
    density_mesh_cutoff=45.0*Hartree,
    k_point_sampling=left_electrode_k_point_sampling,
    )

device_k_point_sampling = MonkhorstPackGrid(
    na=9,
    nb=9,
    nc=231,
    )
device_numerical_accuracy_parameters = NumericalAccuracyParameters(
    density_mesh_cutoff=45.0*Hartree,
    k_point_sampling=device_k_point_sampling,
    )

#----------------------------------------
# Poisson Solver Settings
#----------------------------------------
left_electrode_poisson_solver = FastFourier2DSolver(
    boundary_conditions=[[PeriodicBoundaryCondition(),PeriodicBoundaryCondition()],
                         [PeriodicBoundaryCondition(),PeriodicBoundaryCondition()],
                         [PeriodicBoundaryCondition(),PeriodicBoundaryCondition()]]
    )

#----------------------------------------
# Contour Integral Settings
#----------------------------------------
equilibrium_contour = SemiCircleContour(
    integral_lower_bound=1.17597836274*Hartree,
    circle_points=30,
    )
contour_parameters = ContourParameters(
    equilibrium_contour=equilibrium_contour,
    )

#----------------------------------------
# Electrode Calculators
#----------------------------------------
left_electrode_calculator = LCAOCalculator(
    basis_set=basis_set,
    exchange_correlation=exchange_correlation,
    numerical_accuracy_parameters=left_electrode_numerical_accuracy_parameters,
    poisson_solver=left_electrode_poisson_solver,
    )

#----------------------------------------
# Device Calculator
#----------------------------------------
calculator = DeviceLCAOCalculator(
    basis_set=basis_set,
    exchange_correlation=exchange_correlation,
    numerical_accuracy_parameters=device_numerical_accuracy_parameters,
    contour_parameters=contour_parameters,
    electrode_calculators=
        [left_electrode_calculator, left_electrode_calculator],
    )

surface_configuration.setCalculator(calculator)
nlprint(surface_configuration)
surface_configuration.update()
nlsave('Surface Silver 111 LDA.hdf5', surface_configuration)

# -------------------------------------------------------------
# Surface Bandstructure
# -------------------------------------------------------------
surface_bandstructure = SurfaceBandstructure(
    configuration=surface_configuration,
    route=['X', 'G', 'Y'],
    points_per_segment=150,
    energies=numpy.linspace(-1.5, 1.5, 300)*eV,
    projection_list=ProjectionList(atoms=All),
    contributions=All,
    self_energy_calculator=RecursionSelfEnergy(storage_strategy=NoStorage()),
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1e-03*eV,
    )
nlsave('Surface Silver 111 LDA.hdf5', surface_bandstructure)
nlprint(surface_bandstructure)

surface_bandstructure.py

Projections

When constructing the SurfaceBandstructure it is possible to specify a ProjectionList object with information about which atoms and orbital shells to project on. By default all atoms and orbitals are included. If only the contribution from the two outermost atoms should be included, the following should be set in the above script:

projection_list=ProjectionList(atoms=[9, 10])

when setting up the SurfaceBandstructure object.