SpinTransferTorque

class SpinTransferTorque(configuration, energy=None, kpoints=None, contributions=None, self_energy_calculator=None, energy_zero_parameter=None, infinitesimal=None, z_integration_range=None, adaptive_directions=None, adaptive_method=None)

Class for representing the spin transfer torque for a given device configuration and calculator.

Parameters:
  • configuration (DeviceConfiguration.) – The configuration to calculate the spin torque transfer of.

  • energy (PhysicalQuantity of type energy) – The energy for which the spin transfer torque should be calculated.
    Default: 0.0 * eV

  • kpoints (sequence (size 3) of int | MonkhorstPackGrid | KpointDensity | AdaptiveGrid) – The k-points for which the spin transfer torque should be calculated.
    Default: The Monkhorst-Pack grid used for the self-consistent calculation.

  • contributions (Left | Right) – The density contributions to include in the spin transfer torque.
    Default: Left

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

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

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

  • z_integration_range (list of float) – The range of fractional z-coordinates used for integrating the spin transfer torque. This should be given as a list of two floats. This is only used when the parameter kpoints is an AdaptiveGrid object.
    Default: [0.0, 0.5]

  • adaptive_directions (list of CartesianDirection.{X,Y,Z}) – The directions used for adaptive integration. This should be given as 1, 2, or 3 CartesianDirection flags.
    Default: [CartesianDirection.X, CartesianDirection.Y, CartesianDirection.Z]

  • adaptive_method (RealSpace | OrbitalSpace) – The method used for adaptive integration.
    Default: OrbitalSpace

adaptive()
Returns:

True if adaptive grids are used.

Return type:

bool

adaptiveIntegralValues()
Returns:

The integrated spin-transfer torque in x, y, z direction obtained with adaptive integration. If adaptive grids are not used, this function returns None.

Return type:

numpy.array | None

adaptiveSTTValues()
Returns:

The integrated spin-transfer torque values at each k-point obtained only with adaptive integration. If adaptive grids are not used, this function returns None.

Return type:

numpy.array | None

contributions()
Returns:

The contributions parameter.

Return type:

Left | Right

energy()
Returns:

The energy parameter.

Return type:

PhysicalQuantity of type energy

energyZeroParameter()
Returns:

The energy zero parameter.

Return type:

AverageFermiLevel | AbsoluteEnergy

evaluate(x, y, z, spin=None)

Evaluate in the point x,y,z.

Parameters:
  • x (PhysicalQuantity with type length) – The cartesian x coordinate.

  • y (PhysicalQuantity with type length) – The cartesian y coordinate.

  • z (PhysicalQuantity with type length) – The cartesian z coordinate.

  • spin (Spin.All | Spin.Sum | Spin.X | Spin.Y | Spin.Z) – The spin component to project on.
    Default: The spin that the object was constructed with.

Returns:

The vector grid value at the specified point for the given spin. For Spin.All, a tuple with (Spin.Sum, Spin.X, Spin.Y, Spin.Z) components is returned.

Return type:

PhysicalQuantity

getComponent(component=None)

Get the grid component.

Parameters:

component (int) – The index of the component to extract; 0 for x, 1 for y, 2 for z.
Default: Return all

Returns:

The grid component.

Return type:

RealGrid3D | list of RealGrid3D

gridCoordinate(i, j, k)

Return the grid coordinate for a given grid index.

Parameters:
  • i (int) – The grid index in the A direction.

  • j (int) – The grid index in the B direction.

  • k (int) – The grid index in the C direction.

Returns:

The Cartesian coornate of the given grid index.

Return type:

PhysicalQuantity of type length.

kpoints()
Returns:

A list of k-points obtained with adaptive integration, each as a list of two fractional cartesian coordinates. If adaptive grids are not used, this function returns a MonkhorstPackGrid object.

Return type:

list of list | MonkhorstPackGrid

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

scale(scale)

Scale the field with a float.

Parameters:

scale (float) – The parameter to scale with.

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.

shape()
Returns:

The number of grid points in each direction.

Return type:

tuple of three int.

spinProjection(spin=None)

Query method to get a spin component of the VectorGridValue object.

Parameters:

spin (Spin.Sum | Spin.Up | Spin.Down | Spin.X | Spin.Y | Spin.Z | Spin.RealUpDown | Spin.ImagUpDown) – The spin component for which to return the GridValues object.
Default: Spin.Sum

Returns:

A new VectorGridValues object for the specified spin.

Return type:

VectorGridValues

toArray()
Returns:

The values of the grid as a numpy array slicing off any units.

Return type:

