AtomicOrbital

class sisl.AtomicOrbital(*args, **kwargs)

Bases: sisl.Orbital

A projected atomic orbital consisting of real harmonics

The AtomicOrbital is a specification of the SphericalOrbital by assigning the magnetic quantum number \(m\) to the object.

AtomicOrbital should always be preferred over the SphericalOrbital because it explicitly contains all quantum numbers.

Parameters
  • *args (list of arguments) – list of arguments can be in different input options

  • q0 (float, optional) – initial charge

  • tag (str, optional) – user defined tag

R

maximum radius (in Ang)

Type

float

n

principal quantum number

Type

int

l

azimuthal quantum number

Type

int

m

magnetic quantum number

Type

int

Z

zeta shell

Type

int

P

whether this corresponds to a polarized shell (True)

Type

bool

f

interpolation function that returns f(r) for a given radius

Type

func

q0

initial electronic charge

Type

float

tag

user defined tag

Type

str

Examples

>>> r = np.linspace(0, 5, 50)
>>> f = np.exp(-r)
>>> #                    n, l, m, [Z, [P]]
>>> orb1 = AtomicOrbital(2, 1, 0, 1, (r, f))
>>> orb2 = AtomicOrbital(n=2, l=1, m=0, Z=1, (r, f))
>>> orb3 = AtomicOrbital('2pzZ', (r, f))
>>> orb4 = AtomicOrbital('2pzZ1', (r, f))
>>> orb5 = AtomicOrbital('pz', (r, f))
>>> orb2 == orb3
True
>>> orb2 == orb4
True
>>> orb2 == orb5
True

Attributes

P

R

Z

__doc__

__hash__

__module__

__slots__

f

l

m

n

orb

q0

tag

Methods

__delattr__

Implement delattr(self, name).

__dir__

Default dir() implementation.

__eq__(other)

Return self==value.

__format__

Default object formatter.

__ge__

Return self>=value.

__getattribute__

Return getattr(self, name).

__getstate__()

Return the state of this object

__gt__

Return self>value.

__init__(*args, **kwargs)

Initialize atomic orbital object

__init_subclass__

This method is called when a class is subclassed.

__le__

Return self<=value.

__lt__

Return self<value.

__ne__

Return self!=value.

__new__

Create and return a new object.

__plot__([harmonics, axes])

Plot the orbital radial/spherical harmonics

__reduce__

Helper for pickle.

__reduce_ex__

Helper for pickle.

__repr__()

Return repr(self).

__setattr__

Implement setattr(self, name, value).

__setstate__(d)

Re-create the state of this object

__sizeof__

Size of object in memory, in bytes.

__str__()

A string representation of the object

__subclasshook__

Abstract classes can override this to customize issubclass().

copy()

Create an exact copy of this object

equal(other[, psi, radial])

Compare two orbitals by comparing their radius, and possibly the radial and psi functions

name([tex])

Return named specification of the atomic orbital

psi(r)

Calculate \(\phi(\mathbf r)\) at a given point (or more points)

psi_spher(r, theta, phi[, cos_phi])

Calculate \(\phi(|\mathbf R|, \theta, \phi)\) at a given point (in spherical coordinates)

radial(r[, is_radius])

Calculate the radial part of the wavefunction \(f(\mathbf R)\)

scale(scale)

Scale the orbital by extending R by scale

set_radial(*args)

Update the internal radial function used as a \(f(|\mathbf r|)\)

spher(theta, phi[, cos_phi])

Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)

toGrid([precision, c, R, dtype, Z])

Create a Grid with only this orbital wavefunction on it

toSphere([center])

Return a sphere with radius equal to the orbital size

P
R
Z
copy()[source]

Create an exact copy of this object

equal(other, psi=False, radial=False)[source]

Compare two orbitals by comparing their radius, and possibly the radial and psi functions

Parameters
  • other (Orbital) – comparison orbital

  • psi (bool, optional) – also compare that the full psi are the same

  • radial (bool, optional) – also compare that the radial parts are the same

property f
l
m
n
name(tex=False)[source]

Return named specification of the atomic orbital

orb
psi(r)[source]

Calculate \(\phi(\mathbf r)\) at a given point (or more points)

The position r is a vector from the origin of this orbital.

Parameters

r (array_like) – the vector from the orbital origin

Returns

basis function value at point r

Return type

numpy.ndarray

psi_spher(r, theta, phi, cos_phi=False)[source]

Calculate \(\phi(|\mathbf R|, \theta, \phi)\) at a given point (in spherical coordinates)

This is equivalent to psi however, the input is given in spherical coordinates.

Parameters
  • r (array_like) – the radius from the orbital origin

  • theta (array_like) – azimuthal angle in the \(x-y\) plane (from \(x\))

  • phi (array_like) – polar angle from \(z\) axis

  • cos_phi (bool, optional) – whether phi is actually \(cos(\phi)\) which will be faster because cos is not necessary to call.

Returns

basis function value at point r

Return type

numpy.ndarray

q0
radial(r, is_radius=True)[source]

Calculate the radial part of the wavefunction \(f(\mathbf R)\)

The position r is a vector from the origin of this orbital.

Parameters
  • r (array_like) – radius from the orbital origin, for is_radius=False r must be vectors

  • is_radius (bool, optional) – whether r is a vector or the radius

Returns

radial orbital value at point r

Return type

numpy.ndarray

scale(scale)

Scale the orbital by extending R by scale

set_radial(*args)[source]

Update the internal radial function used as a \(f(|\mathbf r|)\)

See SphericalOrbital.set_radial where these arguments are passed to.

spher(theta, phi, cos_phi=False)[source]

Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)

Parameters
  • theta (array_like) – azimuthal angle in the \(x-y\) plane (from \(x\))

  • phi (array_like) – polar angle from \(z\) axis

  • cos_phi (bool, optional) – whether phi is actually \(cos(\phi)\) which will be faster because cos is not necessary to call.

Returns

spherical harmonics at angles \(\theta\) and \(\phi\)

Return type

numpy.ndarray

tag
toGrid(precision=0.05, c=1.0, R=None, dtype=<class 'numpy.float64'>, Z=1)

Create a Grid with only this orbital wavefunction on it

Parameters
  • precision (float, optional) – used separation in the Grid between voxels (in Ang)

  • c (float or complex, optional) – coefficient for the orbital

  • R (float, optional) – box size of the grid (default to the orbital range)

  • dtype (numpy.dtype, optional) – the used separation in the Grid between voxels

  • Z (int, optional) – atomic number associated with the grid

toSphere(center=None)

Return a sphere with radius equal to the orbital size

Returns

sphere with a radius equal to the radius of this orbital

Return type

Sphere