LALPulsar  6.1.0.1-b72065a
Header TransientCW_utils.h

Detailed Description

Some helper functions useful for "transient CWs", mostly applying transient window functions.

Authors
Reinhard Prix, David Keitel

The approach is described by Prix, Giampanis & Messenger in https://arxiv.org/abs/1104.1704

Prototypes

int XLALParseTransientWindowName (const char *windowName)
 Parse a transient window name string into the corresponding transientWindowType. More...
 
int XLALGetTransientWindowTimespan (UINT4 *t0, UINT4 *t1, transientWindow_t transientWindow)
 Helper-function to determine the total timespan of a transient CW window, ie. More...
 
int XLALApplyTransientWindow (REAL4TimeSeries *series, transientWindow_t transientWindow)
 apply a "transient CW window" described by TransientWindowParams to the given timeseries More...
 
int XLALApplyTransientWindow2NoiseWeights (MultiNoiseWeights *multiNoiseWeights, const MultiLIGOTimeGPSVector *multiTS, transientWindow_t transientWindow)
 apply transient window to give multi noise-weights, associated with given multi timestamps More...
 
CHARXLALPulsarDopplerParams2String (const PulsarDopplerParams *par)
 Turn pulsar doppler-params into a single string that can be used for filenames The format is tRefNNNNNN_RAXXXXX_DECXXXXXX_FreqXXXXX[_f1dotXXXXX][_f2dotXXXXx][_f3dotXXXXX]. More...
 
REAL8 XLALComputeTransientBstat (transientWindowRange_t windowRange, const transientFstatMap_t *FstatMap)
 Compute transient-CW Bayes-factor B_SG = P(x|HypS)/P(x|HypG) (where HypG = Gaussian noise hypothesis), marginalized over start-time and timescale of transient CW signal, using given type and parameters of transient window range. More...
 
pdf1D_t * XLALComputeTransientPosterior_t0 (transientWindowRange_t windowRange, const transientFstatMap_t *FstatMap)
 Compute transient-CW posterior (normalized) on start-time t0, using given type and parameters of transient window range. More...
 
pdf1D_t * XLALComputeTransientPosterior_tau (transientWindowRange_t windowRange, const transientFstatMap_t *FstatMap)
 Compute transient-CW posterior (normalized) on timescale tau, using given type and parameters of transient window range. More...
 
transientFstatMap_tXLALComputeTransientFstatMap (const MultiFstatAtomVector *multiFstatAtoms, transientWindowRange_t windowRange, BOOLEAN useFReg)
 Function to compute transient-window "F-statistic map" over start-time and timescale {t0, tau}. More...
 
FstatAtomVectorXLALmergeMultiFstatAtomsBinned (const MultiFstatAtomVector *multiAtoms, UINT4 deltaT)
 Combine N Fstat-atoms vectors into a single 'canonical' binned and ordered atoms-vector. More...
 
int write_transientCandidate_to_fp (LALFILE *fp, const transientCandidate_t *thisCand, const char timeUnit)
 Write one line for given transient CW candidate into output file. More...
 
int write_transientFstatMap_to_fp (LALFILE *fp, const transientFstatMap_t *FstatMap, const transientWindowRange_t *windowRange, const PulsarDopplerParams *doppler)
 Write full set of t0 and tau grid points (assumed at fixed Doppler parameters) into output file. More...
 
int write_transientCandidateAll_to_fp (LALFILE *fp, const transientCandidate_t *thisCand)
 Write full set of t0 and tau grid points for given transient CW candidate into output file. More...
 
int write_MultiFstatAtoms_to_fp (LALFILE *fp, const MultiFstatAtomVector *multiAtoms)
 Write multi-IFO F-stat atoms 'multiAtoms' into output stream 'fstat'. More...
 
void XLALDestroyTransientFstatMap (transientFstatMap_t *FstatMap)
 Standard destructor for transientFstatMap_t Fully NULL-robust as usual. More...
 
void XLALDestroyTransientCandidate (transientCandidate_t *cand)
 Standard destructor for transientCandidate_t Fully NULL-robust as usual. More...
 
static int XLALCreateExpLUT (void)
 Generate an exponential lookup-table expLUT for e^(-x) over the interval x in [0, xmax], using 'length' points. More...
 
void XLALDestroyExpLUT (void)
 Destructor function for expLUT_t lookup table. More...
 
