SPINspiral
SPINspiral is a parameter-estimation code for gravitational-wave signals detected by LIGO/Virgo
|
Contains interfaces to LAL routines. More...
#include <math.h>
#include <lal/LALStdlib.h>
#include <lal/LALInspiral.h>
#include <lal/GeneratePPNInspiral.h>
#include <lal/GenerateInspiral.h>
#include <lal/DetectorSite.h>
#include <lal/Date.h>
#include <lal/TimeDelay.h>
#include <lal/SkyCoordinates.h>
#include <SPINspiral.h>
#include <stdio.h>
#include <lal/LALStdio.h>
#include <lal/Units.h>
#include <lal/FrameCache.h>
#include <lal/FrameStream.h>
#include <lal/LALDetectors.h>
#include <lal/AVFactories.h>
#include <lal/ResampleTimeSeries.h>
#include <lal/TimeSeries.h>
#include <lal/FrequencySeries.h>
#include <lal/StringInput.h>
#include <lal/VectorOps.h>
#include <lal/Random.h>
#include <lal/XLALError.h>
#include <lal/LIGOLwXMLRead.h>
#include <lal/NRWaveInject.h>
#include <lal/LIGOMetadataUtils.h>
#include <lal/Inject.h>
Functions | |
void | templateLAL12 (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute a 12-parameter (one spin) LAL waveform. | |
void | LALHpHc12 (LALStatus *status, CoherentGW *waveform, SimInspiralTable *injParams, PPNParamStruc *ppnParams, int *l, struct parSet *par, struct interferometer *ifo, int injectionWF, struct runPar run) |
Compute waveform for a 12-parameter (one spin) LAL waveform. | |
void | templateLAL15old (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute waveform for a 15-parameter (two spins) LAL waveform. | |
void | LALHpHc15 (LALStatus *status, CoherentGW *waveform, SimInspiralTable *injParams, PPNParamStruc *ppnParams, int *l, struct parSet *par, struct interferometer *ifo, int injectionWF, struct runPar run) |
Compute a waveform for a 15-parameter (two spins) LAL waveform. | |
void | templateLAL15 (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute waveform for a 15-parameter (two spins) LAL waveform. | |
void | templateLALPhenSpinTaylorRD (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute waveform for a 15-parameter (two spins) LAL PhenSpinTaylorRD waveform. | |
void | templateLALnonSpinning (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute waveform template for a LAL non-spinning inspiral waveform. | |
double | LALFpFc (LALStatus *status, CoherentGW *waveform, SimInspiralTable *injParams, PPNParamStruc *ppnParams, double *wave, int length, struct parSet *par, struct interferometer *ifo, int ifonr) |
Compute detector response for a given detector and given h_+,h_x. | |
void | getWaveformApproximant (const char *familyName, int length, double PNorder, char *waveformApproximant) |
Compose the waveform approximant from the family name and pN order. | |
void | LALfreedomSpin (CoherentGW *waveform) |
Free spinning LAL variables. | |
void | LALfreedomPhenSpinTaylorRD (CoherentGW *waveform) |
Free spinning LAL variables. | |
void | LALfreedomNoSpin (CoherentGW *waveform) |
Free non-spinning LAL variables. |
Contains interfaces to LAL routines.
void templateLAL12 | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute a 12-parameter (one spin) LAL waveform.
Use the LAL <=3.5 PN spinning waveform, with 1 spinning object (12 parameters)
References interferometer::FTin, LALFpFc(), LALfreedomSpin(), LALHpHc12(), and interferometer::samplesize.
Referenced by waveformTemplate().
void LALHpHc12 | ( | LALStatus * | status, |
CoherentGW * | waveform, | ||
SimInspiralTable * | injParams, | ||
PPNParamStruc * | ppnParams, | ||
int * | l, | ||
struct parSet * | par, | ||
struct interferometer * | ifo, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute waveform for a 12-parameter (one spin) LAL waveform.
Compute h_+ and h_x from the parameters in par and interferometer information in ifo. l is a pointer to get the lenght of the waveform computed, this length is also available as waveform->phi->data->length.
References GMST(), runPar::injectionPNorder, runPar::injRevID, longitude(), interferometer::lowCut, runPar::mcmcPNorder, mtpi, parSet::par, runPar::parRevID, interferometer::samplerate, and tpi.
Referenced by templateLAL12().
void templateLAL15old | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute waveform for a 15-parameter (two spins) LAL waveform.
Use the LAL 3.5/2.5 PN spinning waveform, with 2 spinning objects (15 parameters)
References interferometer::FTin, LALFpFc(), LALfreedomSpin(), LALHpHc15(), and interferometer::samplesize.
void LALHpHc15 | ( | LALStatus * | status, |
CoherentGW * | waveform, | ||
SimInspiralTable * | injParams, | ||
PPNParamStruc * | ppnParams, | ||
int * | l, | ||
struct parSet * | par, | ||
struct interferometer * | ifo, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute a waveform for a 15-parameter (two spins) LAL waveform.
Compute h_+ and h_x form the parameters in par and interferometer information in ifo. l is a pointer to get the lenght of the waveform computed, this length is also available in waveform->phi->data->length.
References getWaveformApproximant(), GMST(), interferometer::highCut, runPar::injectionPNorder, runPar::injRevID, longitude(), interferometer::lowCut, McEta2masses(), runPar::mcmcPNorder, mtpi, parSet::par, runPar::parRevID, interferometer::samplerate, and tpi.
Referenced by templateLAL15old().
void templateLAL15 | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute waveform for a 15-parameter (two spins) LAL waveform.
Use the LAL 3.5/2.5 PN spinning waveform, with 2 spinning objects (15 parameters)
References getWaveformApproximant(), GMST(), runPar::injectionPNorder, runPar::injRevID, LALFpFc(), LALfreedomSpin(), longitude(), McEta2masses(), runPar::mcmcPNorder, mtpi, parSet::par, runPar::parRevID, interferometer::samplesize, and tpi.
Referenced by waveformTemplate().
void templateLALPhenSpinTaylorRD | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute waveform for a 15-parameter (two spins) LAL PhenSpinTaylorRD waveform.
Use the LAL 3.5/2.5 PN spinning hybrid PhenSpinTaylorRD waveform, with 2 spinning objects (15 parameters). Phenomenological ring-down attached
References interferometer::FTin, interferometer::FTstart, getWaveformApproximant(), GMST(), runPar::injectionPNorder, runPar::injRevID, LALfreedomPhenSpinTaylorRD(), longitude(), McEta2masses(), runPar::mcmcPNorder, mtpi, parSet::par, runPar::parRevID, interferometer::samplesize, and tpi.
Referenced by waveformTemplate().
void templateLALnonSpinning | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute waveform template for a LAL non-spinning inspiral waveform.
Uses GeneratePPN approximant.
References c3rd, getWaveformApproximant(), GMST(), runPar::injectionPNorder, runPar::injRevID, LALFpFc(), LALfreedomNoSpin(), longitude(), interferometer::lowCut, McEta2masses(), runPar::mcmcPNorder, mtpi, parSet::par, runPar::parRevID, interferometer::samplesize, and tpi.
Referenced by waveformTemplate().
double LALFpFc | ( | LALStatus * | status, |
CoherentGW * | waveform, | ||
SimInspiralTable * | injParams, | ||
PPNParamStruc * | ppnParams, | ||
double * | wave, | ||
int | length, | ||
struct parSet * | par, | ||
struct interferometer * | ifo, | ||
int | ifonr | ||
) |
Compute detector response for a given detector and given h_+,h_x.
Compute the detector response for a given detector (ifonr) and h_+,h_x. the structure waveform must already hold the computed values of h+,x (or just a1, a2, phi and shift as a function of time)
References interferometer::FTstart, and parSet::par.
Referenced by templateLAL12(), templateLAL15(), templateLAL15old(), and templateLALnonSpinning().
void getWaveformApproximant | ( | const char * | familyName, |
int | length, | ||
double | PNorder, | ||
char * | waveformApproximant | ||
) |
Compose the waveform approximant from the family name and pN order.
Referenced by LALHpHc15(), templateLAL15(), templateLALnonSpinning(), and templateLALPhenSpinTaylorRD().
void LALfreedomSpin | ( | CoherentGW * | waveform | ) |
Free spinning LAL variables.
Referenced by templateLAL12(), templateLAL15(), and templateLAL15old().
void LALfreedomPhenSpinTaylorRD | ( | CoherentGW * | waveform | ) |
Free spinning LAL variables.
Referenced by templateLALPhenSpinTaylorRD().
void LALfreedomNoSpin | ( | CoherentGW * | waveform | ) |
Free non-spinning LAL variables.
Referenced by templateLALnonSpinning().