SPINspiral
SPINspiral is a parameter-estimation code for gravitational-wave signals detected by LIGO/Virgo
Functions
src/SPINspiral_routines.c File Reference

(SPINspiral version 0.5.1, revision 275)

Contains various routines. More...

#include <SPINspiral.h>
Include dependency graph for SPINspiral_routines.c:

Functions

double massRatio (double m1, double m2)
 Compute symmetric mass ratio (eta) for given individual masses.
double chirpMass (double m1, double m2)
 Compute chirp mass (mM) for given individual masses.
void McEta2masses (double Mc, double eta, double *m1, double *m2)
 Compute individual masses (m1,m2) for given chirp mass (Mc) and symmetric mass ratio (eta)
void masses2McEta (double m1, double m2, double *Mc, double *eta)
 Compute chirp mass (Mc) and symmetric mass ratio (eta) for given individual masses (m1,m2)
double GMST (double GPSsec)
 Compute Greenwich Mean Sideral time from GPS time.
double rightAscension (double longi, double gmst)
 Compute RA from 'longitude' and GMST.
double longitude (double ra, double gmst)
 Compute 'longitude' from RA and GMST.
double dotProduct (double vec1[3], double vec2[3])
 Compute the dot product of two vectors.
void facVec (double vec1[3], double fac, double vec2[3])
 Multiply a vector by a factor (scalar): vec2 = fac*vec1.
void addVec (double vec1[3], double vec2[3], double result[3])
 Add two vectors result = vec1 + vec2.
void normalise (double vec[3])
 Normalise a vector: nvec = vec / |vec|.
void crossProduct (double vec1[3], double vec2[3], double result[3])
 Compute the cross product of two vectors.
void rotate (double x[3], double ang, double axis[3])
 Rotate the vextor x about angle around the normal vector axis.
int rightHanded (double x[3], double y[3], double z[3])
 Determine whether the vectors x, y and z constitute a right-handed system.
void orthoProject (double x[3], double vec1[3], double vec2[3])
 Determines the orthogonal projection of vector x onto the span of the two ORTHONORMAL vectors vec1 and vec2.
double angle (double x[3], double y[3])
 Determines the angle between vectors x and y.
void coord2vec (double sinlati, double longi, double x[3])
 Convert geographical spherical coordinates to a Cartesian normal vector.
void vec2coord (double x[3], double *sinlati, double *longi)
 Compute geographical spherical coordinates from a Cartesian normal vector.
void setSeed (int *seed)
 Grab a random seed from the system clock.

Detailed Description

Contains various routines.


Function Documentation

double massRatio ( double  m1,
double  m2 
)

Compute symmetric mass ratio (eta) for given individual masses.

double chirpMass ( double  m1,
double  m2 
)

Compute chirp mass (mM) for given individual masses.

void McEta2masses ( double  Mc,
double  eta,
double *  m1,
double *  m2 
)

Compute individual masses (m1,m2) for given chirp mass (Mc) and symmetric mass ratio (eta)

Referenced by dataFT(), LALHpHc15(), templateApostolatos(), templateLAL15(), templateLALnonSpinning(), and templateLALPhenSpinTaylorRD().

void masses2McEta ( double  m1,
double  m2,
double *  Mc,
double *  eta 
)

Compute chirp mass (Mc) and symmetric mass ratio (eta) for given individual masses (m1,m2)

double GMST ( double  GPSsec)

Compute Greenwich Mean Sideral time from GPS time.

Computes the Greenwich Mean Sidereal Time (in radians!) from GPS time (in seconds). See K.R. Lang (1999), p.80sqq.

References pi.

Referenced by dataFT(), LALHpHc12(), LALHpHc15(), localPar(), templateApostolatos(), templateLAL15(), templateLALnonSpinning(), and templateLALPhenSpinTaylorRD().

double rightAscension ( double  longi,
double  gmst 
)

Compute RA from 'longitude' and GMST.

All quantities are in radians

References mtpi, and tpi.

double longitude ( double  ra,
double  gmst 
)

Compute 'longitude' from RA and GMST.

Computes the 'longitude' from the right ascension and GMST. All quantities are in radians. In fact, 'longitude' is something like the Greenwich hour angle of the corresponding RA.

References mtpi, and tpi.

Referenced by LALHpHc12(), LALHpHc15(), localPar(), templateApostolatos(), templateLAL15(), templateLALnonSpinning(), and templateLALPhenSpinTaylorRD().

double dotProduct ( double  vec1[3],
double  vec2[3] 
)

Compute the dot product of two vectors.

Referenced by templateApostolatos().

void facVec ( double  vec1[3],
double  fac,
double  vec2[3] 
)

Multiply a vector by a factor (scalar): vec2 = fac*vec1.

Referenced by templateApostolatos().

void addVec ( double  vec1[3],
double  vec2[3],
double  result[3] 
)

Add two vectors result = vec1 + vec2.

Referenced by templateApostolatos().

void normalise ( double  vec[3])

Normalise a vector: nvec = vec / |vec|.

void crossProduct ( double  vec1[3],
double  vec2[3],
double  result[3] 
)

Compute the cross product of two vectors.

Referenced by templateApostolatos().

void rotate ( double  x[3],
double  ang,
double  axis[3] 
)

Rotate the vextor x about angle around the normal vector axis.

Rotates vector x clockwise around axis (looking along axis while it is pointing towards you). axis must be a UNIT VECTOR

Referenced by IFOinit().

int rightHanded ( double  x[3],
double  y[3],
double  z[3] 
)

Determine whether the vectors x, y and z constitute a right-handed system.

Determines whether vectors x,y & z constitute a right-handed system by checking the sign of the triple product or det(x,y,z).

Referenced by localPar(), and vec2coord().

void orthoProject ( double  x[3],
double  vec1[3],
double  vec2[3] 
)

Determines the orthogonal projection of vector x onto the span of the two ORTHONORMAL vectors vec1 and vec2.

Referenced by localPar(), and vec2coord().

double angle ( double  x[3],
double  y[3] 
)

Determines the angle between vectors x and y.

Referenced by localPar(), and vec2coord().

void coord2vec ( double  sinlati,
double  longi,
double  x[3] 
)

Convert geographical spherical coordinates to a Cartesian normal vector.

Turns geographical coordinates (latitude & longitude) into a vector - Result is a unit (!) vector referring to the (right-handed) coordinate system spanned by the three vectors pointing from geocenter to: x) intersection of greenwich meridian and equator y) intersection of 90 deg. East meridian and equator z) north pole

Referenced by IFOinit(), and localPar().

void vec2coord ( double  x[3],
double *  sinlati,
double *  longi 
)

Compute geographical spherical coordinates from a Cartesian normal vector.

Inverse of coord2vec() (see there for more details)

References angle(), orthoProject(), pi, and rightHanded().

void setSeed ( int *  seed)

Grab a random seed from the system clock.

If seed==0, set (randomise) seed using the system clock (hence a random (random seed) :-)

Referenced by main(), and setRandomInjectionParameters().

 All Data Structures Files Functions Variables Defines