numpy.array.

uniqueString()

Return a unique string representing the state of the object.

unit()
Returns:

The unit of the data in the grid.

Return type:

A physical unit.

unitCell()
Returns:

The unit cell of the grid.

Return type:

PhysicalQuantity of type length.

volumeElement()
Returns:

The volume element of the grid represented by three vectors.

Return type:

PhysicalQuantity of type length.

Usage Examples

Calculate the SpinTransferTorque in the real-space representation of a one-dimensional carbon chain with a small gap in the middle where the electrons have to tunnel through. The spins to the left of the tunnelling gap points in a different direction than the spins in the right part.

# -------------------------------------------------------------
# Left electrode
# -------------------------------------------------------------
# Set up lattice
vector_a = [6.0, 0.0, 0.0]*Angstrom
vector_b = [0.0, 6.0, 0.0]*Angstrom
vector_c = [0.0, 0.0, 5.8]*Angstrom
left_electrode_lattice = UnitCell(vector_a, vector_b, vector_c)
# Define elements
left_electrode_elements = [Carbon, Carbon]
# Define coordinates
left_electrode_coordinates = [[ 3.  ,  3.  ,  1.45],
                              [ 3.  ,  3.  ,  4.35]]*Angstrom
# Set up configuration
left_electrode = BulkConfiguration(
    bravais_lattice=left_electrode_lattice,
    elements=left_electrode_elements,
    cartesian_coordinates=left_electrode_coordinates
    )
# -------------------------------------------------------------
# Right electrode
# -------------------------------------------------------------
# Set up lattice
vector_a = [6.0, 0.0, 0.0]*Angstrom
vector_b = [0.0, 6.0, 0.0]*Angstrom
vector_c = [0.0, 0.0, 5.8]*Angstrom
right_electrode_lattice = UnitCell(vector_a, vector_b, vector_c)
# Define elements
right_electrode_elements = [Carbon, Carbon]
# Define coordinates
right_electrode_coordinates = [[ 3.  ,  3.  ,  1.45],
                               [ 3.  ,  3.  ,  4.35]]*Angstrom
# Set up configuration
right_electrode = BulkConfiguration(
    bravais_lattice=right_electrode_lattice,
    elements=right_electrode_elements,
    cartesian_coordinates=right_electrode_coordinates
    )
# -------------------------------------------------------------
# Central region
# -------------------------------------------------------------
# Set up lattice
vector_a = [6.0, 0.0, 0.0]*Angstrom
vector_b = [0.0, 6.0, 0.0]*Angstrom
vector_c = [0.0, 0.0, 71.5]*Angstrom
central_region_lattice = UnitCell(vector_a, vector_b, vector_c)
# Define elements
central_region_elements = [Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon,
                           Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon,
                           Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon, Carbon]
# Define coordinates
central_region_coordinates = [[  3.  ,   3.  ,   1.45],
                              [  3.  ,   3.  ,   4.35],
                              [  3.  ,   3.  ,   7.25],
                              [  3.  ,   3.  ,  10.15],
                              [  3.  ,   3.  ,  13.05],
                              [  3.  ,   3.  ,  15.95],
                              [  3.  ,   3.  ,  18.85],
                              [  3.  ,   3.  ,  21.75],
                              [  3.  ,   3.  ,  24.65],
                              [  3.  ,   3.  ,  27.55],
                              [  3.  ,   3.  ,  30.45],
                              [  3.  ,   3.  ,  33.35],
                              [  3.  ,   3.  ,  38.15],
                              [  3.  ,   3.  ,  41.05],
                              [  3.  ,   3.  ,  43.95],
                              [  3.  ,   3.  ,  46.85],
                              [  3.  ,   3.  ,  49.75],
                              [  3.  ,   3.  ,  52.65],
                              [  3.  ,   3.  ,  55.55],
                              [  3.  ,   3.  ,  58.45],
                              [  3.  ,   3.  ,  61.35],
                              [  3.  ,   3.  ,  64.25],
                              [  3.  ,   3.  ,  67.15],
                              [  3.  ,   3.  ,  70.05]]*Angstrom
# Set up configuration
central_region = BulkConfiguration(
    bravais_lattice=central_region_lattice,
    elements=central_region_elements,
    cartesian_coordinates=central_region_coordinates
    )
device_configuration = DeviceConfiguration(
    central_region,
    [left_electrode, right_electrode]
    )
# -------------------------------------------------------------
# Calculator
# -------------------------------------------------------------
#----------------------------------------
# Basis Set
#----------------------------------------
basis_set = [
    LDABasis.Carbon_SingleZeta,
    ]

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

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

