SPINspiral
SPINspiral is a parameter-estimation code for gravitational-wave signals detected by LIGO/Virgo
|
Main header file. More...
#include <stdio.h>
#include <math.h>
#include <complex.h>
#include <fftw3.h>
#include <time.h>
#include <remez.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <sys/time.h>
#include <stdlib.h>
#include <lal/LALStdlib.h>
#include <lal/LALInspiral.h>
#include <lal/GeneratePPNInspiral.h>
#include <lal/GenerateInspiral.h>
#include <lal/LALFrameL.h>
#include <lal/LALConstants.h>
Go to the source code of this file.
Data Structures | |
struct | runPar |
struct | MCMCvariables |
struct | parSet |
struct | interferometer |
Defines | |
#define | TRUE (1==1) |
#define | FALSE (!TRUE) |
#define | h2r 0.2617993877991494263 |
#define | r2h 3.819718634205488207 |
#define | d2r 0.01745329251994329509 |
#define | r2d 57.29577951308232311 |
#define | c3rd 0.3333333333333333333 |
#define | M0 4.926e-6 |
#define | max(A, B) ((A)>(B)?(A):(B)) |
#define | min(A, B) ((A)<(B)?(A):(B)) |
#define | USAGE "\n\n\ Usage: lalapps_spinspiral \n\ -i <main input file> \n\ --injXMLfile <injection XML file name> --injXMLnr <injection number in XML file (0-...)> \n\ --mChirp <chirp mass in solar mass> \n\ --eta <eta between 0.03 and 0.25> \n\ --tc <t_c in GPS time, e.g. : 873731234.00000> \n\ --dist <distance in Mpc> \n\ --nIter <number of iterations> \n\ --nSkip <number of skipped steps between stored steps> \n\ --rseed <6 number random seed, e.g. 12345. Use 0 for random random seed ...> \n\ --network <network configuration, e.g. H1=[1], H1L1=[1,2], H1L1V1=[1,2,3], V1H1=[3,1]> \n\ --downsample <downsample factor> \n\ --beforetc <data before t_c being analysed in seconds> \n\ --aftertc <data after t_c being analysed in seconds> \n\ --Flow <low frequency cut off in Hz> \n\ --Fhigh <high frequency cut off in Hz> \n\ --tukey1 <alpha1, the 1st parameter of the template tukey window. default 0.0> \n\ --tukey1 <alpha2, the 2nd parameter of the template tukey window. default 0.0> \n\ --nPSDsegment <number of segments to estimate the PSD> \n\ --lPSDsegment <length of each segment to estimate the PSD> \n\ --outputPath <output path> \n\ --cache <list of cache files, e.g. [cache1,cache2]> \n\ --channel <list of channels, e.g. [H1:LSC-STRAIN,L1:LSC-STRAIN]> \n\ --PSDstart <GPS time for the start of the PSD. default begining of cache file> \n\ --template <waveform template, 3 for spin, 4 for no spin> \n\example: ./lalapps_spinspiral -i ./pipeline/SPINspiral.input --mChirp 1.7 --eta 0.12 --tc 873739311.00000 --dist 15 --nIter 5 --nSkip 1 --downsample 1 --beforetc 5 --aftertc 2 --Flow 45 --Fhigh 1600.0 --nPSDsegment 32 --lPSDsegment 4 --network [1,2] --outputPath ./pipeline/ --cache [./pipeline/H-H1_RDS_C03_L2-873739103-873740831.cache,./pipeline/L-L1_RDS_C03_L2-873739055-873740847.cache] --channel [H1:LSC-STRAIN,L1:LSC-STRAIN] --PSDstart 873740000\n\\n\n" |
Functions | |
void | readCommandLineOptions (int argc, char *argv[], struct runPar *run) |
Read and parse command-line arguments/options. | |
void | readMainInputfile (struct runPar *run) |
Read the main input file. | |
void | readMCMCinputfile (struct runPar *run) |
Read the MCMC input file. | |
void | readDataInputfile (struct runPar *run, struct interferometer ifo[]) |
Read the data input file. | |
void | readInjectionInputfile (struct runPar *run) |
Read the injection input file. | |
void | readParameterInputfile (struct runPar *run) |
Read the (MCMC) parameter input file. | |
void | readSystemInputfile (struct runPar *run) |
Read the input file for system (system-dependent) variables, e.g. SPINspiral.input.system. | |
void | readInjectionXML (struct runPar *run) |
Read an XML injection file if specified. | |
void | setParameterNames (struct runPar *run) |
Set parameters names in the hardcoded parameter database. | |
void | setConstants (void) |
Set global variables. | |
void | setIFOdata (struct runPar *run, struct interferometer ifo[]) |
Set all the data for all IFOs that may be used. | |
void | setRandomInjectionParameters (struct runPar *run) |
Get random values for the injection parameters. | |
void | getInjectionParameters (struct parSet *par, int nInjectionPar, double *parInjectVal) |
Returns the 'injection values' to the parameter set par. | |
void | getStartParameters (struct parSet *par, struct runPar run) |
Returns the 'best-guess values' to the parameter set par. | |
void | startMCMCOffset (struct parSet *par, struct MCMCvariables *mcmc, struct interferometer *ifo[], struct runPar run) |
Choose and print offset starting values for the Markov chain. | |
void | setTemperatureLadder (struct MCMCvariables *mcmc) |
Set up the temperature ladder for parallel tempering. | |
void | setTemperatureLadderOld (struct MCMCvariables *mcmc) |
Set up the temperature ladder for parallel tempering - old routine. | |
void | allocParset (struct parSet *par, int networkSize) |
Allocate memory for the vectors in the struct parSet. | |
void | freeParset (struct parSet *par) |
Deallocate memory for the vectors in the struct parSet. | |
void | copyRun2MCMC (struct runPar run, struct MCMCvariables *mcmc) |
Copy some of the elements of the struct runPar to the struct MCMCvariables. | |
void | setMCMCseed (struct runPar *run) |
void | setSeed (int *seed) |
Grab a random seed from the system clock. | |
void | MCMC (struct runPar run, struct interferometer *ifo[]) |
MCMC routine - forms the MCMC core of the program. | |
void | CholeskyDecompose (double **A, struct MCMCvariables *mcmc) |
Compute Cholesky decomposition matrix. | |
void | par2arr (struct parSet par, double **param, struct MCMCvariables mcmc) |
Put the MCMC parameters in an array. | |
void | arr2par (double **param, struct parSet *par, struct MCMCvariables mcmc) |
Get the MCMC parameters out of their array. | |
double | prior (double *par, int p, struct MCMCvariables mcmc) |
Compute the prior for the given parameter set. | |
double | sigmaPeriodicBoundaries (double sigma, int p, struct MCMCvariables mcmc) |
Bring the adaptation sigma between its periodic boundaries. | |
void | correlatedMCMCupdate (struct interferometer *ifo[], struct parSet *state, struct MCMCvariables *mcmc, struct runPar run) |
Do a correlated block MCMC update. | |
void | uncorrelatedMCMCsingleUpdate (struct interferometer *ifo[], struct parSet *state, struct MCMCvariables *mcmc, struct runPar run) |
Do an uncorrelated, per-parameter MCMC update. | |
void | uncorrelatedMCMCblockUpdate (struct interferometer *ifo[], struct parSet *state, struct MCMCvariables *mcmc, struct runPar run) |
Do an uncorrelated block update. | |
void | writeMCMCheader (struct interferometer *ifo[], struct MCMCvariables mcmc, struct runPar run) |
Write MCMC header to screen and file. | |
void | writeMCMCoutput (struct MCMCvariables mcmc, struct interferometer *ifo[]) |
Write an MCMC iteration as an output line to screen and/or file. | |
void | allocateMCMCvariables (struct MCMCvariables *mcmc) |
Allocate memory for the MCMCvariables struct. | |
void | freeMCMCvariables (struct MCMCvariables *mcmc) |
Deallocate memory for the MCMCvariables struct. | |
void | updateCovarianceMatrix (struct MCMCvariables *mcmc) |
Calculate the new covariance matrix and determine whether the new matrix should be accepted. | |
double | annealTemperature (double temp0, int nburn, int nburn0, int iIter) |
Annealing: set the temperature according to the iteration number and burnin. | |
void | swapChains (struct MCMCvariables *mcmc) |
Parallel tempering: Swap states between two chains. | |
void | writeChainInfo (struct MCMCvariables mcmc) |
Parallel tempering: Print chain and swap info to screen. | |
double | chirpMass (double m1, double m2) |
Compute chirp mass (mM) for given individual masses. | |
double | massRatio (double m1, double m2) |
Compute symmetric mass ratio (eta) 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 rightAscension, 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 angle, 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 | IFOinit (struct interferometer **ifo, int networkSize, struct runPar run) |
Set up IFOs. | |
void | IFOdispose (struct interferometer *ifo, struct runPar run) |
Free allocated IFO variables. | |
double * | filter (int *order, int samplerate, double upperlimit, struct runPar run) |
Compute FIR filter coefficients. | |
double * | downsample (double data[], int *datalength, double coef[], int ncoef, struct runPar run) |
Downsample a time series by a factor downsampleFactor. | |
void | dataFT (struct interferometer *ifo[], int i, int networkSize, struct runPar run) |
Read the data and do a software injection if wanted. | |
double | hannWindow (int j, int N) |
Apply a 'Hann window' to data. | |
double | tukeyWindow (int j, int N, double r) |
Apply a 'Tukey window' to data. | |
double | modifiedTukeyWindow (int j, int N, double r1, double r2) |
Apply a 'Tukey window' to data. | |
void | noisePSDestimate (struct interferometer *ifo[], int ifonr, struct runPar run) |
Returns a (smoothed) estimate of the log- Power Spectral Density. | |
double | logNoisePSD (double f, struct interferometer *ifo) |
double | interpolLogNoisePSD (double f, struct interferometer *ifo) |
Returns a linearly interpolated (log-) noise PSD. | |
void | writeDataToFiles (struct interferometer *ifo[], int networkSize, struct runPar run) |
Write windowed, time-domain data (signal + noise) to disc. | |
void | writeNoiseToFiles (struct interferometer *ifo[], int networkSize, struct runPar run) |
Write noise ASD to disc. | |
void | writeSignalsToFiles (struct interferometer *ifo[], int networkSize, struct runPar run) |
Write signal and its FFT to disc. | |
void | printParameterHeaderToFile (FILE *dump, struct interferometer *ifo, struct runPar run) |
Write an optional header when saving data, ASD or signal to disc. | |
void | waveformTemplate (struct parSet *par, struct interferometer *ifo[], int ifonr, int waveformVersion, int injectionWF, struct runPar run) |
Compute an inspiral waveform. | |
void | templateApostolatos (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute an Apostolatos, 12-parameter, single spin, simple-precession waveform. | |
void | localPar (struct parSet *par, struct interferometer *ifo[], int networkSize, int injectionWF, struct runPar run) |
Calculate the local parameters from the global parameters. | |
void | templateLAL12 (struct parSet *par, struct interferometer *ifo[], int ifonr, int injectionWF, struct runPar run) |
Compute 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 | 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. | |
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 | 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 | LALHpHcNonSpinning (LALStatus *status, CoherentGW *waveform, SimInspiralTable *injParams, PPNParamStruc *ppnParams, int *l, struct parSet *par, struct interferometer *ifo, int injectionWF, struct runPar run) |
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 | LALfreedomNoSpin (CoherentGW *waveform) |
Free non-spinning LAL variables. | |
void | LALfreedomPhenSpinTaylorRD (CoherentGW *waveform) |
Free spinning LAL variables. | |
double | netLogLikelihood (struct parSet *par, int networkSize, struct interferometer *ifo[], int waveformVersion, int injectionWF, struct runPar run) |
Compute the log(Likelihood) for a network of IFOs. | |
double | IFOlogLikelihood (struct parSet *par, struct interferometer *ifo[], int i, int waveformVersion, int injectionWF, struct runPar run) |
Compute the log(Likelihood) for a single IFO. | |
double | logLikelihood_nine (struct parSet *par, int waveformVersion, int injectionWF, struct runPar run) |
Compute a analytical log(Likelihood) | |
double | signalToNoiseRatio (struct parSet *par, struct interferometer *ifo[], int i, int waveformVersion, int injectionWF, struct runPar run) |
Compute the SNR of the waveform with a given parameter set for a single IFO. | |
double | parMatch (struct parSet *par1, int waveformVersion1, int injectionWF1, struct parSet *par2, int waveformVersion2, int injectionWF2, struct interferometer *ifo[], int networkSize, struct runPar run) |
Compute match between waveforms with parameter sets par1 and par2. | |
double | overlapWithData (struct parSet *par, struct interferometer *ifo[], int ifonr, int waveformVersion, int injectionWF, struct runPar run) |
Compute frequency-domain overlap of waveform of given parameters with raw data. | |
double | parOverlap (struct parSet *par1, int waveformVersion1, int injectionWF1, struct parSet *par2, int waveformVersion2, int injectionWF2, struct interferometer *ifo[], int ifonr, struct runPar run) |
Compute the overlap in the frequency domain between two waveforms with parameter sets par1 and par2. | |
double | vecOverlap (fftw_complex *vec1, fftw_complex *vec2, double *noise, int j_1, int j_2, double deltaFT) |
Compute the overlap in the frequency domain between two waveforms with parameter sets par1 and par2. | |
void | signalFFT (fftw_complex *FFTout, struct parSet *par, struct interferometer *ifo[], int ifonr, int waveformVersion, int injectionWF, struct runPar run) |
Compute the FFT of a waveform with given parameter set. | |
double | matchBetweenParameterArrayAndTrueParameters (double *pararray, struct interferometer *ifo[], struct MCMCvariables mcmc, struct runPar run) |
Compute the match between a given parameter array and the injection parameters. | |
void | parseCharacterOptionString (char *input, char **strings[], int *n) |
parse strings from "[one,two,three]" to {"one", "two", "three"} | |
void | readCachefile (struct runPar *run, int ifonr) |
Read a Cache file. Returns an array of what is in the cache file. | |
Variables | |
double | Ms |
double | Mpc |
double | G |
double | c |
double | Mpcs |
double | pi |
double | tpi |
double | mtpi |
Main header file.
#define TRUE (1==1) |
#define h2r 0.2617993877991494263 |
#define r2h 3.819718634205488207 |
#define d2r 0.01745329251994329509 |
#define c3rd 0.3333333333333333333 |
Referenced by readParameterInputfile(), templateApostolatos(), and templateLALnonSpinning().
#define M0 4.926e-6 |
Referenced by templateApostolatos().
#define max | ( | A, | |
B | |||
) | ((A)>(B)?(A):(B)) |
#define min | ( | A, | |
B | |||
) | ((A)<(B)?(A):(B)) |
#define USAGE "\n\n\ Usage: lalapps_spinspiral \n\ -i <main input file> \n\ --injXMLfile <injection XML file name> --injXMLnr <injection number in XML file (0-...)> \n\ --mChirp <chirp mass in solar mass> \n\ --eta <eta between 0.03 and 0.25> \n\ --tc <t_c in GPS time, e.g. : 873731234.00000> \n\ --dist <distance in Mpc> \n\ --nIter <number of iterations> \n\ --nSkip <number of skipped steps between stored steps> \n\ --rseed <6 number random seed, e.g. 12345. Use 0 for random random seed ...> \n\ --network <network configuration, e.g. H1=[1], H1L1=[1,2], H1L1V1=[1,2,3], V1H1=[3,1]> \n\ --downsample <downsample factor> \n\ --beforetc <data before t_c being analysed in seconds> \n\ --aftertc <data after t_c being analysed in seconds> \n\ --Flow <low frequency cut off in Hz> \n\ --Fhigh <high frequency cut off in Hz> \n\ --tukey1 <alpha1, the 1st parameter of the template tukey window. default 0.0> \n\ --tukey1 <alpha2, the 2nd parameter of the template tukey window. default 0.0> \n\ --nPSDsegment <number of segments to estimate the PSD> \n\ --lPSDsegment <length of each segment to estimate the PSD> \n\ --outputPath <output path> \n\ --cache <list of cache files, e.g. [cache1,cache2]> \n\ --channel <list of channels, e.g. [H1:LSC-STRAIN,L1:LSC-STRAIN]> \n\ --PSDstart <GPS time for the start of the PSD. default begining of cache file> \n\ --template <waveform template, 3 for spin, 4 for no spin> \n\example: ./lalapps_spinspiral -i ./pipeline/SPINspiral.input --mChirp 1.7 --eta 0.12 --tc 873739311.00000 --dist 15 --nIter 5 --nSkip 1 --downsample 1 --beforetc 5 --aftertc 2 --Flow 45 --Fhigh 1600.0 --nPSDsegment 32 --lPSDsegment 4 --network [1,2] --outputPath ./pipeline/ --cache [./pipeline/H-H1_RDS_C03_L2-873739103-873740831.cache,./pipeline/L-L1_RDS_C03_L2-873739055-873740847.cache] --channel [H1:LSC-STRAIN,L1:LSC-STRAIN] --PSDstart 873740000\n\\n\n" |
Referenced by readCommandLineOptions().
void readCommandLineOptions | ( | int | argc, |
char * | argv[], | ||
struct runPar * | run | ||
) |
Read and parse command-line arguments/options.
References runPar::cacheFilename, runPar::channelname, runPar::commandSettingsFlag, runPar::dataAfterTc, runPar::dataBeforeTc, runPar::downsampleFactor, runPar::highFrequencyCut, runPar::injXMLfilename, runPar::injXMLnr, runPar::lowFrequencyCut, runPar::mainFilename, runPar::mcmcFilename, runPar::MCMCseed, runPar::mcmcWaveform, runPar::networkSize, runPar::nIter, runPar::nMCMCpar, runPar::outputPath, parseCharacterOptionString(), runPar::PSDsegmentLength, runPar::PSDsegmentNumber, runPar::PSDstart, readCachefile(), runPar::selectifos, runPar::thinOutput, runPar::triggerDist, runPar::triggerEta, runPar::triggerMc, runPar::triggerTc, runPar::tukey1, runPar::tukey2, and USAGE.
Referenced by main().
void readMainInputfile | ( | struct runPar * | run | ) |
Read the main input file.
All parameters that are read here should be(come) members of the runvar struct and lose their global status
References runPar::beVerbose, runPar::dataFilename, runPar::doMatch, runPar::doMCMC, runPar::doSNR, runPar::injectionFilename, runPar::mainFilename, runPar::mcmcFilename, runPar::parameterFilename, runPar::systemFilename, UNUSED, and runPar::writeSignal.
Referenced by main().
void readMCMCinputfile | ( | struct runPar * | run | ) |
Read the MCMC input file.
All parameters that are read here should be(come) members of the runvar struct and lose their global status
References runPar::acceptRateTarget, runPar::adaptiveMCMC, runPar::annealNburn, runPar::annealNburn0, runPar::annealTemp0, runPar::blockFrac, runPar::commandSettingsFlag, runPar::correlatedUpdates, runPar::corrFrac, runPar::matAccFr, runPar::maxTemp, runPar::mcmcFilename, runPar::MCMCseed, runPar::minlogL, runPar::nCorr, runPar::nIter, runPar::nTemps, runPar::parallelTempering, runPar::prMatrixInfo, runPar::prParTempInfo, runPar::saveHotChains, runPar::tempLadder, runPar::thinOutput, and runPar::thinScreenOutput.
Referenced by main().
void readDataInputfile | ( | struct runPar * | run, |
struct interferometer | ifo[] | ||
) |
Read the data input file.
References interferometer::add2channels, interferometer::ch1doubleprecision, interferometer::ch1fileoffset, interferometer::ch1filesize, runPar::channelname, runPar::commandSettingsFlag, runPar::dataAfterTc, runPar::dataBeforeTc, runPar::dataDir, runPar::dataFilename, runPar::datasetName, runPar::downsampleFactor, runPar::highFrequencyCut, interferometer::lati, interferometer::leftArm, interferometer::longi, runPar::lowFrequencyCut, runPar::maxIFOdbaseSize, runPar::networkSize, interferometer::noisedoubleprecision, interferometer::noisefileoffset, interferometer::noisefilesize, interferometer::noiseGPSstart, pi, runPar::PSDsegmentLength, runPar::PSDsegmentNumber, interferometer::rightArm, runPar::selectifos, and runPar::tukeyWin.
Referenced by setIFOdata().
void readInjectionInputfile | ( | struct runPar * | run | ) |
Read the injection input file.
All parameters that are read in here should be(come) members of the runvar struct
References runPar::beVerbose, runPar::geocentricTc, runPar::injBoundLow, runPar::injBoundType, runPar::injBoundUp, runPar::injectionFilename, runPar::injectionPNorder, runPar::injectionSNR, runPar::injectionWaveform, runPar::injectSignal, runPar::injID, runPar::injNumber, runPar::injParUse, runPar::injParVal, runPar::injParValOrig, runPar::injRanPar, runPar::injRanSeed, runPar::injRevID, runPar::injSigma, runPar::nInjectPar, runPar::parAbrev, runPar::parDef, and setRandomInjectionParameters().
Referenced by main().
void readParameterInputfile | ( | struct runPar * | run | ) |
Read the (MCMC) parameter input file.
All parameters that are read in here should be(come) members of the runvar struct
References runPar::beVerbose, c3rd, runPar::commandSettingsFlag, runPar::geocentricTc, runPar::injectionFilename, runPar::injectSignal, runPar::injParVal, runPar::injRevID, runPar::mcmcFilename, runPar::mcmcParUse, runPar::mcmcPNorder, runPar::MCMCseed, runPar::mcmcWaveform, runPar::nMCMCpar, runPar::offsetMCMC, runPar::offsetX, runPar::parAbrev, runPar::parameterFilename, runPar::parBestVal, runPar::parDef, runPar::parFix, runPar::parID, runPar::parNumber, runPar::parRevID, runPar::parSigma, runPar::parStartMCMC, pi, runPar::priorBoundLow, runPar::priorBoundUp, runPar::priorSet, runPar::priorType, tpi, runPar::triggerDist, runPar::triggerEta, runPar::triggerMc, and runPar::triggerTc.
Referenced by main().
void readSystemInputfile | ( | struct runPar * | run | ) |
Read the input file for system (system-dependent) variables, e.g. SPINspiral.input.system.
References runPar::dataDir, runPar::parameterFilename, and runPar::systemFilename.
Referenced by main().
void readInjectionXML | ( | struct runPar * | run | ) |
Read an XML injection file if specified.
References runPar::geocentricTc, runPar::injID, runPar::injNumber, runPar::injParVal, runPar::injRevID, runPar::injXMLfilename, runPar::injXMLnr, runPar::lowFrequencyCutInj, runPar::nInjectPar, and runPar::parAbrev.
Referenced by main().
void setParameterNames | ( | struct runPar * | run | ) |
Set parameters names in the hardcoded parameter database.
References runPar::parAbrev, runPar::parAbrv, and runPar::parDef.
Referenced by main().
void setConstants | ( | void | ) |
void setIFOdata | ( | struct runPar * | run, |
struct interferometer | ifo[] | ||
) |
Set all the data for all IFOs that may be used.
References interferometer::after_tc, interferometer::before_tc, runPar::beVerbose, runPar::dataAfterTc, runPar::dataBeforeTc, runPar::dataFilename, runPar::datasetName, interferometer::highCut, runPar::highFrequencyCut, interferometer::lowCut, runPar::lowFrequencyCut, runPar::mainFilename, runPar::maxIFOdbaseSize, interferometer::radius_eqt, interferometer::radius_pole, readDataInputfile(), and runPar::selectdata.
Referenced by main().
void setRandomInjectionParameters | ( | struct runPar * | run | ) |
Get random values for the injection parameters.
References runPar::injBoundLow, runPar::injBoundUp, runPar::injParVal, runPar::injParValOrig, runPar::injRanPar, runPar::injRanSeed, runPar::injSigma, max, min, runPar::nInjectPar, and setSeed().
Referenced by readInjectionInputfile().
void getInjectionParameters | ( | struct parSet * | par, |
int | nInjectionPar, | ||
double * | injParVal | ||
) |
Returns the 'injection values' to the parameter set par.
References parSet::nPar, and parSet::par.
Referenced by dataFT(), main(), matchBetweenParameterArrayAndTrueParameters(), MCMC(), and writeSignalsToFiles().
void getStartParameters | ( | struct parSet * | par, |
struct runPar | run | ||
) |
Returns the 'best-guess values' to the parameter set par.
References runPar::nMCMCpar, parSet::nPar, parSet::par, and runPar::parBestVal.
void startMCMCOffset | ( | struct parSet * | par, |
struct MCMCvariables * | mcmc, | ||
struct interferometer * | ifo[], | ||
struct runPar | run | ||
) |
Choose and print offset starting values for the Markov chain.
Set MCMC parameters to either the best-guess values or the injection values (where possible). Then, depending on the detailed per-parameter settings, add a random offset from this value. Require a good value for logL (determined by mcmc->minlogL) in order to accept the proposed starting values, unless MCMC parameters do not fully match the injection parameters. This happens independently of whether parameters are fixed or not! Finally, print the selected starting values to screen.
References MCMCvariables::acceptPrior, arr2par(), MCMCvariables::beVerbose, MCMCvariables::injectionWaveform, MCMCvariables::injectSignal, MCMCvariables::injID, MCMCvariables::injParUse, MCMCvariables::injParVal, MCMCvariables::injRevID, MCMCvariables::iTemp, localPar(), MCMCvariables::logL, max, MCMCvariables::mcmcWaveform, netLogLikelihood(), MCMCvariables::networkSize, MCMCvariables::nMCMCpar, MCMCvariables::nParam, MCMCvariables::offsetMCMC, MCMCvariables::offsetX, MCMCvariables::parAbrev, MCMCvariables::param, MCMCvariables::parBestVal, MCMCvariables::parID, MCMCvariables::parSigma, MCMCvariables::parStartMCMC, prior(), MCMCvariables::priorBoundLow, MCMCvariables::priorBoundUp, MCMCvariables::ran, and MCMCvariables::thinScreenOutput.
Referenced by MCMC().
void setTemperatureLadder | ( | struct MCMCvariables * | mcmc | ) |
Set up the temperature ladder for parallel tempering.
References MCMCvariables::beVerbose, MCMCvariables::maxTemp, MCMCvariables::nTemps, MCMCvariables::parallelTempering, MCMCvariables::prParTempInfo, MCMCvariables::tempAmpl, MCMCvariables::tempLadder, and MCMCvariables::tempOverlap.
void setTemperatureLadderOld | ( | struct MCMCvariables * | mcmc | ) |
Set up the temperature ladder for parallel tempering - old routine.
References MCMCvariables::beVerbose, MCMCvariables::maxTemp, min, MCMCvariables::nTemps, MCMCvariables::parallelTempering, MCMCvariables::prParTempInfo, MCMCvariables::tempAmpl, and MCMCvariables::tempLadder.
Referenced by MCMC().
void allocParset | ( | struct parSet * | par, |
int | networkSize | ||
) |
Allocate memory for the vectors in the struct parSet.
References parSet::localti, parSet::locazi, and parSet::loctc.
Referenced by dataFT(), main(), matchBetweenParameterArrayAndTrueParameters(), MCMC(), and writeSignalsToFiles().
void freeParset | ( | struct parSet * | par | ) |
Deallocate memory for the vectors in the struct parSet.
References parSet::localti, parSet::locazi, and parSet::loctc.
Referenced by dataFT(), main(), matchBetweenParameterArrayAndTrueParameters(), MCMC(), and writeSignalsToFiles().
void copyRun2MCMC | ( | struct runPar | run, |
struct MCMCvariables * | mcmc | ||
) |
Copy some of the elements of the struct runPar to the struct MCMCvariables.
References runPar::acceptRateTarget, MCMCvariables::acceptRateTarget, runPar::adaptiveMCMC, MCMCvariables::adaptiveMCMC, runPar::annealNburn, MCMCvariables::annealNburn, runPar::annealNburn0, MCMCvariables::annealNburn0, runPar::annealTemp0, MCMCvariables::annealTemp0, MCMCvariables::baseTime, runPar::beVerbose, MCMCvariables::beVerbose, runPar::blockFrac, MCMCvariables::blockFrac, MCMCvariables::chTemp, runPar::correlatedUpdates, MCMCvariables::correlatedUpdates, runPar::corrFrac, MCMCvariables::corrFrac, runPar::geocentricTc, MCMCvariables::geocentricTc, runPar::injectionPNorder, MCMCvariables::injectionPNorder, runPar::injectionWaveform, MCMCvariables::injectionWaveform, runPar::injectSignal, MCMCvariables::injectSignal, runPar::injID, MCMCvariables::injID, runPar::injParUse, MCMCvariables::injParUse, runPar::injParVal, MCMCvariables::injParVal, runPar::injRevID, MCMCvariables::injRevID, runPar::matAccFr, MCMCvariables::matAccFr, max, runPar::maxnPar, MCMCvariables::maxnPar, runPar::maxTemp, MCMCvariables::maxTemp, runPar::mcmcParUse, MCMCvariables::mcmcParUse, runPar::mcmcPNorder, MCMCvariables::mcmcPNorder, runPar::MCMCseed, runPar::mcmcWaveform, MCMCvariables::mcmcWaveform, runPar::minlogL, MCMCvariables::minlogL, runPar::nCorr, MCMCvariables::nCorr, runPar::networkSize, MCMCvariables::networkSize, runPar::nInjectPar, MCMCvariables::nInjectPar, runPar::nIter, MCMCvariables::nIter, runPar::nMCMCpar, MCMCvariables::nMCMCpar, runPar::nTemps, MCMCvariables::nTemps, runPar::offsetMCMC, MCMCvariables::offsetMCMC, runPar::offsetX, MCMCvariables::offsetX, runPar::parAbrev, MCMCvariables::parAbrev, runPar::parAbrv, MCMCvariables::parAbrv, runPar::parallelTempering, MCMCvariables::parallelTempering, runPar::parBestVal, MCMCvariables::parBestVal, runPar::parDBn, MCMCvariables::parDBn, runPar::parFix, MCMCvariables::parFix, runPar::parID, MCMCvariables::parID, runPar::parName, MCMCvariables::parName, runPar::parNumber, MCMCvariables::parNumber, runPar::parRevID, MCMCvariables::parRevID, runPar::parSigma, MCMCvariables::parSigma, runPar::parStartMCMC, MCMCvariables::parStartMCMC, runPar::priorBoundLow, MCMCvariables::priorBoundLow, runPar::priorBoundUp, MCMCvariables::priorBoundUp, runPar::priorType, MCMCvariables::priorType, runPar::prMatrixInfo, MCMCvariables::prMatrixInfo, runPar::prParTempInfo, MCMCvariables::prParTempInfo, runPar::saveHotChains, MCMCvariables::saveHotChains, MCMCvariables::seed, runPar::tempLadder, MCMCvariables::tempLadder, runPar::thinOutput, MCMCvariables::thinOutput, runPar::thinScreenOutput, and MCMCvariables::thinScreenOutput.
Referenced by MCMC().
void setMCMCseed | ( | struct runPar * | run | ) |
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().
void MCMC | ( | struct runPar | run, |
struct interferometer * | ifo[] | ||
) |
MCMC routine - forms the MCMC core of the program.
Initialise and build a Markov chain
References MCMCvariables::acceptPrior, MCMCvariables::acceptRateTarget, MCMCvariables::adaptiveMCMC, MCMCvariables::adaptScale, MCMCvariables::adaptSigma, allocateMCMCvariables(), allocParset(), MCMCvariables::annealNburn, MCMCvariables::annealNburn0, MCMCvariables::annealTemp0, annealTemperature(), arr2par(), MCMCvariables::baseTime, MCMCvariables::beVerbose, MCMCvariables::blockFrac, MCMCvariables::chTemp, copyRun2MCMC(), correlatedMCMCupdate(), MCMCvariables::correlatedUpdates, MCMCvariables::corrFrac, MCMCvariables::corrUpdate, MCMCvariables::covar, MCMCvariables::decreaseSigma, MCMCvariables::dlogL, MCMCvariables::fouts, freeMCMCvariables(), freeParset(), getInjectionParameters(), getStartParameters(), MCMCvariables::hist, MCMCvariables::iHist, MCMCvariables::iIter, MCMCvariables::increaseSigma, MCMCvariables::injectionWaveform, MCMCvariables::injectSignal, MCMCvariables::injID, MCMCvariables::injParUse, MCMCvariables::injParVal, MCMCvariables::injRevID, MCMCvariables::iTemp, localPar(), MCMCvariables::logL, max, MCMCvariables::maxdlogL, MCMCvariables::maxLparam, MCMCvariables::mcmcWaveform, MCMCvariables::minlogL, MCMCvariables::nCorr, netLogLikelihood(), MCMCvariables::networkSize, MCMCvariables::nInjectPar, MCMCvariables::nIter, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, MCMCvariables::nParFit, MCMCvariables::nTemps, runPar::outputPath, par2arr(), MCMCvariables::parallelTempering, MCMCvariables::param, MCMCvariables::parBestVal, MCMCvariables::parFix, MCMCvariables::parID, MCMCvariables::parSigma, MCMCvariables::prParTempInfo, MCMCvariables::ran, MCMCvariables::saveHotChains, MCMCvariables::seed, setTemperatureLadderOld(), sigmaPeriodicBoundaries(), startMCMCOffset(), swapChains(), MCMCvariables::tempAmpl, MCMCvariables::tempLadder, MCMCvariables::tempOverlap, tpi, uncorrelatedMCMCblockUpdate(), uncorrelatedMCMCsingleUpdate(), updateCovarianceMatrix(), writeChainInfo(), writeMCMCheader(), and writeMCMCoutput().
Referenced by main().
void CholeskyDecompose | ( | double ** | A, |
struct MCMCvariables * | mcmc | ||
) |
Compute Cholesky decomposition matrix.
Performs Cholesky decompositon of matrix A and returns result in the same matrix - adapted from PJG Fortran function If matrix is not positive definite, return zeroes
References MCMCvariables::nMCMCpar, and MCMCvariables::parFix.
Referenced by updateCovarianceMatrix().
void par2arr | ( | struct parSet | par, |
double ** | param, | ||
struct MCMCvariables | mcmc | ||
) |
Put the MCMC parameters in an array.
References MCMCvariables::iTemp, MCMCvariables::nMCMCpar, and parSet::par.
Referenced by correlatedMCMCupdate(), MCMC(), uncorrelatedMCMCblockUpdate(), and uncorrelatedMCMCsingleUpdate().
void arr2par | ( | double ** | param, |
struct parSet * | par, | ||
struct MCMCvariables | mcmc | ||
) |
Get the MCMC parameters out of their array.
References MCMCvariables::iTemp, MCMCvariables::nMCMCpar, and parSet::par.
Referenced by correlatedMCMCupdate(), MCMC(), startMCMCOffset(), uncorrelatedMCMCblockUpdate(), and uncorrelatedMCMCsingleUpdate().
double prior | ( | double * | par, |
int | p, | ||
struct MCMCvariables | mcmc | ||
) |
Compute the prior for the given parameter set.
Contains boundary conditions and prior information for the MCMC. Try to avoid returning 0, to increase jump sizes
References mtpi, pi, MCMCvariables::priorBoundLow, MCMCvariables::priorBoundUp, MCMCvariables::priorType, and tpi.
Referenced by correlatedMCMCupdate(), startMCMCOffset(), uncorrelatedMCMCblockUpdate(), and uncorrelatedMCMCsingleUpdate().
double sigmaPeriodicBoundaries | ( | double | sigma, |
int | p, | ||
struct MCMCvariables | mcmc | ||
) |
Bring the adaptation sigma between its periodic boundaries.
References min, pi, MCMCvariables::priorType, and tpi.
Referenced by MCMC(), and uncorrelatedMCMCsingleUpdate().
void correlatedMCMCupdate | ( | struct interferometer * | ifo[], |
struct parSet * | state, | ||
struct MCMCvariables * | mcmc, | ||
struct runPar | run | ||
) |
Do a correlated block MCMC update.
Do an update for all non-fixed MCMC parameters. Use the covariance matrix to take into account correlations. The covariance matrix has been constructed from previous iterations. Use adaptation to scale the whole matrix. Some experiments with larger jumps using the 'hotter' covariance matrix.
References MCMCvariables::accepted, MCMCvariables::acceptPrior, MCMCvariables::adaptiveMCMC, MCMCvariables::adaptSigmaOut, arr2par(), MCMCvariables::chTemp, MCMCvariables::corrSig, MCMCvariables::covar, MCMCvariables::decreaseSigma, MCMCvariables::increaseSigma, MCMCvariables::iTemp, localPar(), MCMCvariables::logL, max, MCMCvariables::mcmcWaveform, min, MCMCvariables::minlogL, netLogLikelihood(), MCMCvariables::networkSize, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, par2arr(), MCMCvariables::param, MCMCvariables::parFix, prior(), and MCMCvariables::ran.
Referenced by MCMC().
void uncorrelatedMCMCsingleUpdate | ( | struct interferometer * | ifo[], |
struct parSet * | state, | ||
struct MCMCvariables * | mcmc, | ||
struct runPar | run | ||
) |
Do an uncorrelated, per-parameter MCMC update.
Do an update for all non-fixed MCMC parameters in the current T chain. Propose a jump and decide whether to accept or not on a per-parameter basis. Use adaptation.
References MCMCvariables::accepted, MCMCvariables::acceptPrior, MCMCvariables::acceptRateTarget, MCMCvariables::adaptiveMCMC, MCMCvariables::adaptScale, MCMCvariables::adaptSigma, MCMCvariables::adaptSigmaOut, arr2par(), MCMCvariables::chTemp, MCMCvariables::iIter, MCMCvariables::iTemp, localPar(), MCMCvariables::logL, max, MCMCvariables::mcmcWaveform, min, MCMCvariables::minlogL, netLogLikelihood(), MCMCvariables::networkSize, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, par2arr(), MCMCvariables::param, MCMCvariables::parFix, prior(), MCMCvariables::ran, and sigmaPeriodicBoundaries().
Referenced by MCMC().
void uncorrelatedMCMCblockUpdate | ( | struct interferometer * | ifo[], |
struct parSet * | state, | ||
struct MCMCvariables * | mcmc, | ||
struct runPar | run | ||
) |
Do an uncorrelated block update.
Do an update for all non-fixed MCMC parameters in the current T chain. Propose a jump and decide whether to accept or not for all parameters at once. No adaptation here, some experimenting with larger jumps every now and then.
References MCMCvariables::accepted, MCMCvariables::acceptPrior, MCMCvariables::adaptSigma, arr2par(), MCMCvariables::chTemp, MCMCvariables::iTemp, localPar(), MCMCvariables::logL, max, MCMCvariables::mcmcWaveform, min, MCMCvariables::minlogL, netLogLikelihood(), MCMCvariables::networkSize, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, par2arr(), MCMCvariables::param, MCMCvariables::parFix, prior(), and MCMCvariables::ran.
Referenced by MCMC().
void writeMCMCheader | ( | struct interferometer * | ifo[], |
struct MCMCvariables | mcmc, | ||
struct runPar | run | ||
) |
Write MCMC header to screen and file.
References interferometer::after_tc, MCMCvariables::annealNburn, interferometer::before_tc, interferometer::deltaFT, MCMCvariables::fouts, interferometer::FTsize, interferometer::FTstart, interferometer::highCut, interferometer::lowCut, MCMCvariables::maxTemp, runPar::mcmcPNorder, runPar::mcmcWaveform, interferometer::name, MCMCvariables::nCorr, runPar::netsnr, runPar::networkSize, MCMCvariables::nIter, runPar::nMCMCpar, MCMCvariables::nMCMCpar, MCMCvariables::nTemps, MCMCvariables::offsetMCMC, MCMCvariables::parAbrev, MCMCvariables::parID, interferometer::samplerate, interferometer::samplesize, MCMCvariables::saveHotChains, MCMCvariables::seed, interferometer::snr, and MCMCvariables::tempLadder.
Referenced by MCMC().
void writeMCMCoutput | ( | struct MCMCvariables | mcmc, |
struct interferometer * | ifo[] | ||
) |
Write an MCMC iteration as an output line to screen and/or file.
References MCMCvariables::accepted, MCMCvariables::fouts, MCMCvariables::iIter, interferometer::index, MCMCvariables::iTemp, MCMCvariables::logL, MCMCvariables::nMCMCpar, MCMCvariables::parAbrev, MCMCvariables::param, MCMCvariables::parID, MCMCvariables::saveHotChains, MCMCvariables::thinOutput, and MCMCvariables::thinScreenOutput.
Referenced by MCMC().
void allocateMCMCvariables | ( | struct MCMCvariables * | mcmc | ) |
Allocate memory for the MCMCvariables struct.
Allocate memory for the MCMCvariables struct. Don't forget to deallocate whatever you put here in freeMCMCvariables()
References MCMCvariables::accepted, MCMCvariables::acceptElems, MCMCvariables::acceptPrior, MCMCvariables::adaptScale, MCMCvariables::adaptSigma, MCMCvariables::adaptSigmaOut, MCMCvariables::corrSig, MCMCvariables::corrUpdate, MCMCvariables::covar, MCMCvariables::dlogL, MCMCvariables::hist, MCMCvariables::histDev, MCMCvariables::histMean, MCMCvariables::iHist, MCMCvariables::logL, max, MCMCvariables::maxdlogL, MCMCvariables::maxLparam, MCMCvariables::nCorr, MCMCvariables::nInjectPar, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, MCMCvariables::nParFit, MCMCvariables::nTemps, MCMCvariables::param, MCMCvariables::swapTs1, MCMCvariables::swapTs2, MCMCvariables::swapTss, and MCMCvariables::tempAmpl.
Referenced by MCMC().
void freeMCMCvariables | ( | struct MCMCvariables * | mcmc | ) |
Deallocate memory for the MCMCvariables struct.
References MCMCvariables::accepted, MCMCvariables::acceptElems, MCMCvariables::acceptPrior, MCMCvariables::adaptScale, MCMCvariables::adaptSigma, MCMCvariables::adaptSigmaOut, MCMCvariables::corrSig, MCMCvariables::corrUpdate, MCMCvariables::covar, MCMCvariables::dlogL, MCMCvariables::hist, MCMCvariables::histDev, MCMCvariables::histMean, MCMCvariables::iHist, MCMCvariables::logL, MCMCvariables::maxdlogL, MCMCvariables::maxLparam, MCMCvariables::nlogL, MCMCvariables::nMCMCpar, MCMCvariables::nParam, MCMCvariables::nTemps, MCMCvariables::param, MCMCvariables::ran, MCMCvariables::swapTs1, MCMCvariables::swapTs2, MCMCvariables::swapTss, and MCMCvariables::tempAmpl.
Referenced by MCMC().
void updateCovarianceMatrix | ( | struct MCMCvariables * | mcmc | ) |
Calculate the new covariance matrix and determine whether the new matrix should be accepted.
Compute the new covariance matrix for the current T chain. Determine by the 'improvement' of the new matrix whether it should be accepted.
References MCMCvariables::acceptElems, CholeskyDecompose(), MCMCvariables::corrUpdate, MCMCvariables::covar, MCMCvariables::hist, MCMCvariables::histDev, MCMCvariables::histMean, MCMCvariables::iIter, MCMCvariables::iTemp, MCMCvariables::matAccFr, max, MCMCvariables::nCorr, MCMCvariables::nMCMCpar, MCMCvariables::nParFit, MCMCvariables::nTemps, MCMCvariables::parFix, and MCMCvariables::prMatrixInfo.
Referenced by MCMC().
double annealTemperature | ( | double | temp0, |
int | nburn, | ||
int | nburn0, | ||
int | iIter | ||
) |
void swapChains | ( | struct MCMCvariables * | mcmc | ) |
Parallel tempering: Swap states between two chains.
References MCMCvariables::logL, max, min, MCMCvariables::nMCMCpar, MCMCvariables::nTemps, MCMCvariables::param, MCMCvariables::ran, MCMCvariables::swapTs1, MCMCvariables::swapTs2, MCMCvariables::swapTss, and MCMCvariables::tempLadder.
Referenced by MCMC().
void writeChainInfo | ( | struct MCMCvariables | mcmc | ) |
Parallel tempering: Print chain and swap info to screen.
References MCMCvariables::accepted, MCMCvariables::acceptElems, MCMCvariables::chTemp, MCMCvariables::corrUpdate, MCMCvariables::histDev, MCMCvariables::iIter, MCMCvariables::iTemp, MCMCvariables::nMCMCpar, MCMCvariables::nTemps, MCMCvariables::parAbrv, MCMCvariables::parID, MCMCvariables::prParTempInfo, MCMCvariables::swapTs1, MCMCvariables::swapTs2, and MCMCvariables::swapTss.
Referenced by MCMC().
double chirpMass | ( | double | m1, |
double | m2 | ||
) |
Compute chirp mass (mM) for given individual masses.
double massRatio | ( | double | m1, |
double | m2 | ||
) |
Compute symmetric mass ratio (eta) 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 | ||
) |
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.
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 IFOinit | ( | struct interferometer ** | ifo, |
int | networkSize, | ||
struct runPar | run | ||
) |
Set up IFOs.
Determines interferometer arm (unit-) vectors (and others), given position (lat/long) and arm angles. Vectors refer to the (right-handed) Earth coordinate system spanned by the three vectors: x) from geocenter to intersection of greenwich meridian with equator plane; y) from geocenter to intersection of 90E meridian with equator plane; z) from geocenter to north pole.
References interferometer::after_tc, interferometer::before_tc, runPar::beVerbose, coord2vec(), dataFT(), interferometer::dataTrafo, interferometer::deltaFT, interferometer::FTsize, runPar::geocentricTc, interferometer::highCut, interferometer::highIndex, interferometer::index, interferometer::indexRange, interpolLogNoisePSD(), interferometer::leftvec, interferometer::lowCut, interferometer::lowIndex, interferometer::noisePSD, noisePSDestimate(), interferometer::normalvec, interferometer::orthoArm, pi, interferometer::positionvec, runPar::PSDsegmentLength, runPar::PSDsegmentNumber, r2d, interferometer::radius_eqt, interferometer::radius_pole, interferometer::raw_dataTrafo, interferometer::rightvec, rotate(), and runPar::tukeyWin.
Referenced by main().
void IFOdispose | ( | struct interferometer * | ifo, |
struct runPar | run | ||
) |
Free allocated IFO variables.
References runPar::beVerbose, interferometer::dataTrafo, interferometer::FTin, interferometer::FTout, interferometer::FTplan, interferometer::FTwindow, interferometer::index, interferometer::name, interferometer::noisePSD, interferometer::raw_dataTrafo, interferometer::raw_noisePSD, and interferometer::rawDownsampledWindowedData.
Referenced by main().
Compute FIR filter coefficients.
The filter has a total of N=(order+order-1) coefficients that are symmetric, i.e. coef[k] = coef[N-k]. For details & filter properties see the function downsampling() below.
References BANDPASS, runPar::beVerbose, runPar::downsampleFactor, and remez().
Referenced by dataFT(), and noisePSDestimate().
double* downsample | ( | double | data[], |
int * | datalength, | ||
double | filtercoef[], | ||
int | ncoef, | ||
struct runPar | run | ||
) |
Downsample a time series by a factor downsampleFactor.
Downsamples a time series by factor downsampleFactor by first low-pass filtering it using a finite-impulse-response (FIR) filter and then thinning the data. Filter coefficients are determined using the 'Parks-McClellan' or 'Remez exchange' algorithm. The resulting data vector is shorter than original. Returned vector is allocated using fftw_malloc() and thus must be freed again using fftw_free().
References runPar::downsampleFactor.
Referenced by dataFT(), and noisePSDestimate().
void dataFT | ( | struct interferometer * | ifo[], |
int | ifonr, | ||
int | networkSize, | ||
struct runPar | run | ||
) |
Read the data and do a software injection if wanted.
Computes the Fourier Transform for the specified range of the specified Frame (".gwf") file, after adding up the two (signal & noise) channels, or injecting a waveform template into the noise. Also takes care of preparing FT stuff (ifo[ifonr]->FTplan, ->FTin, ->FTout, ...).
References interferometer::after_tc, allocParset(), interferometer::before_tc, runPar::beVerbose, interferometer::ch1fileoffset, interferometer::ch1filesize, interferometer::ch2fileoffset, interferometer::ch2filesize, runPar::commandSettingsFlag, interferometer::deltaFT, downsample(), runPar::downsampleFactor, filter(), runPar::FrameGPSstart, runPar::FrameLength, runPar::FrameName, freeParset(), interferometer::FTin, interferometer::FTout, interferometer::FTplan, interferometer::FTsize, interferometer::FTstart, interferometer::FTwindow, runPar::geocentricTc, getInjectionParameters(), GMST(), runPar::injectionWaveform, runPar::injectSignal, runPar::injParVal, runPar::injRevID, localPar(), parSet::localti, parSet::locazi, parSet::loctc, McEta2masses(), runPar::nFrame, runPar::nInjectPar, parSet::par, pi, interferometer::raw_dataTrafo, interferometer::rawDownsampledWindowedData, interferometer::samplerate, interferometer::samplesize, runPar::tukeyWin, tukeyWindow(), and waveformTemplate().
Referenced by IFOinit().
double hannWindow | ( | int | j, |
int | N | ||
) |
Apply a 'Hann window' to data.
Window data using a Hann window. j = 0, ..., N-1
References pi.
Referenced by noisePSDestimate().
double tukeyWindow | ( | int | j, |
int | N, | ||
double | r | ||
) |
double modifiedTukeyWindow | ( | int | j, |
int | N, | ||
double | r1, | ||
double | r2 | ||
) |
Apply a 'Tukey window' to data.
Window data using a Tukey window. r1 for lower frequency windowing, r2 for upper frequency windowing. For r=0 equal to rectangular window, for r=1 equal to Hann window (0<r<1 denotes the fraction of the window in which it behaves sinusoidal). j = 0, ..., N-1
References pi.
Referenced by waveformTemplate().
void noisePSDestimate | ( | struct interferometer * | ifo[], |
int | ifonr, | ||
struct runPar | run | ||
) |
Returns a (smoothed) estimate of the log- Power Spectral Density.
Data is split into K segments of M seconds, and K-1 overlapping segments of length 2M are eventually windowed and transformed.
References runPar::beVerbose, runPar::commandSettingsFlag, downsample(), runPar::downsampleFactor, filter(), runPar::FrameGPSstart, runPar::FrameLength, runPar::FrameName, hannWindow(), runPar::nFrame, interferometer::noisefileoffset, interferometer::noisefilesize, interferometer::noiseGPSstart, runPar::PSDsegmentLength, runPar::PSDsegmentNumber, interferometer::PSDsize, runPar::PSDstart, and interferometer::raw_noisePSD.
Referenced by IFOinit().
double logNoisePSD | ( | double | f, |
struct interferometer * | ifo | ||
) |
double interpolLogNoisePSD | ( | double | f, |
struct interferometer * | ifo | ||
) |
Returns a linearly interpolated (log-) noise PSD.
References interferometer::highCut, interferometer::lowCut, interferometer::PSDsize, and interferometer::raw_noisePSD.
Referenced by IFOinit().
void writeDataToFiles | ( | struct interferometer * | ifo[], |
int | networkSize, | ||
struct runPar | run | ||
) |
Write windowed, time-domain data (signal + noise) to disc.
References runPar::beVerbose, interferometer::deltaFT, interferometer::indexRange, runPar::MCMCseed, printParameterHeaderToFile(), interferometer::samplesize, and runPar::writeSignal.
Referenced by main().
void writeNoiseToFiles | ( | struct interferometer * | ifo[], |
int | networkSize, | ||
struct runPar | run | ||
) |
Write noise ASD to disc.
The noise ASD is the square root of the estimated noise PSD (no injected signal).
References runPar::beVerbose, interferometer::deltaFT, interferometer::indexRange, runPar::MCMCseed, printParameterHeaderToFile(), and runPar::writeSignal.
Referenced by main().
void writeSignalsToFiles | ( | struct interferometer * | ifo[], |
int | networkSize, | ||
struct runPar | run | ||
) |
Write signal and its FFT to disc.
Write a signal with the injection parameters and its FFT to disc, as *-signal.dat.* and *-signalFFT.dat.*.
References allocParset(), runPar::beVerbose, freeParset(), interferometer::FTout, getInjectionParameters(), interferometer::highIndex, runPar::injectionWaveform, runPar::injParVal, localPar(), runPar::MCMCseed, runPar::nInjectPar, printParameterHeaderToFile(), interferometer::samplesize, waveformTemplate(), and runPar::writeSignal.
Referenced by main().
void printParameterHeaderToFile | ( | FILE * | dump, |
struct interferometer * | ifo, | ||
struct runPar | run | ||
) |
Write an optional header when saving data, ASD or signal to disc.
The header idendtifies the (waveform) parameters.
References runPar::injID, runPar::injParVal, runPar::nInjectPar, runPar::parAbrev, and interferometer::snr.
Referenced by writeDataToFiles(), writeNoiseToFiles(), and writeSignalsToFiles().
void waveformTemplate | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute an inspiral waveform.
Compute an inspiral waveform template. waveformVersion determines the template to use. injectionWF indicates whether this is an injection waveform (1) or not (0).
References modifiedTukeyWindow(), interferometer::samplesize, templateApostolatos(), templateLAL12(), templateLAL15(), templateLALnonSpinning(), templateLALPhenSpinTaylorRD(), runPar::tukey1, and runPar::tukey2.
Referenced by dataFT(), IFOlogLikelihood(), signalFFT(), signalToNoiseRatio(), and writeSignalsToFiles().
void templateApostolatos | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute an Apostolatos, 12-parameter, single spin, simple-precession waveform.
Compute a spinning, 'simple-precession' template in restricted 1.5PN order with 1 spin (Apostolatos et al., 1994, PhRvD..49.6274A). The output vector ifo[ifonr]->FTin is of length ifo[ifonr]->samplesize, starting at 'tstart'(?) and with resolution ifo[ifonr]->samplerate.
References addVec(), c3rd, crossProduct(), dotProduct(), facVec(), interferometer::FTin, GMST(), interferometer::highCut, runPar::injParUse, runPar::injRevID, parSet::localti, parSet::locazi, parSet::loctc, longitude(), interferometer::lowCut, M0, McEta2masses(), runPar::mcmcParUse, Mpcs, mtpi, parSet::NdJ, parSet::par, runPar::parRevID, pi, interferometer::samplesize, and tpi.
Referenced by waveformTemplate().
void localPar | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | networkSize, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Calculate the local parameters from the global parameters.
Calculate the local (i.e. in the detector frame) parameters from the global parameters. par : pointer to parameter set (struct) ifo : pointer to interferometer data (struct)
References angle(), c, coord2vec(), interferometer::FTstart, GMST(), runPar::injRevID, parSet::localti, parSet::locazi, parSet::loctc, longitude(), mtpi, orthoProject(), parSet::par, runPar::parRevID, pi, interferometer::positionvec, rightHanded(), and tpi.
Referenced by correlatedMCMCupdate(), dataFT(), main(), matchBetweenParameterArrayAndTrueParameters(), MCMC(), startMCMCOffset(), uncorrelatedMCMCblockUpdate(), uncorrelatedMCMCsingleUpdate(), and writeSignalsToFiles().
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 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 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().
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 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 LALHpHcNonSpinning | ( | LALStatus * | status, |
CoherentGW * | waveform, | ||
SimInspiralTable * | injParams, | ||
PPNParamStruc * | ppnParams, | ||
int * | l, | ||
struct parSet * | par, | ||
struct interferometer * | ifo, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
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 LALfreedomNoSpin | ( | CoherentGW * | waveform | ) |
Free non-spinning LAL variables.
Referenced by templateLALnonSpinning().
void LALfreedomPhenSpinTaylorRD | ( | CoherentGW * | waveform | ) |
Free spinning LAL variables.
Referenced by templateLALPhenSpinTaylorRD().
double netLogLikelihood | ( | struct parSet * | par, |
int | networkSize, | ||
struct interferometer * | ifo[], | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute the log(Likelihood) for a network of IFOs.
References IFOlogLikelihood(), and logLikelihood_nine().
Referenced by correlatedMCMCupdate(), MCMC(), startMCMCOffset(), uncorrelatedMCMCblockUpdate(), and uncorrelatedMCMCsingleUpdate().
double IFOlogLikelihood | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | i, | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute the log(Likelihood) for a single IFO.
References vecOverlap(), and waveformTemplate().
Referenced by netLogLikelihood().
double logLikelihood_nine | ( | struct parSet * | par, |
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute a analytical log(Likelihood)
References runPar::nMCMCpar, parSet::par, and runPar::parRevID.
Referenced by netLogLikelihood().
double signalToNoiseRatio | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | i, | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute the SNR of the waveform with a given parameter set for a single IFO.
References vecOverlap(), and waveformTemplate().
Referenced by main().
double parMatch | ( | struct parSet * | par1, |
int | waveformVersion1, | ||
int | injectionWF1, | ||
struct parSet * | par2, | ||
int | waveformVersion2, | ||
int | injectionWF2, | ||
struct interferometer * | ifo[], | ||
int | networkSize, | ||
struct runPar | run | ||
) |
Compute match between waveforms with parameter sets par1 and par2.
References signalFFT(), and vecOverlap().
Referenced by main(), and matchBetweenParameterArrayAndTrueParameters().
double overlapWithData | ( | struct parSet * | par, |
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute frequency-domain overlap of waveform of given parameters with raw data.
References signalFFT(), and vecOverlap().
double parOverlap | ( | struct parSet * | par1, |
int | waveformVersion1, | ||
int | injectionWF1, | ||
struct parSet * | par2, | ||
int | waveformVersion2, | ||
int | injectionWF2, | ||
struct interferometer * | ifo[], | ||
int | ifonr, | ||
struct runPar | run | ||
) |
Compute the overlap in the frequency domain between two waveforms with parameter sets par1 and par2.
References signalFFT(), and vecOverlap().
Referenced by main().
double vecOverlap | ( | fftw_complex * | vec1, |
fftw_complex * | vec2, | ||
double * | noise, | ||
int | j_1, | ||
int | j_2, | ||
double | deltaFT | ||
) |
Compute the overlap in the frequency domain between two waveforms with parameter sets par1 and par2.
Referenced by IFOlogLikelihood(), overlapWithData(), parMatch(), parOverlap(), and signalToNoiseRatio().
void signalFFT | ( | fftw_complex * | FFTout, |
struct parSet * | par, | ||
struct interferometer * | ifo[], | ||
int | ifonr, | ||
int | waveformVersion, | ||
int | injectionWF, | ||
struct runPar | run | ||
) |
Compute the FFT of a waveform with given parameter set.
References interferometer::FTsize, and waveformTemplate().
Referenced by overlapWithData(), parMatch(), and parOverlap().
double matchBetweenParameterArrayAndTrueParameters | ( | double * | pararray, |
struct interferometer * | ifo[], | ||
struct MCMCvariables | mcmc, | ||
struct runPar | run | ||
) |
Compute the match between a given parameter array and the injection parameters.
References allocParset(), freeParset(), getInjectionParameters(), MCMCvariables::injectionWaveform, MCMCvariables::injParVal, localPar(), MCMCvariables::mcmcWaveform, MCMCvariables::networkSize, MCMCvariables::nInjectPar, MCMCvariables::nMCMCpar, parSet::par, and parMatch().
void parseCharacterOptionString | ( | char * | input, |
char ** | strings[], | ||
int * | n | ||
) |
parse strings from "[one,two,three]" to {"one", "two", "three"}
Read and parse command-line arguments/options Parses a character string (passed as one of the options) and decomposes it into individual parameter character strings. Input is of the form input : "[one,two,three]" and the resulting output is strings : {"one", "two", "three"} length of parameter names is by now limited to 512 characters. (should 'theoretically' (untested) be able to digest white space as well. Irrelevant for command line options, though.)
Referenced by readCommandLineOptions().
void readCachefile | ( | struct runPar * | run, |
int | ifonr | ||
) |
Read a Cache file. Returns an array of what is in the cache file.
References runPar::cacheFilename, runPar::FrameDetector, runPar::FrameGPSstart, runPar::FrameLength, runPar::FrameName, runPar::FramePrefix, and runPar::nFrame.
Referenced by readCommandLineOptions().
double Ms |
Referenced by setConstants().
double Mpc |
Referenced by setConstants().
double G |
Referenced by setConstants().
double c |
Referenced by ComputeA(), localPar(), remez(), and setConstants().
double Mpcs |
Referenced by setConstants(), and templateApostolatos().
double pi |
double tpi |
double mtpi |