REAL8 XLALFastNegExp (REAL8 mx)
 Fast exponential function e^-x using lookup-table (LUT). More...
 
static REAL8 XLALGetRectangularTransientWindowValue (UINT4 timestamp, UINT4 t0, UINT4 t1)
 Function to compute the value of a rectangular transient-window at a given timestamp. More...
 
static REAL8 XLALGetExponentialTransientWindowValue (UINT4 timestamp, UINT4 t0, UINT4 t1, UINT4 tau)
 Function to compute the value of an exponential transient-window at a given timestamp. More...
 
static REAL8 XLALGetTransientWindowValue (UINT4 timestamp, UINT4 t0, UINT4 t1, UINT4 tau, transientWindowType_t type)
 Function to compute the value of a given transient-window function at a given timestamp. More...
 

Data Structures

struct  transientWindowRange_t
 Struct defining a range of transient windows. More...
 
struct  transientFstatMap_t
 Struct holding a transient-window "F-statistic map" over start-time and timescale {t0, tau}. More...
 
struct  transientCandidate_t
 Struct holding a transient CW candidate. More...
 

Macros

#define DAY24   (24 * 3600)
 standard 24h day = 86400 seconds ==> this is what's used in the definition of 'tauDays' More...
 
#define TRANSIENT_EXP_EFOLDING   3.0
 e-folding parameter for exponential window, after which we truncate the window for efficiency. More...
 

Function Documentation

◆ XLALParseTransientWindowName()

int XLALParseTransientWindowName ( const char windowName)

Parse a transient window name string into the corresponding transientWindowType.

Definition at line 74 of file TransientCW_utils.c.

◆ XLALGetTransientWindowTimespan()

int XLALGetTransientWindowTimespan ( UINT4 t0,
UINT4 t1,
transientWindow_t  transientWindow 
)

Helper-function to determine the total timespan of a transient CW window, ie.

the earliest and latest timestamps of non-zero window function.

Parameters
[out]t0window start-time
[out]t1window end-time
[in]transientWindowwindow-parameters

Definition at line 110 of file TransientCW_utils.c.

◆ XLALApplyTransientWindow()

int XLALApplyTransientWindow ( REAL4TimeSeries series,
transientWindow_t  transientWindow 
)

apply a "transient CW window" described by TransientWindowParams to the given timeseries

Parameters
seriesinput timeseries to apply window to
transientWindowtransient-CW window to apply

Definition at line 159 of file TransientCW_utils.c.

◆ XLALApplyTransientWindow2NoiseWeights()

int XLALApplyTransientWindow2NoiseWeights ( MultiNoiseWeights multiNoiseWeights,
const MultiLIGOTimeGPSVector multiTS,
transientWindow_t  transientWindow 
)

apply transient window to give multi noise-weights, associated with given multi timestamps

Parameters
multiNoiseWeights[in/out] noise weights to apply transient window to
[in]multiTSassociated timestamps of noise-weights
[in]transientWindowtransient window parameters

Definition at line 222 of file TransientCW_utils.c.

◆ XLALPulsarDopplerParams2String()

CHAR * XLALPulsarDopplerParams2String ( const PulsarDopplerParams par)

Turn pulsar doppler-params into a single string that can be used for filenames The format is tRefNNNNNN_RAXXXXX_DECXXXXXX_FreqXXXXX[_f1dotXXXXX][_f2dotXXXXx][_f3dotXXXXX].

Definition at line 305 of file TransientCW_utils.c.

◆ XLALComputeTransientBstat()

REAL8 XLALComputeTransientBstat ( transientWindowRange_t  windowRange,
const transientFstatMap_t FstatMap 
)

Compute transient-CW Bayes-factor B_SG = P(x|HypS)/P(x|HypG) (where HypG = Gaussian noise hypothesis), marginalized over start-time and timescale of transient CW signal, using given type and parameters of transient window range.

Note: this function is a C-implemention, partially based-on/inspired-by Stefanos Giampanis' original matlab implementation of this search function.

Note2: if window->type == none, uses a single rectangular window covering all the data.

Parameters
[in]windowRangetype and parameters specifying transient window range
[in]FstatMappre-computed transient-Fstat map F_mn over {t0, tau} ranges

Definition at line 367 of file TransientCW_utils.c.

◆ XLALComputeTransientPosterior_t0()