# Use the special noncollinear mixing scheme
iteration_control_parameters = IterationControlParameters(
    algorithm=PulayMixer(noncollinear_mixing=True)
    )

#----------------------------------------
# Electrode Calculators
#----------------------------------------
left_electrode_calculator = LCAOCalculator(
    basis_set=basis_set,
    exchange_correlation=exchange_correlation,
    poisson_solver=left_electrode_poisson_solver,
    )
right_electrode_calculator = LCAOCalculator(
    basis_set=basis_set,
    exchange_correlation=exchange_correlation,
    poisson_solver=right_electrode_poisson_solver,
    )

#----------------------------------------
# Device Calculator
#----------------------------------------
calculator = DeviceLCAOCalculator(
    basis_set=basis_set,
    exchange_correlation=exchange_correlation,
    iteration_control_parameters = iteration_control_parameters,
    electrode_calculators=
        [left_electrode_calculator, right_electrode_calculator],
    )

# Define the spin rotation
theta = 90*Degrees
left_spins = [(i, 1, 0*Degrees, 0*Degrees) for i in range(12)]
right_spins = [(i, 1, theta, 0*Degrees ) for i in range(12,24)]
spin_list = left_spins+right_spins
initial_spin = InitialSpin(scaled_spins=spin_list)

# Setup the initial state
device_configuration.setCalculator(
calculator,
    initial_spin=initial_spin,
)

# Calculate and save
device_configuration.update()
nlsave("STT.nc", device_configuration)

# -------------------------------------------------------------
# Mulliken population
# -------------------------------------------------------------
mulliken_population = MullikenPopulation(device_configuration)
nlsave("STT.nc", mulliken_population)
nlprint(mulliken_population)

#----------------------------------------
# Spin Transfer Torque
#----------------------------------------
spin_transfer_torque = SpinTransferTorque(
    configuration=device_configuration,
    energy=0.0*eV,
    kpoints=MonkhorstPackGrid(1,1),
    contributions=Left,
    self_energy_calculator=RecursionSelfEnergy(),
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1.0e-6*eV
)
nlsave("STT.nc", spin_transfer_torque)

spin_transfer_torque.py

Spin transfer torque calculations can require a very fine k-point sampling in order to get converged results which can be very computationally demanding. An alternative to the MonkhorstPack grid, which is the default, is to use an AdaptiveGrid which is automatically refined in the the regions in k-space where a fine resolution is required. A calculation, with adaptive k-point grid, can be set up as follows:

#----------------------------------------
# Adaptive Grid.
#----------------------------------------

adaptive_grid = AdaptiveGrid(
	kA_range=[-0.5, 0.5],
	kB_range=[-0.5, 0.5],
	tolerance=1e-2,
	error_measure=Relative,
	number_of_initial_levels=2,
	maximum_number_of_levels=7)

#----------------------------------------
# Spin Transfer Torque
#----------------------------------------
spin_transfer_torque = SpinTransferTorque(
    configuration=device_configuration,
    energy=0.0*eV,
    kpoints=adaptive_grid,
    contributions=Left,
    self_energy_calculator=RecursionSelfEnergy(),
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1.0e-6*eV,
    z_integration_range=[0.0, 0.5],
    adaptive_directions=[0,1,2],
    adaptive_method=RealSpace,
)
nlsave("STT.nc", spin_transfer_torque)

spin_transfer_torque_adaptive.py

The last three input arguments to the SpinTransferTorque object are special for the use of adaptive grids. They are used to define a single function value from the three-dimensional vector grid, which represents the spin transfer torque (STT). The STT is integrated over the transverse directions, x and y, and along the z-axis in an interval defined by the z_integration_range such that \(z_{min}=\) z_integration_range[0] and \(z_{max}=\) z_integration_range[1]. In fractional coordinates, this results in:

\[S_i = \int_{-0.5}^{0.5}dx\int_{-0.5}^{0.5}dy\int_{z_{min}}^{z_{max}}dz \;STT_i(x,y,z),\]

where \(STT_i\) is the \(i\)’th cartesian component of the STT. (0,1,2) correspond to x-, y-, and z-direcrions respectively. Which of the cartesian directions is included, is determined by the input argument adaptive_directions.

The last input argument, ‘adaptive_method’, can either be RealSpace or OrbitalSpace and refers to either a real-space or orbital-space evaluation of the integrated STT \(S_i\) . More information about the two approaches can be found in the tutorial Spin Transfer Torque.

Notes

More information is available in a technical note: Noncollinear spins and spin transfer torque in devices.