sisl.SphericalOrbital
- class sisl.SphericalOrbital(l, rf_or_func, q0=0.0, tag='', **kwargs)
Bases:
Orbital
An arbitrary orbital class which only contains the harmonical part of the wavefunction where \(\phi(\mathbf r)=f(|\mathbf r|)Y_l^m(\theta,\varphi)\)
Note that in this case the used spherical harmonics is:
\[Y^m_l(\theta,\varphi) = (-1)^m\sqrt{\frac{2l+1}{4\pi} \frac{(l-m)!}{(l+m)!}} e^{i m \theta} P^m_l(\cos(\varphi))\]The resulting orbital is
\[\phi_{lmn}(\mathbf r) = f(|\mathbf r|) Y^m_l(\theta, \varphi)\]where typically \(f(|\mathbf r|)\equiv\phi_{ln}(|\mathbf r|)\). The above equation clarifies that this class is only intended for each \(l\), and that subsequent \(m\) orders may be extracted by altering the spherical harmonic. Also, the quantum number \(n\) is not necessary as that value is implicit in the \(\phi_{ln}(|\mathbf r|)\) function.
- Parameters
- f
interpolation function that returns f(r) for a given radius
- Type
func
Examples
>>> from scipy.interpolate import interp1d >>> orb = SphericalOrbital(1, (np.arange(10.), np.arange(10.))) >>> orb.equal(SphericalOrbital(1, interp1d(np.arange(10.), np.arange(10.), ... fill_value=(0., 0.), kind="cubic", bounds_error=False))) True
Methods
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 a named specification of the orbital (tag)
psi
(r[, m])Calculate \(\phi(\mathbf R)\) at a given point (or more points)
psi_spher
(r, theta, phi[, m, cos_phi])Calculate \(\phi(|\mathbf R|, \theta, \phi)\) at a given point (in spherical coordinates)
radial
(r, *args, **kwargs)Calculate the radial part of spherical orbital \(R(\mathbf r)\)
scale
(scale)Scale the orbital by extending R by scale
set_radial
(*args, **kwargs)Update the internal radial function used as a \(f(|\mathbf r|)\)
spher
(theta, phi[, m, cos_phi])Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)
toAtomicOrbital
([m, n, zeta, P, q0])Create a list of AtomicOrbital objects
toGrid
([precision, c, R, dtype, atom])Create a Grid with only this orbital wavefunction on it
toSphere
([center])Return a sphere with radius equal to the orbital size
Maxmimum radius of orbital
\(l\) quantum number
Initial charge
Named tag of orbital
- property R
Maxmimum radius of orbital
- equal(other, psi=False, radial=False)[source]
Compare two orbitals by comparing their radius, and possibly the radial and psi functions
- property l
\(l\) quantum number
- name(tex=False)
Return a named specification of the orbital (tag)
- psi(r, m=0)[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 of (:, 3)) – vector from the orbital origin
m (int, optional) – magnetic quantum number, must be in range
-self.l <= m <= self.l
- Returns
basis function value at point r
- Return type
- psi_spher(r, theta, phi, m=0, 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
m (int, optional) – magnetic quantum number, must be in range
-self.l <= m <= self.l
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
- property q0
Initial charge
- radial(r, *args, **kwargs)
Calculate the radial part of spherical orbital \(R(\mathbf r)\)
The position r is a vector from the origin of this orbital.
- Parameters
r (array_like) – radius from the orbital origin
*args – arguments passed to the radial function
**args – keyword arguments passed to the radial function
- Returns
radial orbital value at point r
- Return type
- scale(scale)
Scale the orbital by extending R by scale
- set_radial(*args, **kwargs)
Update the internal radial function used as a \(f(|\mathbf r|)\)
This can be called in several ways:
- set_radial(r, f)
which uses
scipy.interpolate.UnivariateSpline(r, f, k=3, s=0, ext=1, check_finite=False)
to define the interpolation function (see interp keyword). Here the maximum radius of the orbital is the maximum r value, regardless off(r)
is zero for smaller r.- set_radial(func)
which sets the interpolation function directly. The maximum orbital range is determined automatically to a precision of 0.0001 AA.
- Parameters
r (numpy.ndarray) – the radial positions and the radial function values at r.
f (numpy.ndarray) – the radial positions and the radial function values at r.
func (callable) – a function which enables evaluation of the radial function. The function should accept a single array and return a single array.
interp (callable, optional) – When two non-keyword arguments are passed this keyword will be used. It is the interpolation function which should return the equivalent of func. By using this one can define a custom interpolation routine. It should accept two arguments,
interp(r, f)
and return a callable that returns interpolation values. See examples for different interpolation routines.
- Return type
Examples
>>> from scipy import interpolate as interp >>> o = SphericalOrbital(1, lambda x:x) >>> r = np.linspace(0, 4, 300) >>> f = np.exp(-r) >>> def i_univariate(r, f): ... return interp.UnivariateSpline(r, f, k=3, s=0, ext=1, check_finite=False) >>> def i_interp1d(r, f): ... return interp.interp1d(r, f, kind="cubic", fill_value=(f[0], 0.), bounds_error=False) >>> def i_spline(r, f): ... from functools import partial ... tck = interp.splrep(r, f, k=3, s=0) ... return partial(interp.splev, tck=tck, der=0, ext=1) >>> R = np.linspace(0, 4, 400) >>> o.set_radial(r, f, interp=i_univariate) >>> f_univariate = o.radial(R) >>> o.set_radial(r, f, interp=i_interp1d) >>> f_interp1d = o.radial(R) >>> o.set_radial(r, f, interp=i_spline) >>> f_spline = o.radial(R) >>> np.allclose(f_univariate, f_interp1d) True >>> np.allclose(f_univariate, f_spline) True
- spher(theta, phi, m=0, 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
m (int, optional) – magnetic quantum number, must be in range
-self.l <= m <= self.l
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\) and given quantum number m
- Return type
- property tag
Named tag of orbital
- toAtomicOrbital(m=None, n=None, zeta=1, P=False, q0=None)[source]
Create a list of AtomicOrbital objects
This defaults to create a list of AtomicOrbital objects for every m (for m in -l:l). One may optionally specify the sub-set of m to retrieve.
- Parameters
m (int or list or None) – if
None
it defaults to-l:l
, else only for the requested mzeta (int, optional) – the specified zeta-shell
n (int, optional) – specify the \(n\) quantum number
P (bool, optional) – whether the orbitals are polarized.
q0 (float, optional) – the initial charge per orbital, initially \(q_0 / (2l+1)\) with \(q_0\) from this object
- Returns
AtomicOrbital (for passed m an atomic orbital will be returned)
list of AtomicOrbital (for each \(m\in[-l;l]\) an atomic orbital will be returned in the list)
- toGrid(precision=0.05, c=1.0, R=None, dtype=<class 'numpy.float64'>, atom=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
atom (optional) – atom associated with the grid; either an atom instance or something that
Atom(atom)
would convert to a proper atom.