Generate samples of various statistics (F-stat, F-atoms, B-stat,...) drawn from their respective distributions, assuming Gaussian noise, and drawing signal params from their (given) priors.
This is based on synthesizeBstat, and is mostly meant to be used for efficient Monte-Carlos studies, ROC curves etc
Prototypes | |
int | XLALDrawCorrelatedNoise (PulsarAmplitudeVect n_mu, const gsl_matrix *L, gsl_rng *rng) |
Generate 4 random-noise draws n_mu = {n_1, n_2, n_3, n_4} with correlations according to the matrix M = L L^T, which is passed in as input. More... | |
FstatAtomVector * | XLALGenerateFstatAtomVector (const DetectorStateSeries *detStates, const AMCoeffs *amcoeffs) |
Generate an FstatAtomVector for given antenna-pattern functions. More... | |
MultiFstatAtomVector * | XLALGenerateMultiFstatAtomVector (const MultiDetectorStateSeries *multiDetStates, const MultiAMCoeffs *multiAM) |
Generate a multi-FstatAtomVector for given antenna-pattern functions. More... | |
int | XLALAddNoiseToFstatAtomVector (FstatAtomVector *atoms, gsl_rng *rng) |
Add Gaussian-noise components to given FstatAtomVector. More... | |
int | XLALAddNoiseToMultiFstatAtomVector (MultiFstatAtomVector *multiAtoms, gsl_rng *rng) |
Add Gaussian-noise components to given multi-FstatAtomVector. More... | |
REAL8 | XLALAddSignalToFstatAtomVector (FstatAtomVector *atoms, AntennaPatternMatrix *M_mu_nu, const PulsarAmplitudeVect A_Mu, transientWindow_t transientWindow) |
Add signal s_mu = M_mu_nu A^nu within the given transient-window to given atoms. More... | |
REAL8 | XLALAddSignalToMultiFstatAtomVector (MultiFstatAtomVector *multiAtoms, AntennaPatternMatrix *M_mu_nu, const PulsarAmplitudeVect A_Mu, transientWindow_t transientWindow, INT4 lineX) |
Add given signal s_mu = M_mu_nu A^nu within the given transient-window to multi-IFO noise-atoms. More... | |
MultiFstatAtomVector * | XLALSynthesizeTransientAtoms (InjParams_t *injParamsOut, SkyPosition skypos, AmplitudePrior_t AmpPrior, transientWindowRange_t transientInjectRange, const MultiDetectorStateSeries *multiDetStates, BOOLEAN SignalOnly, multiAMBuffer_t *multiAMBuffer, gsl_rng *rng, INT4 lineX, const MultiNoiseWeights *multiNoiseWeights) |
Generates a multi-Fstat atoms vector for given parameters, drawing random parameters wherever required. More... | |
int | XLALRescaleMultiFstatAtomVector (MultiFstatAtomVector *multiAtoms, REAL8 rescale) |
Rescale a given multi-Fstat atoms vector {Fa,Fb} by given scalar factor. More... | |
int | write_InjParams_to_fp (FILE *fp, const InjParams_t *par, const UINT4 dataStartGPS, const BOOLEAN outputMmunuX, const UINT4 numDetectors) |
Write an injection-parameters structure to the given file-pointer, adding one line with the injection parameters. More... | |
Data Structures | |
struct | AmplitudePrior_t |
Signal (amplitude) parameter ranges. More... | |
struct | multiAMBuffer_t |
struct for buffering of AM-coeffs, if signal for same sky-position is injected More... | |
struct | InjParams_t |
Hold all (generally) randomly drawn injection parameters: skypos, amplitude-params, M_mu_nu, transient-window, SNR. More... | |
Enumerations | |
enum | AmpPriorType_t { AMP_PRIOR_TYPE_PHYSICAL = 0 , AMP_PRIOR_TYPE_CANONICAL , AMP_PRIOR_TYPE_LAST } |
Enumeration of allowed amplitude-prior types. More... | |
int XLALDrawCorrelatedNoise | ( | PulsarAmplitudeVect | n_mu, |
const gsl_matrix * | L, | ||
gsl_rng * | rng | ||
) |
Generate 4 random-noise draws n_mu = {n_1, n_2, n_3, n_4} with correlations according to the matrix M = L L^T, which is passed in as input.
Note: you need to pass a pre-allocated 4-vector n_mu. Note2: this function is meant as a lower-level noise-generation utility, called from a higher-level function to translate the antenna-pattern functions into pre-factorized Lcor
[out] | n_mu | 4d vector of noise-components {n_mu}, with correlation L * L^T |
[in] | L | correlator matrix to get n_mu = L_mu_nu * norm_nu from 4 uncorr. unit variates norm_nu |
rng | gsl random-number generator |
Definition at line 55 of file SynthesizeCWDraws.c.
FstatAtomVector * XLALGenerateFstatAtomVector | ( | const DetectorStateSeries * | detStates, |
const AMCoeffs * | amcoeffs | ||
) |
Generate an FstatAtomVector for given antenna-pattern functions.
Simply creates FstatAtomVector and initializes with antenna-pattern function.
detStates | input detector-state series, only used for timestamps |
amcoeffs | input antenna-pattern functions {a_i, b_i} |
Definition at line 104 of file SynthesizeCWDraws.c.
MultiFstatAtomVector * XLALGenerateMultiFstatAtomVector | ( | const MultiDetectorStateSeries * | multiDetStates, |
const MultiAMCoeffs * | multiAM | ||
) |
Generate a multi-FstatAtomVector for given antenna-pattern functions.
Simply creates MultiFstatAtomVector and initializes with antenna-pattern function.
[in] | multiDetStates | multi-detector state series, only used for timestamps |
multiAM | input antenna-pattern functions {a_i, b_i} |
Definition at line 156 of file SynthesizeCWDraws.c.
int XLALAddNoiseToFstatAtomVector | ( | FstatAtomVector * | atoms, |
gsl_rng * | rng | ||
) |
Add Gaussian-noise components to given FstatAtomVector.
atoms | input atoms-vector, noise will be added to this |
rng | random-number generator |
Definition at line 209 of file SynthesizeCWDraws.c.
int XLALAddNoiseToMultiFstatAtomVector | ( | MultiFstatAtomVector * | multiAtoms, |
gsl_rng * | rng | ||
) |
Add Gaussian-noise components to given multi-FstatAtomVector.
multiAtoms | input multi atoms-vector, noise will be added to this |
rng | random-number generator |
Definition at line 287 of file SynthesizeCWDraws.c.
REAL8 XLALAddSignalToFstatAtomVector | ( | FstatAtomVector * | atoms, |
AntennaPatternMatrix * | M_mu_nu, | ||
const PulsarAmplitudeVect | A_Mu, | ||
transientWindow_t | transientWindow | ||
) |
Add signal s_mu = M_mu_nu A^nu within the given transient-window to given atoms.
RETURN: SNR^2 of the injected signal and the effective AntennaPatternMatrix M_mu_nu for this signal.
atoms | [in/out] atoms vectors containing antenna-functions and possibly noise {Fa,Fb} | |
[out] | M_mu_nu | effective antenna-pattern matrix for the injected signal |
[in] | A_Mu | input canonical amplitude vector A^mu = {A1,A2,A3,A4} |
transientWindow | transient signal window |
Definition at line 325 of file SynthesizeCWDraws.c.
REAL8 XLALAddSignalToMultiFstatAtomVector | ( | MultiFstatAtomVector * | multiAtoms, |
AntennaPatternMatrix * | M_mu_nu, | ||
const PulsarAmplitudeVect | A_Mu, | ||
transientWindow_t | transientWindow, | ||
INT4 | lineX | ||
) |
Add given signal s_mu = M_mu_nu A^nu within the given transient-window to multi-IFO noise-atoms.
RETURN: SNR^2 of the injected signal and the effective AntennaPatternMatrix M_mu_nu for this signal.
multiAtoms | [in/out] multi atoms vectors containing antenna-functions and possibly noise {Fa,Fb} | |
[out] | M_mu_nu | effective multi-IFO antenna-pattern matrix for the injected signal |
[in] | A_Mu | input canonical amplitude vector A^mu = {A1,A2,A3,A4} |
[in] | transientWindow | transient signal window |
[in] | lineX | if >= 0: generate signal only for detector 'lineX': must be within 0,...(Ndet-1) |
Definition at line 461 of file SynthesizeCWDraws.c.
MultiFstatAtomVector * XLALSynthesizeTransientAtoms | ( | InjParams_t * | injParamsOut, |
SkyPosition | skypos, | ||
AmplitudePrior_t | AmpPrior, | ||
transientWindowRange_t | transientInjectRange, | ||
const MultiDetectorStateSeries * | multiDetStates, | ||
BOOLEAN | SignalOnly, | ||
multiAMBuffer_t * | multiAMBuffer, | ||
gsl_rng * | rng, | ||
INT4 | lineX, | ||
const MultiNoiseWeights * | multiNoiseWeights | ||
) |
Generates a multi-Fstat atoms vector for given parameters, drawing random parameters wherever required.
Input: detector states, signal sky-pos (or allsky), amplitudes (or range), transient window range
[out] | injParamsOut | return summary of injected signal parameters (can be NULL) |
skypos | (Alpha,Delta,system). Use Alpha < 0 to signal 'allsky' | |
[in] | AmpPrior | amplitude-parameter priors to draw signals from |
transientInjectRange | transient-window range for injections (flat priors) | |
[in] | multiDetStates | multi-detector state series covering observation time |
[in] | SignalOnly | switch to generate signal draws without noise |
multiAMBuffer | [in/out] buffer for AM-coefficients if re-using same skyposition (must be !=NULL) | |
rng | [in/out] gsl random-number generator | |
[in] | lineX | if >= 0: generate signal only for detector 'lineX': must be within 0,...(Ndet-1) |
[in] | multiNoiseWeights | per-detector noise weights SX^-1/S^-1, no per-SFT variation (can be NULL for unit weights) |
Definition at line 521 of file SynthesizeCWDraws.c.
int XLALRescaleMultiFstatAtomVector | ( | MultiFstatAtomVector * | multiAtoms, |
REAL8 | rescale | ||
) |
Rescale a given multi-Fstat atoms vector {Fa,Fb} by given scalar factor.
This is used to rescale signal atoms to desired fixed SNR.
multiAtoms | [in/out] multi atoms vectors containing a signal in {Fa,Fb} to be rescaled |
rescale | rescale factor: Fa' = rescale * Fa, and Fb'= rescale * Fb |
Definition at line 679 of file SynthesizeCWDraws.c.
int write_InjParams_to_fp | ( | FILE * | fp, |
const InjParams_t * | par, | ||
const UINT4 | dataStartGPS, | ||
const BOOLEAN | outputMmunuX, | ||
const UINT4 | numDetectors | ||
) |
Write an injection-parameters structure to the given file-pointer, adding one line with the injection parameters.
[in] | fp | file-pointer to output file |
[in] | par | injection params to write. NULL means write header-comment line |
[in] | dataStartGPS | data start-time in GPS seconds (used to turn window 't0' into offset from dataStartGPS) |
[in] | outputMmunuX | write per-IFO antenna pattern matrices? |
[in] | numDetectors | number of detectors, only needed to construct M_mu_nu_X_header_string |
Definition at line 717 of file SynthesizeCWDraws.c.
enum AmpPriorType_t |
Enumeration of allowed amplitude-prior types.
Enumerator | |
---|---|
AMP_PRIOR_TYPE_PHYSICAL | 'physical' priors: isotropic pdf{cosi,psi,phi0} AND flat pdf(h0) |
AMP_PRIOR_TYPE_CANONICAL | 'canonical' priors: uniform in A^mu up to h_max |
AMP_PRIOR_TYPE_LAST |
Definition at line 62 of file SynthesizeCWDraws.h.