Prototypes | |
static UNUSED REAL8 | XLALSimIMRSpinPrecEOBHamiltonianDeltaR (SpinEOBHCoeffs *coeffs, const REAL8 r, const REAL8 eta, const REAL8 a) |
This function calculates the DeltaR potential function in the spin EOB Hamiltonian. More... | |
static REAL8 | XLALSimIMRSpinPrecEOBHamiltonian (const REAL8 eta, REAL8Vector *restrict x, REAL8Vector *restrict p, REAL8Vector *restrict s1Vec, REAL8Vector *restrict s2Vec, REAL8Vector *restrict sigmaKerr, REAL8Vector *restrict sigmaStar, int tortoise, SpinEOBHCoeffs *coeffs) |
static REAL8 | XLALSimIMRSpinPrecEOBHamiltonianDeltaT (SpinEOBHCoeffs *coeffs, const REAL8 r, const REAL8 eta, const REAL8 a) |
This function calculates the function \(\Delta_t(r)\) which appears in the spinning EOB potential function. More... | |
static REAL8 | inner_product (const REAL8 values1[], const REAL8 values2[]) |
Functions to compute the inner product and cross products between vectors. More... | |
static REAL8 * | cross_product (const REAL8 values1[], const REAL8 values2[], REAL8 result[]) |
static REAL8 * | rotate_vector (const REAL8 v[], const REAL8 k[], const REAL8 theta, REAL8 result[]) |
Rotate the vector v around the axis k by an angle theta, counterclockwise Note that this assumes that the rotation axis k is normalized. More... | |
static UNUSED REAL8 | XLALSimIMRSpinPrecEOBNonKeplerCoeff (const REAL8 values[], SpinEOBParams *funcParams) |
Function to calculate the non-Keplerian coefficient for the PRECESSING EOB model. More... | |
static REAL8 | XLALSimIMRSpinPrecEOBCalcOmega (const REAL8 values[], SpinEOBParams *funcParams) |
Function to calculate the value of omega for the PRECESSING EOB waveform. More... | |
static int | XLALSpinPrecHcapRvecDerivative (double UNUSED t, const REAL8 values[], REAL8 dvalues[], void *funcParams) |
Function to calculate numerical derivatives of the spin EOB Hamiltonian, to get \(dr/dt\), as decribed in Eqs. More... | |
static double | GSLSpinPrecHamiltonianWrapperForRvecDerivs (double x, void *params) |
Wrapper for GSL to call the Hamiltonian function. More... | |
static double | GSLSpinPrecHamiltonianWrapperFordHdpphi (double x, void *params) |
Wrapper for GSL to call the Hamiltonian function. More... | |
static REAL8 | XLALSimIMRSpinPrecEOBHamiltonian (const REAL8 eta, REAL8Vector *restrict x, REAL8Vector *restrict p, REAL8Vector *restrict s1Vec, REAL8Vector *restrict s2Vec, REAL8Vector *restrict sigmaKerr, REAL8Vector *restrict sigmaStar, INT4 tortoise, SpinEOBHCoeffs *coeffs) |
Function to calculate the value of the spinning Hamiltonian for given values of the dynamical variables (in a Cartesian co-ordinate system). More... | |
Go to the source code of this file.
Macros | |
#define | _LALSIMIMRSPINEOBHAMILTONIANPREC_C |
#define _LALSIMIMRSPINEOBHAMILTONIANPREC_C |
Functions for calculating the effective one-body Hamiltonian for spinning binaries, as described in Taracchini et al. ( PRD 86, 024011 (2012), arXiv 1202.0790 ). All equation numbers in this file refer to equations of this paper, unless otherwise specified. This code borrows hugely from a C implementation originally written by Enrico Barausse, following Barausse and Buonanno PRD 81, 084024 (2010) and PRD 84, 104027 (2011), henceforth BB1 and BB2
Definition at line 15 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
This function calculates the DeltaR potential function in the spin EOB Hamiltonian.
This function calculates the function \(\Delta_r(r)\) which appears in the spinning EOB potential function.
Eqs. 5.83 of PRD 81, 084024 (2010)
coeffs | < Pre-computed coefficients which appear in the function |
r | < Current orbital radius (in units of total mass) |
eta | < Symmetric mass ratio |
a | < Normalized deformed Kerr spin |
Definition at line 649 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
|
static |
This function calculates the function \(\Delta_t(r)\) which appears in the spinning EOB potential function.
Eqs. 5.73 of PRD 81, 084024 (2010) augmented by 4PN, linear-in-eta corrections: see also section "New 4PN term in the radial potential" of https://dcc.ligo.org/T1400476
coeffs | < Pre-computed coefficients which appear in the function |
r | < Current orbital radius (in units of total mass) |
eta | < Symmetric mass ratio |
a | < Normalized deformed Kerr spin |
Definition at line 605 of file LALSimIMRSpinEOBHamiltonianPrec.c.
Functions to compute the inner product and cross products between vectors.
Definition at line 684 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Definition at line 694 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Rotate the vector v around the axis k by an angle theta, counterclockwise Note that this assumes that the rotation axis k is normalized.
Definition at line 709 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Function to calculate the non-Keplerian coefficient for the PRECESSING EOB model.
radius \(r\) times the cuberoot of the returned number is \(r_\Omega\) defined in Eq. A2, i.e. the function returns \((r_{\Omega} / r)^3\) = \(1/(r^3 (\partial Hreal/\partial p_\phi |p_r=0)^2)\).
values | < Dynamical variables |
funcParams | < EOB parameters |
Definition at line 988 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Function to calculate the value of omega for the PRECESSING EOB waveform.
Needs the dynamics in Cartesian coordinates.
First, the frame is rotated so that L is along the y-axis. this rotation includes the spins.
Second, \(\vec{r}\) and \(\vec{p}\) are converted to polar coordinates (and not the spins). As L is along the y-axis, \(\theta\) defined as the angle between L and the y-axis is 0, which is a cyclic coordinate now and that fixes \(p_\theta = 0\).
Third, \(p_r\) is set to 0.
Fourth, the polar \((r,\phi,p_r=0,p_\phi)\) and the Cartesian spin vectors are used to compute the derivative \(\partial Hreal/\partial p_\phi |p_r=0\).
the polarvalues, respectively, are \({r, \theta, \phi, p_r, p_\theta, p_\phi}\)
values | < Dynamical variables |
funcParams | < EOB parameters |
Definition at line 738 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Function to calculate numerical derivatives of the spin EOB Hamiltonian, to get \(dr/dt\), as decribed in Eqs.
A4 of PRD 81, 084041 (2010) This function is not used by the spin-aligned SEOBNRv1 model.
t | < UNUSED |
values | < Dynamical variables |
dvalues | < Time derivatives of variables (returned) |
funcParams | < EOB parameters |
Definition at line 1026 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Wrapper for GSL to call the Hamiltonian function.
This is simply the function GSLSpinPrecHamiltonianWrapper copied over. The alternative was to make it non-static which increases runtime as static functions can be better optimized.
Definition at line 1634 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Wrapper for GSL to call the Hamiltonian function.
This is simply the function GSLSpinPrecHamiltonianWrapper copied over. The alternative was to make it non-static which increases runtime as static functions can be better optimized.
Definition at line 1742 of file LALSimIMRSpinEOBHamiltonianPrec.c.
|
static |
Function to calculate the value of the spinning Hamiltonian for given values of the dynamical variables (in a Cartesian co-ordinate system).
The inputs are as follows:
x - the separation vector r expressed in Cartesian co-ordinates p - the momentum vector (with the radial component tortoise pr*) sigmaKerr - spin of the effective Kerr background (a combination of the individual spin vectors) sigmaStar - spin of the effective particle (a different combination of the individual spins). coeffs - coefficients which crop up in the Hamiltonian. These can be calculated using the XLALCalculateSpinEOBParams() function.
The function returns a REAL8, which will be the value of the Hamiltonian if all goes well; otherwise, it will return the XLAL REAL8 failure NaN. The Hamiltonian function is described in PRD 81, 084024 (2010) and PRD 84, 104027 (2011)
eta | < Symmetric mass ratio |
x | < Position vector |
p | < Momentum vector (tortoise radial component pr*) |
s1Vec | < Spin vector 1 |
s2Vec | < Spin vector 2 |
sigmaKerr | < Spin vector sigma_kerr |
sigmaStar | < Spin vector sigma_star |
tortoise | < flag to state whether the momentum is the tortoise co-ord |
coeffs | < Structure containing various coefficients |
Definition at line 126 of file LALSimIMRSpinEOBHamiltonianPrec.c.