LALSimulation  5.4.0.1-fe68b98
LALSimIMRSpinEOBInitialConditionsPrec.c File Reference

Prototypes

static REAL8 CalculateDotProductPrec (const REAL8 a[], const REAL8 b[])
 Calculate the dot product of two vectors. More...
 
static REAL8 CalculateCrossProductPrec (const INT4 i, const REAL8 a[], const REAL8 b[])
 Calculates the ith component of the cross product of two vectors, where i is in the range 0-2. More...
 
static int NormalizeVectorPrec (REAL8 a[])
 Normalizes the given vector. More...
 
static int CalculateRotationMatrixPrec (gsl_matrix *rotMatrix, gsl_matrix *rotInverse, REAL8 r[], REAL8 v[], REAL8 L[])
 Calculate the rotation matrix and its inverse. More...
 
static int ApplyRotationMatrixPrec (gsl_matrix *rotMatrix, REAL8 a[])
 Applies a rotation matrix to a given vector. More...
 
static int SphericalToCartesianPrec (REAL8 qCart[], REAL8 pCart[], const REAL8 qSph[], const REAL8 pSph[])
 Performs a co-ordinate transformation from spherical to Cartesian co-ordinates. More...
 
static int CartesianToSphericalPrec (REAL8 qSph[], REAL8 pSph[], const REAL8 qCart[], const REAL8 pCart[])
 Perform a co-ordinate transformation from Cartesian to spherical co-ordinates. More...
 
static int XLALFindSphericalOrbitPrec (const gsl_vector *x, void *params, gsl_vector *f)
 Root function for gsl root finders. More...
 
static double GSLSpinHamiltonianDerivWrapperPrec (double x, void *params)
 Wrapper for calculating specific derivative of the Hamiltonian in spherical co-ordinates, dH/dr, dH/dptheta and dH/dpphi. More...
 
static int XLALRobustDerivative (const gsl_function *F, double x, double h, double *result, double *absErr)
 
static REAL8 XLALCalculateSphHamiltonianDeriv2Prec (const int idx1, const int idx2, const REAL8 values[], SpinEOBParams *params, INT4 use_optimized)
 Function to calculate the second derivative of the Hamiltonian. More...
 
static int XLALSimIMRSpinEOBInitialConditionsPrec (REAL8Vector *initConds, const REAL8 mass1, const REAL8 mass2, const REAL8 fMin, const REAL8 inc, const REAL8 spin1[], const REAL8 spin2[], SpinEOBParams *params, INT4 use_optimized)
 Main function for calculating the spinning EOB initial conditions, following the quasi-spherical initial conditions described in Sec. More...
 

Go to the source code of this file.

Macros

#define _LALSIMIMRSPINPRECEOBINITIALCONDITIONS_C
 

Variables

static REAL8 scale1 = 1
 
static REAL8 scale2 = 2
 
static REAL8 scale3 = 200
 

Macro Definition Documentation

◆ _LALSIMIMRSPINPRECEOBINITIALCONDITIONS_C

#define _LALSIMIMRSPINPRECEOBINITIALCONDITIONS_C

Definition at line 21 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

Function Documentation

◆ CalculateDotProductPrec()

static REAL8 CalculateDotProductPrec ( const REAL8  a[],
const REAL8  b[] 
)
inlinestatic

Calculate the dot product of two vectors.

Parameters
a< first vector
b< second vector

Definition at line 34 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ CalculateCrossProductPrec()

static REAL8 CalculateCrossProductPrec ( const INT4  i,
const REAL8  a[],
const REAL8  b[] 
)
inlinestatic

Calculates the ith component of the cross product of two vectors, where i is in the range 0-2.

Parameters
i< i-th component of cross product
a< first vector
b< second vector

Definition at line 48 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ NormalizeVectorPrec()

static int NormalizeVectorPrec ( REAL8  a[])
inlinestatic

Normalizes the given vector.

Definition at line 62 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ CalculateRotationMatrixPrec()

static int CalculateRotationMatrixPrec ( gsl_matrix *  rotMatrix,
gsl_matrix *  rotInverse,
REAL8  r[],
REAL8  v[],
REAL8  L[] 
)
static

Calculate the rotation matrix and its inverse.

Rotate the ex-ey-ez frame to the r-v-L frame. This static function is called only by XLALSimIMRSpinEOBInitialConditionsPrec

a, b, g are the angles alpha, beta and gamma

