SpECTRE  v2024.04.12
ylm Namespace Reference

Items related to spherical harmonics. More...

Namespaces

namespace  Tags
 Holds tags and ComputeItems associated with a ylm::Strahlkorper.
 

Classes

class  Spherepack
 Defines the C++ interface to SPHEREPACK. More...
 
class  SpherepackIterator
 Iterates over spectral coefficients stored in SPHEREPACK format. More...
 
class  Strahlkorper
 A star-shaped surface expanded in spherical harmonics. More...
 

Functions

template<typename Frame >
void change_expansion_center_of_strahlkorper (gsl::not_null< Strahlkorper< Frame > * > strahlkorper, const std::array< double, 3 > &new_center)
 Changes the expansion center of a Strahlkorper, where the expansion center is defined as the point about which the spectral basis of the Strahlkorper is expanded, which is the quantity returned by Strahlkorper::expansion_center().
 
template<typename Frame >
void change_expansion_center_of_strahlkorper_to_physical (gsl::not_null< Strahlkorper< Frame > * > strahlkorper)
 Changes the expansion center of a Strahlkorper to the physical center. Because Strahlkorper::physical_center() returns only an approximate quantity, change_expansion_center_of_strahlkorper_to_physical is iterative, and does not return exactly the same result as passing Strahlkorper::physical_center() to change_expansion_center_of_strahlkorper.
 
template<typename Frame >
void fill_ylm_legend_and_data (gsl::not_null< std::vector< std::string > * > legend, gsl::not_null< std::vector< double > * > data, const ylm::Strahlkorper< Frame > &strahlkorper, double time, size_t max_l)
 Fills the legend and row of spherical harmonic data to write to disk. More...
 
template<typename Frame >
std::vector< ylm::Strahlkorper< Frame > > read_surface_ylm (const std::string &file_name, const std::string &surface_subfile_name, size_t requested_number_of_times_from_end)
 Returns a list of ylm::Strahlkorpers constructed from reading in spherical harmonic data for a surface at a requested list of times. More...
 
template<typename Frame >
ylm::Strahlkorper< Frame > read_surface_ylm_single_time (const std::string &file_name, const std::string &surface_subfile_name, double time, double relative_epsilon)
 Similar to ylm::read_surface_ylm, this reads in spherical harmonic data for a surface and constructs a ylm::Strahlkorper. However, this function only does it at a specific time and returns a single ylm::Strahlkorper. More...
 
bool operator== (const Spherepack &lhs, const Spherepack &rhs)
 
bool operator!= (const Spherepack &lhs, const Spherepack &rhs)
 
bool operator== (const SpherepackIterator &lhs, const SpherepackIterator &rhs)
 
bool operator!= (const SpherepackIterator &lhs, const SpherepackIterator &rhs)
 
std::ostreamoperator<< (std::ostream &os, const SpherepackIterator::CoefficientArray &coefficient_array)
 
template<typename Frame >
bool operator== (const Strahlkorper< Frame > &lhs, const Strahlkorper< Frame > &rhs)
 
template<typename Frame >
bool operator!= (const Strahlkorper< Frame > &lhs, const Strahlkorper< Frame > &rhs)
 
void real_spherical_harmonic (gsl::not_null< DataVector * > spherical_harmonic, const DataVector &theta, const DataVector &phi, size_t l, int m)
 Evaluates a real spherical harmonic of order l at the requested angles \(\theta\) and \(\phi\). More...
 
DataVector real_spherical_harmonic (const DataVector &theta, const DataVector &phi, size_t l, int m)
 Evaluates a real spherical harmonic of order l at the requested angles \(\theta\) and \(\phi\). More...
 
template<typename Fr >
tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > theta_phi (const Strahlkorper< Fr > &strahlkorper)
 
template<typename Fr >
void theta_phi (const gsl::not_null< tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > * > theta_phi, const Strahlkorper< Fr > &strahlkorper)
 
