SpECTRE
v2023.01.13
|
A star-shaped surface expanded in spherical harmonics. More...
#include <Strahlkorper.hpp>
Classes | |
struct | Center |
struct | Lmax |
struct | Radius |
Public Types | |
using | options = tmpl::list< Lmax, Radius, Center > |
Public Member Functions | |
Strahlkorper (size_t l_max, size_t m_max, double radius, std::array< double, 3 > center) | |
Construct a sphere of radius radius with a given center. | |
Strahlkorper (size_t l_max, double radius, std::array< double, 3 > center) | |
Construct a sphere of radius radius , setting m_max =l_max . | |
Strahlkorper (size_t l_max, size_t m_max, const DataVector &radius_at_collocation_points, std::array< double, 3 > center) | |
Construct a Strahlkorper from a DataVector containing the radius at the collocation points. More... | |
Strahlkorper (size_t l_max, size_t m_max, const Strahlkorper &another_strahlkorper) | |
Prolong or restrict another surface to the given l_max and m_max . | |
Strahlkorper (DataVector coefs, const Strahlkorper &another_strahlkorper) | |
Construct a Strahlkorper from another Strahlkorper, but explicitly specifying the coefficients. Here coefficients are in the same storage scheme as the coefficients() member function returns. | |
Strahlkorper (DataVector coefs, Strahlkorper &&another_strahlkorper) | |
Move-construct a Strahlkorper from another Strahlkorper, explicitly specifying the coefficients. | |
void | pup (PUP::er &p) |
Serialization for Charm++. | |
const DataVector & | coefficients () const |
DataVector & | coefficients () |
const std::array< double, 3 > & | expansion_center () const |
Point about which the spectral basis of the Strahlkorper is expanded. The center is given in the frame in which the Strahlkorper is defined. This center must be somewhere inside the Strahlkorper, but in principle it can be anywhere. See physical_center() for a different measure. | |
std::array< double, 3 > | physical_center () const |
Approximate physical center (determined by \(l=1\) coefficients) Implementation of Eqs. (38)-(40) in [73]. | |
double | average_radius () const |
Average radius of the surface (determined by \(Y_{00}\) coefficient) | |
size_t | l_max () const |
Maximum \(l\) in \(Y_{lm}\) decomposition. | |
size_t | m_max () const |
Maximum \(m\) in \(Y_{lm}\) decomposition. | |
double | radius (double theta, double phi) const |
Radius at a particular angle \((\theta,\phi)\). This is inefficient if done at multiple points many times. See YlmSpherepack for alternative ways of computing this. | |
bool | point_is_contained (const std::array< double, 3 > &x) const |
Determine if a point x is contained inside the surface. The point must be given in Cartesian coordinates in the frame in which the Strahlkorper is defined. This is inefficient if done at multiple points many times. | |
const YlmSpherepack & | ylm_spherepack () const |
Static Public Attributes | |
static constexpr Options::String | help |
A star-shaped surface expanded in spherical harmonics.
Strahlkorper< Frame >::Strahlkorper | ( | size_t | l_max, |
size_t | m_max, | ||
const DataVector & | radius_at_collocation_points, | ||
std::array< double, 3 > | center | ||
) |
Construct a Strahlkorper from a DataVector containing the radius at the collocation points.
radius_at_collocation_points
. Instead, the constructed Strahlkorper will match the shape given by radius_at_collocation_points
only to order (l_max
,m_max
). This is because the YlmSpherepack representation of the Strahlkorper has more collocation points than spectral coefficients. Specifically, radius_at_collocation_points
has \((l_{\rm max} + 1) (2 m_{\rm max} + 1)\) degrees of freedom, but because there are only \(m_{\rm max}^2+(l_{\rm max}-m_{\rm max})(2m_{\rm max}+1)\) spectral coefficients, it is not possible to choose spectral coefficients to exactly match all points in radius_at_collocation_points
.
|
inline |
These coefficients are stored as SPHEREPACK coefficients. Suppose you represent a set of coefficients \(F^{lm}\) in the expansion
\[ f(\theta,\phi) = \sum_{l=0}^{l_{max}} \sum_{m=-l}^{l} F^{lm} Y^{lm}(\theta,\phi) \]
Here the \(Y^{lm}(\theta,\phi)\) are the usual complex-valued scalar spherical harmonics, so \(F^{lm}\) are also complex-valued. But here we assume that \(f(\theta,\phi)\) is real, so therefore the \(F^{lm}\) obey \(F^{l-m} = (-1)^m (F^{lm})^\star\). So one does not need to store both real and imaginary parts for both positive and negative \(m\), and the stored coefficients can all be real.
So the stored coefficients are:
\begin{align} \text{coefficients()(l,m)} &= (-1)^m \sqrt{\frac{2}{\pi}} \Re(F^{lm}) \quad \text{for} \quad m\ge 0, \\ \text{coefficients()(l,m)} &= (-1)^m \sqrt{\frac{2}{\pi}} \Im(F^{lm}) \quad \text{for} \quad m<0 \end{align}
|
staticconstexpr |