${}_{s}Y_{\ell,m}$ functions
The spin-weighted spherical harmonics are an important set of functions defined on the rotation group $๐๐(3)$, or more generally, the spin group $๐๐ฉ๐ข๐ง(3)$ that covers it. They are eigenfunctions of the left- and right-Lie derivatives, and are particularly useful in describing the angular dependence of polarized fields, like the electromagnetic field and gravitational-wave field. Originally introduced by Newman and Penrose [7], they are essentially components of Wigner's $\frak{D}$ matrices:
\[{}_{s}Y_{\ell,m}(\mathbf{R}) = (-1)^s \sqrt{\frac{2\ell+1}{4\pi}} \, \frak{D}^{(\ell)}_{m, -s}(\mathbf{R}).\]
As such, they can be computed using the same $H$ recursion algorithm as the Wigner $\frak{D}^{(\ell)}_{m, -s}$ matrices. But because not all values of $s \in -\ell:\ell$ are used, we can be much more efficient in both storage and computation time.
The user interface is very similar to the one for Wigner's $๐$ and $d$ matrices:
using Quaternionic
using SphericalFunctions
R = randn(RotorF64)
โโโโ = 8
s = -2
Y = sYlm_values(R, โโโโ, s)
Again, the results can take up a lot of memory, so for maximum efficiency when calling this function repeatedly with different R
values, it is best to pre-allocate the necessary memory with the sYlm_prep
function, and the pass that in as an argument to sYlm_values!
:
Y_storage = sYlm_prep(โโโโ, s)
Y = sYlm_values!(Y_storage, R, s)
(Beware that, as noted in the documentation for sYlm_values!
, the output Y
is just a reference to part of the Y_storage
object, so you should not reuse Y_storage
until you have copied or otherwise finished using Y
.)
The output Y
is a single vector of Complex
numbers with the same base type as R
. The ordering of the elements is described in the documentation for sYlm_values!
. It is also possible to efficiently view slices of this vector as a series of individual vectors using a sYlm_iterator
:
for (โ, Yหก) in zip(0:โโโโ, sYlm_iterator(Y, โโโโ))
# Do something with the matrix Yหก[โ+mโฒ+1, โ+m+1]
end
Docstrings
SphericalFunctions.sYlm_values
โ FunctionsYlm_values(R, โโโโ, spin)
sYlm_values(ฮธ, ฯ, โโโโ, spin)
Compute values of the spin-weighted spherical harmonic ${}_{s}Y_{\ell, m}(R)$ for all $\ell \leq \ell_\mathrm{max}$.
See sYlm_values!
for details about the input and output values.
This function only appropriate when you need to evaluate the ${}_{s}Y_{\ell, m}$ for a single value of R
or ฮธ, ฯ
because it allocates large arrays and performs many calculations that could be reused. If you need to evaluate the matrices for many values of R
or ฮธ, ฯ
, you should pre-allocate the storage with sYlm_prep
, and then call sYlm_values!
with the result instead.
SphericalFunctions.sYlm_values!
โ FunctionsYlm_values!(sYlm_storage, R, spin)
sYlm_values!(sYlm_storage, ฮธ, ฯ, spin)
sYlm_values!(sYlm, R, โโโโ, spin)
sYlm_values!(sYlm, ฮธ, ฯ, โโโโ, spin)
Compute values of the spin-weighted spherical harmonic ${}_{s}Y_{\ell, m}(R)$ for all $\ell \leq \ell_\mathrm{max}$.
The spherical harmonics of spin weight $s$ are related to Wigner's $\mathfrak{D}$ matrix as
\[\begin{aligned} {}_{s}Y_{\ell, m}(R) &= (-1)^s \sqrt{\frac{2\ell+1}{4\pi}} \mathfrak{D}^{(\ell)}_{m, -s}(R) \\ &= (-1)^s \sqrt{\frac{2\ell+1}{4\pi}} \bar{\mathfrak{D}}^{(\ell)}_{-s, m}(\bar{R}). \end{aligned}\]
In all cases, the result is returned in a 1-dimensional array ordered as
[
โYโโ(R)
for โ โ 0:โโโโ
for m โ -โ:โ
]
When the first argument is Y
, it will be modified, so it must be at least as large as that array. When the first argument is sYlm_storage
, it should be the quantity returned by sYlm_prep
, and the result will be written into the Y
field of that tuple. Both of these options โ especially the latter โ reduce the number of allocations needed on each call to the corresponding functions, which should increase the speed significantly. Note that the Y
or sYlm_storage
arguments must have types compatible with the type of R
or ฮธ, ฯ
.
When using the sYlm_storage
argument (which is recommended), the returned quantity sYlm
will be an alias of sYlm_storage[1]
. If you want to retain that data after the next call to sYlm_values!
, you should copy it with copy(sYlm)
.
The ฮธ, ฯ
arguments are spherical coordinates as described in the documentation of Quaternionic.from_spherical_coordinates
.
See also sYlm_values
for a simpler function call when you only need to evaluate the ${}_{s}Y_{\ell, m}$ for a single value of R
or ฮธ, ฯ
.
Examples
using Quaternionic, SphericalFunctions
spin = -2
โโโโ = 8
T = Float64
R = Rotor{T}(1, 2, 3, 4) # Will be normalized automatically
sYlm_storage = sYlm_prep(โโโโ, spin, T)
sYlm = sYlm_values!(sYlm_storage, R, spin)
SphericalFunctions.sYlm_prep
โ FunctionsYlm_prep(โโโโ, sโโโ, [T=Float64, [โโแตขโ=0]])
Construct storage space and pre-compute recursion coefficients to compute spin-weighted spherical-harmonic values ${}_{s}Y_{\ell, m}$ in place.
This returns the sYlm_storage
arguments needed by sYlm_values!
.
Note that the result of this function can be passed to sYlm_values!
, even if the value of spin
passed to that function is smaller (in absolute value) than the sโโโ
passed to this function. That is, the sYlm_storage
returned by this function can be used to compute ${}_{s}Y_{\ell, m}$ values for numerous values of the spin.
SphericalFunctions.sYlm_iterator
โ TypesYlm_iterator(Y, โโโโ, [โโแตขโ, [iโแตขโ]])
Construct an Iterator that returns sub-vectors of Y
, each of which consists of elements $(โ,-โ)$ through $(โ,โ)$, for $โ$ from โโแตขโ
through โโโโ
.
Note that the returned objects are views into the original Y
data โ meaning that you may alter their values.
Because the result is a vector restricted to a particular $โ$ value, you can index the $(โ, m)$ element as [โ+m+1]
. For example, you might use this as something like
for (โ, Yหก) in zip(โโแตขโ:โโโโ, sYlm_iterator(Y, โโโโ))
for m in -โ:โ
Yหก[โ+m+1] = <...>
end
end
By default, Y
is assumed to contain all possible values, beginning with (0,0)
. However, if โโแตขโ
is not 0, this can be ambiguous: do we mean that Y
really starts with the (0,0)
element and we are just asking to begin the iteration higher? Or do we mean that Y
doesn't even contain data for lower โ
values? We can resolve this using iโแตขโ
, which gives the index of โโแตขโ
in Y
. By default, we assume the first case, and set iโแตขโ=Ysize(โโแตขโ-1)+1
. However, if Y
doesn't contain data below โโแตขโ
, we could use iโแตขโ=1
to indicate the index in Y
at which to find $(โโแตขโ,-โโแตขโ)$.
Also note that no bounds checking is done, either at instantiation time or during iteration. You are responsible for ensuring that the size of Y
and the values of โโโโ
, โโแตขโ
, and iโแตขโ
make sense.