template<typename Fr >
tnsr::i< DataVector, 3, Fr > rhat (const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
void rhat (const gsl::not_null< tnsr::i< DataVector, 3, Fr > * > r_hat, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
ylm::Tags::aliases::Jacobian< Fr > jacobian (const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
void jacobian (const gsl::not_null< ylm::Tags::aliases::Jacobian< Fr > * > jac, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
ylm::Tags::aliases::InvJacobian< Fr > inv_jacobian (const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
void inv_jacobian (const gsl::not_null< ylm::Tags::aliases::InvJacobian< Fr > * > inv_jac, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
ylm::Tags::aliases::InvHessian< Fr > inv_hessian (const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
void inv_hessian (const gsl::not_null< ylm::Tags::aliases::InvHessian< Fr > * > inv_hess, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
Scalar< DataVectorradius (const Strahlkorper< Fr > &strahlkorper)
 
template<typename Fr >
void radius (const gsl::not_null< Scalar< DataVector > * > result, const Strahlkorper< Fr > &strahlkorper)
 
template<typename Fr >
tnsr::I< DataVector, 3, Fr > cartesian_coords (const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius, const tnsr::i< DataVector, 3, Fr > &r_hat)
 
template<typename Fr >
void cartesian_coords (const gsl::not_null< tnsr::I< DataVector, 3, Fr > * > coords, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius, const tnsr::i< DataVector, 3, Fr > &r_hat)
 
template<typename Fr >
tnsr::I< DataVector, 3, Fr > cartesian_coords (const Strahlkorper< Fr > &strahlkorper)
 
template<typename Fr >
tnsr::i< DataVector, 3, Fr > cartesian_derivs_of_scalar (const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius_of_strahlkorper, const ylm::Tags::aliases::InvJacobian< Fr > &inv_jac)
 
template<typename Fr >
void cartesian_derivs_of_scalar (const gsl::not_null< tnsr::i< DataVector, 3, Fr > * > dx_scalar, const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius_of_strahlkorper, const ylm::Tags::aliases::InvJacobian< Fr > &inv_jac)
 
template<typename Fr >
tnsr::ii< DataVector, 3, Fr > cartesian_second_derivs_of_scalar (const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius_of_strahlkorper, const ylm::Tags::aliases::InvJacobian< Fr > &inv_jac, const ylm::Tags::aliases::InvHessian< Fr > &inv_hess)
 
template<typename Fr >
void cartesian_second_derivs_of_scalar (const gsl::not_null< tnsr::ii< DataVector, 3, Fr > * > d2x_scalar, const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius_of_strahlkorper, const ylm::Tags::aliases::InvJacobian< Fr > &inv_jac, const ylm::Tags::aliases::InvHessian< Fr > &inv_hess)
 
template<typename Fr >
Scalar< DataVectorlaplacian_of_scalar (const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
void laplacian_of_scalar (const gsl::not_null< Scalar< DataVector > * > laplacian, const Scalar< DataVector > &scalar, const Strahlkorper< Fr > &strahlkorper, const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &theta_phi)
 
template<typename Fr >
ylm::Tags::aliases::Jacobian< Fr > tangents (const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius, const tnsr::i< DataVector, 3, Fr > &r_hat, const ylm::Tags::aliases::Jacobian< Fr > &jac)
 
template<typename Fr >
void tangents (const gsl::not_null< ylm::Tags::aliases::Jacobian< Fr > * > result, const Strahlkorper< Fr > &strahlkorper, const Scalar< DataVector > &radius, const tnsr::i< DataVector, 3, Fr > &r_hat, const ylm::Tags::aliases::Jacobian< Fr > &jac)
 
template<typename Fr >
tnsr::i< DataVector, 3, Fr > normal_one_form (const tnsr::i< DataVector, 3, Fr > &dx_radius, const tnsr::i< DataVector, 3, Fr > &r_hat)
 
template<typename Fr >
void normal_one_form (const gsl::not_null< tnsr::i< DataVector, 3, Fr > * > one_form, const tnsr::i< DataVector, 3, Fr > &dx_radius, const tnsr::i< DataVector, 3, Fr > &r_hat)
 
template<typename Fr >
std::vector< std::array< double, 4 > > fit_ylm_coeffs (const DataVector &times, const std::vector< Strahlkorper< Fr > > &strahlkorpers)
 
template<typename Frame >
void time_deriv_of_strahlkorper (gsl::not_null< Strahlkorper< Frame > * > time_deriv, const std::deque< std::pair< double, Strahlkorper< Frame > > > &previous_strahlkorpers)
 Compute the time derivative of a Strahlkorper from a number of previous Strahlkorpers. More...
 
Scalar< double > ylm_to_stf_0 (const ModalVector &l0_coefs)
 Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l. More...
 
template<typename Frame >
tnsr::i< double, 3, Frame > ylm_to_stf_1 (const ModalVector &l1_coefs)
 Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l. More...
 
template<typename Frame >
tnsr::ii< double, 3, Frame > ylm_to_stf_2 (const ModalVector &l2_coefs)
 Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l. More...
 

Detailed Description

Items related to spherical harmonics.

Contains functions that depend on a Strahlkorper but not on a metric.

Function Documentation

◆ cartesian_coords() [1/3]

template<typename Fr >
void ylm::cartesian_coords ( const gsl::not_null< tnsr::I< DataVector, 3, Fr > * >  coords,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat 
)
Parameters
coordsThe returned Cartesian coordinates.
strahlkorperThe Strahlkorper surface.
radiusThe radius as a function of angle, as returned by ylm::radius.
r_hatThe Euclidean radial unit vector as returned by ylm::rhat.

◆ cartesian_coords() [2/3]

template<typename Fr >
tnsr::I< DataVector, 3, Fr > ylm::cartesian_coords ( const Strahlkorper< Fr > &  strahlkorper)

This overload computes radius, theta_phi, and r_hat internally. Use the other overloads if you already have these quantities.

Parameters
strahlkorperThe Strahlkorper surface.

◆ cartesian_coords() [3/3]

template<typename Fr >
tnsr::I< DataVector, 3, Fr > ylm::cartesian_coords ( const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat 
)

cartesian_coords(i) is \(x_{\rm surf}^i\), the vector of \((x,y,z)\) coordinates of each point on the Strahlkorper surface.

Parameters
strahlkorperThe Strahlkorper surface.
radiusThe radius as a function of angle, as returned by ylm::radius.
r_hatThe Euclidean radial unit vector as returned by ylm::rhat.

◆ cartesian_derivs_of_scalar() [1/2]

template<typename Fr >
void ylm::cartesian_derivs_of_scalar ( const gsl::not_null< tnsr::i< DataVector, 3, Fr > * >  dx_scalar,
const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius_of_strahlkorper,
const ylm::Tags::aliases::InvJacobian< Fr > &  inv_jac 
)
Parameters
dx_scalarThe returned derivatives of the scalar.
scalarThe scalar to be differentiated.
strahlkorperThe Strahlkorper surface.
radius_of_strahlkorperThe radius of the Strahlkorper at each point, as returned by ylm::radius.
inv_jacThe inverse Jacobian as returned by ylm::inv_jacobian

◆ cartesian_derivs_of_scalar() [2/2]

template<typename Fr >
tnsr::i< DataVector, 3, Fr > ylm::cartesian_derivs_of_scalar ( const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius_of_strahlkorper,
const ylm::Tags::aliases::InvJacobian< Fr > &  inv_jac 
)

dx_scalar(i) is \(\partial f/\partial x^i\) evaluated on the surface. Here \(f=f(r,\theta,\phi)=f(\theta,\phi)\) is some scalar function independent of the radial coordinate. \(f\) is considered a function of Cartesian coordinates \(f=f(\theta(x,y,z),\phi(x,y,z))\) for this operation.

Parameters
scalarThe scalar to be differentiated.
strahlkorperThe Strahlkorper surface.
radius_of_strahlkorperThe radius of the Strahlkorper at each point, as returned by ylm::radius.
inv_jacThe inverse Jacobian as returned by ylm::inv_jacobian

◆ cartesian_second_derivs_of_scalar() [1/2]

template<typename Fr >
void ylm::cartesian_second_derivs_of_scalar ( const gsl::not_null< tnsr::ii< DataVector, 3, Fr > * >  d2x_scalar,
const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius_of_strahlkorper,
const ylm::Tags::aliases::InvJacobian< Fr > &  inv_jac,
const ylm::Tags::aliases::InvHessian< Fr > &  inv_hess 
)
Parameters
d2x_scalarThe returned 2nd derivatives of the scalar.
scalarThe scalar to be differentiated.
strahlkorperThe Strahlkorper surface.
radius_of_strahlkorperThe radius of the Strahlkorper at each point, as returned by ylm::radius.
inv_jacThe inverse Jacobian as returned by ylm::inv_jacobian
inv_hessThe inverse Hessian as returned by `ylminv_hessian.

◆ cartesian_second_derivs_of_scalar() [2/2]

template<typename Fr >
tnsr::ii< DataVector, 3, Fr > ylm::cartesian_second_derivs_of_scalar ( const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius_of_strahlkorper,
const ylm::Tags::aliases::InvJacobian< Fr > &  inv_jac,
const ylm::Tags::aliases::InvHessian< Fr > &  inv_hess 
)

d2x_scalar(i,j) is \(\partial^2 f/\partial x^i\partial x^j\) evaluated on the surface. Here \(f=f(r,\theta,\phi)=f(\theta,\phi)\) is some scalar function independent of the radial coordinate. \(f\) is considered a function of Cartesian coordinates \(f=f(\theta(x,y,z),\phi(x,y,z))\) for this operation.

Parameters
scalarThe scalar to be differentiated.
strahlkorperThe Strahlkorper surface.
radius_of_strahlkorperThe radius of the Strahlkorper at each point, as returned by ylm::radius.
inv_jacThe inverse Jacobian as returned by ylm::inv_jacobian
inv_hessThe inverse Hessian as returned by `ylminv_hessian.

◆ fill_ylm_legend_and_data()

template<typename Frame >
void ylm::fill_ylm_legend_and_data ( gsl::not_null< std::vector< std::string > * >  legend,
gsl::not_null< std::vector< double > * >  data,
const ylm::Strahlkorper< Frame > &  strahlkorper,
double  time,
size_t  max_l 
)

Fills the legend and row of spherical harmonic data to write to disk.

The number of coefficients to write is based on max_l, the maximum value that the input strahlkorper could possibly have. When strahlkorper.l_max() < max_l, coefficients with \(l\) higher than strahlkorper.l_max() will simply be zero. Assuming the same max_l is always used for a given surface, we will always write the same number of columns for each row, as max_l sets the number of columns to write

◆ fit_ylm_coeffs()

template<typename Fr >
std::vector< std::array< double, 4 > > ylm::fit_ylm_coeffs ( const DataVector times,
const std::vector< Strahlkorper< Fr > > &  strahlkorpers 
)

The linear least squares fit of the polynomial of order 3 given a std::vector of Strahlkorpers to their \(Y_l^m\) coefficients. Assumes the the \(l_{\max}\) and \(m_{\max}\) of each Strahlkorper are the same, and the returned vector consists of \(2l_{\max}m_{\max}\) (the number of \(Y_l^m\) coefficients) std::array<double, 4>s, each of which consists of the four coefficients that define the best fit cubic to each \(Y_l^m\) coefficient of the Strahlkorper as a function of time.

Parameters
timesThe time corresponding to each Strahlkorper to be fit to.
strahlkorpersThe Strahlkorper surfaces which consists of a set of \(Y_l^m\) coefficients corresponding to the shape of the Strahlkorper at a particular time.

◆ inv_hessian() [1/2]

template<typename Fr >
void ylm::inv_hessian ( const gsl::not_null< ylm::Tags::aliases::InvHessian< Fr > * >  inv_hess,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

InvHessian(k,i,j) is \(r\partial (J^{-1}){}^k_j/\partial x^i\), where \((J^{-1}){}^k_j\) is the inverse Jacobian. Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). InvHessian is not symmetric because the Jacobians are Pfaffian. InvHessian doesn't depend on the shape of the surface.

◆ inv_hessian() [2/2]

template<typename Fr >
ylm::Tags::aliases::InvHessian< Fr > ylm::inv_hessian ( const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi)

InvHessian(k,i,j) is \(r\partial (J^{-1}){}^k_j/\partial x^i\), where \((J^{-1}){}^k_j\) is the inverse Jacobian. Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). InvHessian is not symmetric because the Jacobians are Pfaffian. InvHessian doesn't depend on the shape of the surface.

◆ inv_jacobian() [1/2]

template<typename Fr >
void ylm::inv_jacobian ( const gsl::not_null< ylm::Tags::aliases::InvJacobian< Fr > * >  inv_jac,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

InvJacobian(0,i) is \(r\partial\theta/\partial x^i\), and InvJacobian(1,i) is \(r\sin\theta\partial\phi/\partial x^i\). Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). InvJacobian doesn't depend on the shape of the surface.

◆ inv_jacobian() [2/2]

template<typename Fr >
ylm::Tags::aliases::InvJacobian< Fr > ylm::inv_jacobian ( const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi)

InvJacobian(0,i) is \(r\partial\theta/\partial x^i\), and InvJacobian(1,i) is \(r\sin\theta\partial\phi/\partial x^i\). Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). InvJacobian doesn't depend on the shape of the surface.

◆ jacobian() [1/2]

template<typename Fr >
void ylm::jacobian ( const gsl::not_null< ylm::Tags::aliases::Jacobian< Fr > * >  jac,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

Jacobian(i,0) is \(\frac{1}{r}\partial x^i/\partial\theta\), and Jacobian(i,1) is \(\frac{1}{r\sin\theta}\partial x^i/\partial\phi\). Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). Jacobian doesn't depend on the shape of the surface.

◆ jacobian() [2/2]

template<typename Fr >
ylm::Tags::aliases::Jacobian< Fr > ylm::jacobian ( const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi)

Jacobian(i,0) is \(\frac{1}{r}\partial x^i/\partial\theta\), and Jacobian(i,1) is \(\frac{1}{r\sin\theta}\partial x^i/\partial\phi\). Here \(r\) means \(\sqrt{x^2+y^2+z^2}\). Jacobian doesn't depend on the shape of the surface.

◆ laplacian_of_scalar() [1/2]

template<typename Fr >
void ylm::laplacian_of_scalar ( const gsl::not_null< Scalar< DataVector > * >  laplacian,
const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

\(\nabla^2 f\), the flat Laplacian of a scalar \(f\) on the surface. This is \(\eta^{ij}\partial^2 f/\partial x^i\partial x^j\), where \(f=f(r,\theta,\phi)=f(\theta,\phi)\) is some scalar function independent of the radial coordinate. \(f\) is considered a function of Cartesian coordinates \(f=f(\theta(x,y,z),\phi(x,y,z))\) for this operation.

◆ laplacian_of_scalar() [2/2]

template<typename Fr >
Scalar< DataVector > ylm::laplacian_of_scalar ( const Scalar< DataVector > &  scalar,
const Strahlkorper< Fr > &  strahlkorper,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

\(\nabla^2 f\), the flat Laplacian of a scalar \(f\) on the surface. This is \(\eta^{ij}\partial^2 f/\partial x^i\partial x^j\), where \(f=f(r,\theta,\phi)=f(\theta,\phi)\) is some scalar function independent of the radial coordinate. \(f\) is considered a function of Cartesian coordinates \(f=f(\theta(x,y,z),\phi(x,y,z))\) for this operation.

◆ normal_one_form() [1/2]

template<typename Fr >
void ylm::normal_one_form ( const gsl::not_null< tnsr::i< DataVector, 3, Fr > * >  one_form,
const tnsr::i< DataVector, 3, Fr > &  dx_radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat 
)
Parameters
one_formThe returned normal one form.
dx_radiusThe Cartesian derivatives of the radius, as returned by ylm::cartesian_derivs_of_scalar with ylm::radius passed in as the scalar.
r_hatThe radial unit vector as returned by ylm::rhat.

◆ normal_one_form() [2/2]

template<typename Fr >
tnsr::i< DataVector, 3, Fr > ylm::normal_one_form ( const tnsr::i< DataVector, 3, Fr > &  dx_radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat 
)

normal_one_form(i) is \(s_i\), the (unnormalized) normal one-form to the surface, expressed in Cartesian components. This is computed by \(x_i/r-\partial r_{\rm surf}/\partial x^i\), where \(x_i/r\) is Rhat and \(\partial r_{\rm surf}/\partial x^i\) is DxRadius. See Eq. (8) of [14]. Note on the word "normal": \(s_i\) points in the correct direction (it is "normal" to the surface), but it does not have unit length (it is not "normalized"; normalization requires a metric).

Parameters
dx_radiusThe Cartesian derivatives of the radius, as returned by ylm::cartesian_derivs_of_scalar with ylm::radius passed in as the scalar.
r_hatThe radial unit vector as returned by ylm::rhat.

◆ radius() [1/2]

template<typename Fr >
void ylm::radius ( const gsl::not_null< Scalar< DataVector > * >  result,
const Strahlkorper< Fr > &  strahlkorper 
)

(Euclidean) distance \(r_{\rm surf}(\theta,\phi)\) from the expansion center to each point of the Strahlkorper surface.

◆ radius() [2/2]

template<typename Fr >
Scalar< DataVector > ylm::radius ( const Strahlkorper< Fr > &  strahlkorper)

(Euclidean) distance \(r_{\rm surf}(\theta,\phi)\) from the expansion center to each point of the Strahlkorper surface.

◆ read_surface_ylm_single_time()

template<typename Frame >
ylm::Strahlkorper< Frame > ylm::read_surface_ylm_single_time ( const std::string file_name,
const std::string surface_subfile_name,
double  time,
double  relative_epsilon 
)

Similar to ylm::read_surface_ylm, this reads in spherical harmonic data for a surface and constructs a ylm::Strahlkorper. However, this function only does it at a specific time and returns a single ylm::Strahlkorper.

Note
If two times are found within epsilon of the time, then an error will occur. Similarly, if no time is found within the epsilon, then an error will occur as well.
Parameters
file_namename of the H5 file containing the surface's spherical harmonic data
surface_subfile_namename of the subfile (with no leading slash nor the .dat extension) within file_name that contains the surface's spherical harmonic data to read in
timeTime to read the coefficients at.
relative_epsilonHow much error is allowed when looking for a specific time. This is useful so users don't have to know the specific time to machine precision.

◆ rhat() [1/2]

template<typename Fr >
void ylm::rhat ( const gsl::not_null< tnsr::i< DataVector, 3, Fr > * >  r_hat,
const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi 
)

r_hat(i) is \(\hat{r}_i = x_i/\sqrt{x^2+y^2+z^2}\) on the Strahlkorper surface. Doesn't depend on the shape of the surface.

We need to choose upper vs lower indices for rhat; it doesn't matter because rhat is a quantity defined with a Euclidean metric, so we choose the lower index arbitrarily.

◆ rhat() [2/2]

template<typename Fr >
tnsr::i< DataVector, 3, Fr > ylm::rhat ( const tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > &  theta_phi)

r_hat(i) is \(\hat{r}_i = x_i/\sqrt{x^2+y^2+z^2}\) on the Strahlkorper surface. Doesn't depend on the shape of the surface.

We need to choose upper vs lower indices for rhat; it doesn't matter because rhat is a quantity defined with a Euclidean metric, so we choose the lower index arbitrarily.

◆ tangents() [1/2]

template<typename Fr >
void ylm::tangents ( const gsl::not_null< ylm::Tags::aliases::Jacobian< Fr > * >  result,
const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat,
const ylm::Tags::aliases::Jacobian< Fr > &  jac 
)
Parameters
resultThe computed tangent vectors.
strahlkorperThe Strahlkorper surface.
radiusThe radius of the Strahlkorper at each point, as returned by ylm::radius.
r_hatThe radial unit vector as returned by ylm::rhat.
jacThe jacobian as returned by ylm::jacobian.

◆ tangents() [2/2]

template<typename Fr >
ylm::Tags::aliases::Jacobian< Fr > ylm::tangents ( const Strahlkorper< Fr > &  strahlkorper,
const Scalar< DataVector > &  radius,
const tnsr::i< DataVector, 3, Fr > &  r_hat,
const ylm::Tags::aliases::Jacobian< Fr > &  jac 
)

tangents(i,j) is \(\partial x_{\rm surf}^i/\partial q^j\), where \(x_{\rm surf}^i\) are the Cartesian coordinates of the surface (i.e. cartesian_coords) and are considered functions of \((\theta,\phi)\).

\(\partial/\partial q^0\) means \(\partial/\partial\theta\); and \(\partial/\partial q^1\) means \(\csc\theta\,\,\partial/\partial\phi\). Note that the vectors tangents(i,0) and tangents(i,1) are orthogonal to the normal_one_form \(s_i\), i.e. \(s_i \partial x_{\rm surf}^i/\partial q^j = 0\); this statement is independent of a metric. Also, tangents(i,0) and tangents(i,1) are not necessarily orthogonal to each other, since orthogonality between 2 vectors (as opposed to a vector and a one-form) is metric-dependent.

Parameters
strahlkorperThe Strahlkorper surface.
radiusThe radius of the Strahlkorper at each point, as returned by ylm::radius.
r_hatThe radial unit vector as returned by ylm::rhat.
jacThe jacobian as returned by ylm::jacobian.

◆ theta_phi() [1/2]

template<typename Fr >
void ylm::theta_phi ( const gsl::not_null< tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > * >  theta_phi,
const Strahlkorper< Fr > &  strahlkorper 
)

\((\theta,\phi)\) on the Strahlkorper surface. Doesn't depend on the shape of the surface.

We need to choose upper vs lower indices for theta_phi; it doesn't matter because these are coordinates and not geometric objects, so we choose lower indices arbitrarily.

◆ theta_phi() [2/2]

template<typename Fr >
tnsr::i< DataVector, 2, ::Frame::Spherical< Fr > > ylm::theta_phi ( const Strahlkorper< Fr > &  strahlkorper)

\((\theta,\phi)\) on the Strahlkorper surface. Doesn't depend on the shape of the surface.

We need to choose upper vs lower indices for theta_phi; it doesn't matter because these are coordinates and not geometric objects, so we choose lower indices arbitrarily.

◆ time_deriv_of_strahlkorper()

template<typename Frame >
void ylm::time_deriv_of_strahlkorper ( gsl::not_null< Strahlkorper< Frame > * >  time_deriv,
const std::deque< std::pair< double, Strahlkorper< Frame > > > &  previous_strahlkorpers 
)

Compute the time derivative of a Strahlkorper from a number of previous Strahlkorpers.

Details

Does simple 1D FD with non-uniform spacing using fd::non_uniform_1d_weights.

Parameters
time_derivStrahlkorper whose coefficients are the time derivative of previous_strahlkorpers' coefficients.
previous_strahlkorpersAll previous Strahlkorpers and the times they are at. They are expected to have the most recent Strahlkorper in the front and the Strahlkorper furthest in the past in the back of the deque.

◆ ylm_to_stf_0()

Scalar< double > ylm::ylm_to_stf_0 ( const ModalVector l0_coefs)

Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l.

Details

Spherical harmonics of degree l are equivalent to symmetric trace-free tensors of rank l. This equivalence and the transformation is given e.g. in [181], Eqs. (2.10) - (2.14). The conversion coefficients are hard-coded to numerical precision and implemented up to order l=2.

The spherical harmonic coefficients are expected to be sorted with ascending m, i.e. (-m, -m+1, ... , m)

◆ ylm_to_stf_1()

template<typename Frame >
tnsr::i< double, 3, Frame > ylm::ylm_to_stf_1 ( const ModalVector l1_coefs)

Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l.

Details

Spherical harmonics of degree l are equivalent to symmetric trace-free tensors of rank l. This equivalence and the transformation is given e.g. in [181], Eqs. (2.10) - (2.14). The conversion coefficients are hard-coded to numerical precision and implemented up to order l=2.

The spherical harmonic coefficients are expected to be sorted with ascending m, i.e. (-m, -m+1, ... , m)

◆ ylm_to_stf_2()

template<typename Frame >
tnsr::ii< double, 3, Frame > ylm::ylm_to_stf_2 ( const ModalVector l2_coefs)

Converts real spherical harmonic coefficients of degree l into a symmetric trace-free tensor of rank l.

Details

Spherical harmonics of degree l are equivalent to symmetric trace-free tensors of rank l. This equivalence and the transformation is given e.g. in [181], Eqs. (2.10) - (2.14). The conversion coefficients are hard-coded to numerical precision and implemented up to order l=2.

The spherical harmonic coefficients are expected to be sorted with ascending m, i.e. (-m, -m+1, ... , m)