pdf1D_t * XLALComputeTransientPosterior_t0 ( transientWindowRange_t  windowRange,
const transientFstatMap_t FstatMap 
)

Compute transient-CW posterior (normalized) on start-time t0, using given type and parameters of transient window range.

NOTE: the returned pdf has a number of sample-points N_t0Range given by the size of the input matrix FstatMap (namely N_t0Range = t0Band / dt0)

Parameters
[in]windowRangetype and parameters specifying transient window range
[in]FstatMappre-computed transient-Fstat map F_mn over {t0, tau} ranges

Definition at line 432 of file TransientCW_utils.c.

◆ XLALComputeTransientPosterior_tau()

pdf1D_t * XLALComputeTransientPosterior_tau ( transientWindowRange_t  windowRange,
const transientFstatMap_t FstatMap 
)

Compute transient-CW posterior (normalized) on timescale tau, using given type and parameters of transient window range.

NOTE: the returned pdf has a number of sample-points N_tauRange given by the size of the input matrix FstatMap (namely N_tauRange = tauBand / dtau)

Parameters
[in]windowRangetype and parameters specifying transient window range
[in]FstatMappre-computed transient-Fstat map F_mn over {t0, tau} ranges

Definition at line 519 of file TransientCW_utils.c.

◆ XLALComputeTransientFstatMap()

transientFstatMap_t * XLALComputeTransientFstatMap ( const MultiFstatAtomVector multiFstatAtoms,
transientWindowRange_t  windowRange,
BOOLEAN  useFReg 
)

Function to compute transient-window "F-statistic map" over start-time and timescale {t0, tau}.

Returns a 2D matrix F_mn, with m = index over start-times t0, and n = index over timescales tau, in steps of dt0 in [t0, t0+t0Band], and dtau in [tau, tau+tauBand] as defined in transientWindowRange

Note: if window->type == none, we compute a single rectangular window covering all the data.

Note2: if the experimental switch useFReg is true, returns FReg=F - log(D) instead of F. This option is of little practical interest, except for demonstrating that marginalizing (1/D)e^F is less sensitive than marginalizing e^F (see transient methods-paper [in prepartion])

Parameters
[in]multiFstatAtomsmulti-IFO F-statistic atoms
[in]windowRangetype and parameters specifying transient window range to search
[in]useFRegexperimental switch: compute FReg = F - log(D) instead of F

Definition at line 612 of file TransientCW_utils.c.

◆ XLALmergeMultiFstatAtomsBinned()

FstatAtomVector * XLALmergeMultiFstatAtomsBinned ( const MultiFstatAtomVector multiAtoms,
UINT4  deltaT 
)

Combine N Fstat-atoms vectors into a single 'canonical' binned and ordered atoms-vector.

The function pre-sums all atoms on a regular 'grid' of timestep bins deltaT covering the full data-span. Atoms with timestamps falling into the bin i : [t_i, t_{i+1} ) are pre-summed and returned as atoms[i], where t_i = t_0 + i * deltaT.

Note: this pre-binning is equivalent to using a rectangular transient window on the deltaT timescale, which is OK even with a different transient window, provided deltaT << transient-window timescale!

Bins containing no atoms are returned with all values set to zero.

Definition at line 860 of file TransientCW_utils.c.

◆ write_transientCandidate_to_fp()

int write_transientCandidate_to_fp ( LALFILE *  fp,
const transientCandidate_t thisCand,
const char  timeUnit 
)

Write one line for given transient CW candidate into output file.

NOTE: if input thisCand == NULL, we write a header comment-line explaining the fields

Definition at line 943 of file TransientCW_utils.c.

◆ write_transientFstatMap_to_fp()

int write_transientFstatMap_to_fp ( LALFILE *  fp,
const transientFstatMap_t FstatMap,
const transientWindowRange_t windowRange,
const PulsarDopplerParams doppler 
)

Write full set of t0 and tau grid points (assumed at fixed Doppler parameters) into output file.

NOTE: if input FstatMap == NULL, we write a header comment-line explaining the fields

if doppler = NULL, we skip those columns

Definition at line 1001 of file TransientCW_utils.c.

◆ write_transientCandidateAll_to_fp()

int write_transientCandidateAll_to_fp ( LALFILE *  fp,
const transientCandidate_t thisCand 
)