Implement the Rotation Matrix following the "x-convention"

  1. rotate about the z-axis by an angle a, rotation matrix Rz(a);
  2. rotate about the former x-axis (now x') by an angle b, rotation matrix Rx(b);
  3. rotate about the former z-axis (now z') by an algle g, rotation matrix Rz(g);
Parameters
rotMatrixOUTPUT, rotation matrix
rotInverseOUTPUT, rotation matrix inversed
rposition vector
vvelocity vector
Lorbital angular momentum

Definition at line 79 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ ApplyRotationMatrixPrec()

static int ApplyRotationMatrixPrec ( gsl_matrix *  rotMatrix,
REAL8  a[] 
)
static

Applies a rotation matrix to a given vector.

Parameters
rotMatrixrotation matrix
aOUTPUT, vector rotated

Definition at line 178 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ SphericalToCartesianPrec()

static int SphericalToCartesianPrec ( REAL8  qCart[],
REAL8  pCart[],
const REAL8  qSph[],
const REAL8  pSph[] 
)
static

Performs a co-ordinate transformation from spherical to Cartesian co-ordinates.

In the code from Tyson Littenberg, this was only applicable to the special theta=pi/2, phi=0 case. This special transformation is a static function called only by GSLSpinHamiltonianDerivWrapperPrec and XLALSimIMRSpinEOBInitialConditionsPrec

Parameters
qCart< OUTPUT, position vector in Cartesean coordinates
pCart< OUTPUT, momentum vector in Cartesean coordinates
qSph< position vector in spherical coordinates
pSph< momentum vector in spherical coordinates

Definition at line 204 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ CartesianToSphericalPrec()

static int CartesianToSphericalPrec ( REAL8  qSph[],
REAL8  pSph[],
const REAL8  qCart[],
const REAL8  pCart[] 
)
static

Perform a co-ordinate transformation from Cartesian to spherical co-ordinates.

In the code from Tyson Littenberg, this was only applicable to the special theta=pi/2, phi=0 case. This special transformation is a static function called only by XLALSimIMRSpinEOBInitialConditionsPrec

Parameters
qSph< OUTPUT, position vector in spherical coordinates
pSph< OUTPUT, momentum vector in Cartesean coordinates
qCart< position vector in spherical coordinates
pCart< momentum vector in Cartesean coordinates

Definition at line 249 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ XLALFindSphericalOrbitPrec()

static int XLALFindSphericalOrbitPrec ( const gsl_vector *  x,
void *  params,
gsl_vector *  f 
)
static

Root function for gsl root finders.

The gsl root finder will look for gsl_vector *x position in parameter space where the returned values in gsl_vector *f are zero. The functions on which we search for roots are: dH/dr, dH/dptheta and dH/dpphi-omega, namely, Eqs. 4.8 and 4.9 of Buonanno et al. PRD 74, 104005 (2006).

Parameters
x< Parameters requested by gsl root finder
params< Spin EOB parameters
f< Function values for the given parameters

Definition at line 297 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ GSLSpinHamiltonianDerivWrapperPrec()

static double GSLSpinHamiltonianDerivWrapperPrec ( double  x,
void *  params 
)
static

Wrapper for calculating specific derivative of the Hamiltonian in spherical co-ordinates, dH/dr, dH/dptheta and dH/dpphi.

It only works for the specific co-ord system we use here

Parameters
x< Derivative at x
params< Function parameters

Definition at line 449 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ XLALRobustDerivative()

static int XLALRobustDerivative ( const gsl_function *  F,
double  x,
double  h,
double *  result,
double *  absErr 
)
static

Definition at line 565 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ XLALCalculateSphHamiltonianDeriv2Prec()

static REAL8 XLALCalculateSphHamiltonianDeriv2Prec ( const int  idx1,
const int  idx2,
const REAL8  values[],
SpinEOBParams params,
INT4  use_optimized 
)
static

Function to calculate the second derivative of the Hamiltonian.

The derivatives are taken with respect to indices idx1, idx2

Parameters
idx1< Derivative w.r.t. index 1
idx2< Derivative w.r.t. index 2
values< Dynamical variables in spherical coordinates
params< Spin EOB Parameters
use_optimized< use_optimized=1 -> use EXACT instead of finite difference derivatives... UNSUPPORTED AT MOMENT.

Definition at line 610 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ XLALSimIMRSpinEOBInitialConditionsPrec()

static int XLALSimIMRSpinEOBInitialConditionsPrec ( REAL8Vector initConds,
const REAL8  mass1,
const REAL8  mass2,
const REAL8  fMin,
const REAL8  inc,
const REAL8  spin1[],
const REAL8  spin2[],
SpinEOBParams params,
INT4  use_optimized 
)
static

Main function for calculating the spinning EOB initial conditions, following the quasi-spherical initial conditions described in Sec.

IV A of Buonanno, Chen & Damour PRD 74, 104005 (2006). All equation numbers in the comments of this function refer to this paper.

Inputs are binary parameters (masses, spin vectors and inclination), EOB model parameters and initial frequency.

Outputs are initial dynamical variables in a vector (x, y, z, px, py, pz, S1x, S1y, S1z, S2x, S2y, S2z).

In the paper, the initial conditions are solved for a given initial radius, while in this function, they are solved for a given inital orbital frequency.

This difference is reflected in solving Eq. (4.8). The calculation is done in 5 steps: STEP 1) Rotate to LNhat0 along z-axis and N0 along x-axis frame, where LNhat0 and N0 are initial normal to orbital plane and initial orbital separation; STEP 2) After rotation in STEP 1, in spherical coordinates, phi0 and theta0 are given directly in Eq. (4.7), r0, pr0, ptheta0 and pphi0 are obtained by solving Eqs. (4.8) and (4.9) (using gsl_multiroot_fsolver). At this step, we find initial conditions for a spherical orbit without radiation reaction. STEP 3) Rotate to L0 along z-axis and N0 along x-axis frame, where L0 is the initial orbital angular momentum and L0 is calculated using initial position and linear momentum obtained in STEP 2. STEP 4) In the L0-N0 frame, we calculate (dE/dr)|sph using Eq. (4.14), then initial dr/dt using Eq. (4.10), and finally pr0 using Eq. (4.15). STEP 5) Rotate back to the original inertial frame by inverting the rotation of STEP 3 and then inverting the rotation of STEP 1.

Parameters
initConds< OUTPUT, Initial dynamical variables
mass1< mass 1
mass2< mass 2
fMin< Initial frequency (given)
inc< Inclination
spin1< Initial spin vector 1
spin2< Initial spin vector 2
params< Spin EOB parameters
use_optimized< use_optimized=1 -> use EXACT instead of finite difference derivatives, where needed

Definition at line 705 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

Variable Documentation

◆ scale1

REAL8 scale1 = 1
static

Definition at line 27 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ scale2

REAL8 scale2 = 2
static

Definition at line 27 of file LALSimIMRSpinEOBInitialConditionsPrec.c.

◆ scale3

REAL8 scale3 = 200
static

Definition at line 27 of file LALSimIMRSpinEOBInitialConditionsPrec.c.