Write full set of t0 and tau grid points for given transient CW candidate into output file.

NOTE: if input thisCand == NULL, we write a header comment-line explaining the fields

Definition at line 1061 of file TransientCW_utils.c.

◆ write_MultiFstatAtoms_to_fp()

int write_MultiFstatAtoms_to_fp ( LALFILE *  fp,
const MultiFstatAtomVector multiAtoms 
)

Write multi-IFO F-stat atoms 'multiAtoms' into output stream 'fstat'.

Definition at line 1089 of file TransientCW_utils.c.

◆ XLALDestroyTransientFstatMap()

void XLALDestroyTransientFstatMap ( transientFstatMap_t FstatMap)

Standard destructor for transientFstatMap_t Fully NULL-robust as usual.

Definition at line 1127 of file TransientCW_utils.c.

◆ XLALDestroyTransientCandidate()

void XLALDestroyTransientCandidate ( transientCandidate_t cand)

Standard destructor for transientCandidate_t Fully NULL-robust as usual.

Definition at line 1148 of file TransientCW_utils.c.

◆ XLALCreateExpLUT()

int XLALCreateExpLUT ( void  )
static

Generate an exponential lookup-table expLUT for e^(-x) over the interval x in [0, xmax], using 'length' points.

Definition at line 1171 of file TransientCW_utils.c.

◆ XLALDestroyExpLUT()

void XLALDestroyExpLUT ( void  )

Destructor function for expLUT_t lookup table.

Definition at line 1202 of file TransientCW_utils.c.

◆ XLALFastNegExp()

REAL8 XLALFastNegExp ( REAL8  mx)

Fast exponential function e^-x using lookup-table (LUT).

We need to compute exp(-x) for x >= 0, typically in a B-stat integral of the form int e^-x dx: this means that small values e^(-x) will not contribute much to the integral and are less important than values close to 1. Therefore we pre-compute a LUT of e^(-x) for x in [0, xmax], in Npoints points, and set e^(-x) = 0 for x < xmax.

NOTE: if module-global expLUT=NULL, we create it here NOTE: if argument is negative, we use math-lib exp(-x) instead of LUT

Definition at line 1229 of file TransientCW_utils.c.

◆ XLALGetRectangularTransientWindowValue()

static REAL8 XLALGetRectangularTransientWindowValue ( UINT4  timestamp,
UINT4  t0,
UINT4  t1 
)
inlinestatic

Function to compute the value of a rectangular transient-window at a given timestamp.

This is the central function defining the rectangular window properties.

Parameters
timestamptimestamp for which to compute window-value
t0start-time of rectangular window
t1end-time of rectangular window

Definition at line 156 of file TransientCW_utils.h.

◆ XLALGetExponentialTransientWindowValue()

static REAL8 XLALGetExponentialTransientWindowValue ( UINT4  timestamp,
UINT4  t0,
UINT4  t1,
UINT4  tau 
)
inlinestatic

Function to compute the value of an exponential transient-window at a given timestamp.

This is the central function defining the exponential window properties.

Parameters
timestamptimestamp for which to compute window-value
t0start-time of exponential window
t1end-time of exponential window
taucharacteristic time of the exponential window

Definition at line 175 of file TransientCW_utils.h.

◆ XLALGetTransientWindowValue()

static REAL8 XLALGetTransientWindowValue ( UINT4  timestamp,
UINT4  t0,
UINT4  t1,
UINT4  tau,
transientWindowType_t  type 
)
inlinestatic

Function to compute the value of a given transient-window function at a given timestamp.

This is a simple wrapper to the actual window-defining functions

Parameters
timestamptimestamp for which to compute window-value
t0start-time of window
t1end-time of window
taucharacteristic time of window
typewindow type

Definition at line 200 of file TransientCW_utils.h.

Macro Definition Documentation

◆ DAY24

#define DAY24   (24 * 3600)

standard 24h day = 86400 seconds ==> this is what's used in the definition of 'tauDays'

Definition at line 63 of file TransientCW_utils.h.

◆ TRANSIENT_EXP_EFOLDING

#define TRANSIENT_EXP_EFOLDING   3.0

e-folding parameter for exponential window, after which we truncate the window for efficiency.

3 e-foldings means we lose only about e^(-2x3) ~1e-8 of signal power!

Definition at line 70 of file TransientCW_utils.h.