LALPulsar  6.1.0.1-89842e6
PulsarCrossCorr_v2.h File Reference

Prototypes

int XLALGetDopplerShiftedFrequencyInfo (REAL8Vector *shiftedFreqs, UINT4Vector *lowestBins, COMPLEX8Vector *expSignalPhases, REAL8VectorSequence *sincList, UINT4 numBins, PulsarDopplerParams *dopp, SFTIndexList *sftIndices, MultiSFTVector *inputSFTs, MultiSSBtimes *multiTimes, MultiUINT4Vector *badBins, REAL8 Tsft)
 Calculate the Doppler-shifted frequency associated with each SFT in a list. More...
 
int XLALCreateSFTIndexListFromMultiSFTVect (SFTIndexList **indexList, MultiSFTVector *sfts)
 Construct flat SFTIndexList out of a MultiSFTVector. More...
 
int XLALCreateSFTPairIndexList (SFTPairIndexList **pairIndexList, SFTIndexList *indexList, MultiSFTVector *sfts, REAL8 maxLag, BOOLEAN inclAutoCorr)
 Construct list of SFT pairs for inclusion in statistic. More...
 
int XLALCreateSFTPairIndexListResamp (MultiResampSFTPairMultiIndexList **resampMultiPairIndexList, SFTPairIndexList **pairIndexList, SFTIndexList *indexList, MultiSFTVector *sfts, REAL8 maxLag, BOOLEAN inclAutoCorr, BOOLEAN inclSameDetector, REAL8 Tsft, REAL8 Tshort)
 Resampling-modified: construct list of SFT pairs for inclusion in statistic. More...
 
int XLALCreateSFTPairIndexListShortResamp (MultiResampSFTPairMultiIndexList **resampPairIndexList, const REAL8 maxLag, const BOOLEAN inclAutoCorr, const BOOLEAN inclSameDetector, const REAL8 Tsft, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes)
 
int XLALTestResampPairIndexList (MultiResampSFTPairMultiIndexList *resampMultiPairIndexList)
 Check that the contents of a resampling multi-pair index list are sensible by inspection. More...
 
int XLALCalculateCrossCorrGammas (REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, SFTPairIndexList *pairIndexList, SFTIndexList *indexList, MultiAMCoeffs *multiCoeffs)
 Construct vector of G_alpha amplitudes for each SFT pair. More...
 
int XLALCalculateCrossCorrGammasResamp (REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, MultiResampSFTPairMultiIndexList *resampMultiPairIndexList, MultiAMCoeffs *multiCoeffs)
 test function for RESAMPLING More...
 
int XLALCalculateCrossCorrGammasResampShort (REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, MultiResampSFTPairMultiIndexList *resampMultiPairIndexList, MultiAMCoeffs *multiCoeffs)
 test function for RESAMPLING with tShort More...
 
int XLALCalculatePulsarCrossCorrStatistic (REAL8 *ccStat, REAL8 *evSquared, REAL8Vector *curlyGAmp, COMPLEX8Vector *expSignalPhases, UINT4Vector *lowestBins, REAL8VectorSequence *sincList, SFTPairIndexList *sftPairs, SFTIndexList *sftIndices, MultiSFTVector *inputSFTs, MultiNoiseWeights *multiWeights, UINT4 numBins)
 Calculate multi-bin cross-correlation statistic. More...
 
int XLALCalculatePulsarCrossCorrStatisticResamp (REAL8Vector *_LAL_RESTRICT_ ccStatVector, REAL8Vector *_LAL_RESTRICT_ evSquaredVector, REAL8Vector *_LAL_RESTRICT_ numeEquivAve, REAL8Vector *_LAL_RESTRICT_ numeEquivCirc, const REAL8Vector *_LAL_RESTRICT_ resampCurlyGAmp, const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampPairs, const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights, const PulsarDopplerParams *_LAL_RESTRICT_ binaryTemplateSpacings, const PulsarDopplerParams *_LAL_RESTRICT_ dopplerpos, const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_a, const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_b, ResampCrossCorrWorkspace *_LAL_RESTRICT_ ws, COMPLEX8 *_LAL_RESTRICT_ ws1KFaX_k, COMPLEX8 *_LAL_RESTRICT_ ws1KFbX_k, COMPLEX8 *_LAL_RESTRICT_ ws2LFaX_k, COMPLEX8 *_LAL_RESTRICT_ ws2LFbX_k)
 
int XLALCalculateCrossCorrPhaseDerivatives (REAL8VectorSequence **phaseDerivs, const PulsarDopplerParams *dopplerPoint, const EphemerisData *edat, SFTIndexList *indexList, MultiSSBtimes *multiTimes, const DopplerCoordinateSystem *coordSys)
 calculate signal phase derivatives wrt Doppler coords, for each SFT More...
 
int XLALCalculateCrossCorrPhaseDerivativesShort (REAL8VectorSequence **resampPhaseDerivs, const PulsarDopplerParams *dopplerPoint, const EphemerisData *edat, SFTIndexList *indexList, MultiResampSFTPairMultiIndexList *resampMultiPairs, MultiSSBtimes *multiTimes, const DopplerCoordinateSystem *coordSys)
 (test function) MODIFIED for Tshort More...
 
int XLALCalculateCrossCorrPhaseMetric (gsl_matrix **g_ij, gsl_vector **eps_i, REAL8 *sumGammaSq, const REAL8VectorSequence *phaseDerivs, const SFTPairIndexList *pairIndexList, const REAL8Vector *Gamma_ave, const REAL8Vector *Gamma_circ, const DopplerCoordinateSystem *coordSys)
 calculate phase metric for CW cross-correlation search, as well as vector used for parameter offsets More...
 
int XLALCalculateCrossCorrPhaseMetricShort (gsl_matrix **g_ij, gsl_vector **eps_i, REAL8 *sumGammaSq, const REAL8VectorSequence *phaseDerivs, const MultiResampSFTPairMultiIndexList *resampMultiPairs, const REAL8Vector *Gamma_ave, const REAL8Vector *Gamma_circ, const DopplerCoordinateSystem *coordSys)
 (test function) Redesigning to use Tshort instead More...
 
int XLALCalculateLMXBCrossCorrDiagMetric (REAL8 *hSens, REAL8 *g_ff, REAL8 *g_aa, REAL8 *g_TT, REAL8 *g_pp, REAL8 *weightedMuTAve, PulsarDopplerParams DopplerParams, REAL8Vector *G_alpha, SFTPairIndexList *pairIndexList, SFTIndexList *indexList, MultiSFTVector *sfts, MultiNoiseWeights *multiWeights)
 
int XLALCalculateLMXBCrossCorrDiagMetricShort (REAL8 *hSens, REAL8 *g_ff, REAL8 *g_aa, REAL8 *g_TT, REAL8 *g_pp, const PulsarDopplerParams DopplerParams, const REAL8Vector *_LAL_RESTRICT_ G_alpha, const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampMultiPairIndexList, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights)
 
LIGOTimeGPSVectorXLALModifyTimestampsFromSFTsShort (REAL8TimeSeries **sciFlag, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times, const REAL8 tShort, const UINT4 numShortPerDet)
 
LIGOTimeGPSVectorXLALExtractTimestampsFromSFTsShort (REAL8TimeSeries **sciFlag, const SFTVector *sfts, REAL8 tShort, UINT4 numShortPerDet)
 ‍** (possible future function) Wrapper for Bessel Orbital Space stepper in the manner of V. Dergachev's loosely-coherent search *‍/ More...
 
MultiLIGOTimeGPSVectorXLALModifyMultiTimestampsFromSFTs (MultiREAL8TimeSeries **scienceFlagVect, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet)
 
MultiLIGOTimeGPSVectorXLALExtractMultiTimestampsFromSFTsShort (MultiREAL8TimeSeries **scienceFlagVect, const MultiSFTVector *multiSFTs, REAL8 tShort, UINT4 numShortPerDet)
 Given a multi-SFT vector, return a MultiLIGOTimeGPSVector holding the SFT timestamps (test function using tShort) More...
 
int XLALFillDetectorTensorShort (DetectorState *detState, const LALDetector *detector)
 (test function) fill detector state with tShort, importing various slightly-modified LALPulsar functions for testing More...
 
DetectorStateSeriesXLALGetDetectorStatesShort (const LIGOTimeGPSVector *timestamps, const LALDetector *detector, const EphemerisData *edat, REAL8 tOffset, REAL8 tShort, UINT4 numShortPerDet)
 (test function) get detector states for tShort More...
 
MultiDetectorStateSeriesXLALGetMultiDetectorStatesShort (const MultiLIGOTimeGPSVector *multiTS, const MultiLALDetector *multiIFO, const EphemerisData *edat, REAL8 tOffset, REAL8 tShort, UINT4 numShortPerDet)
 (test function) get multi detector states for tShort More...
 
int XLALModifyAMCoeffsWeights (REAL8Vector **resampMultiWeightsX, const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const UINT4 maxNumStepsOldIfGapless, const UINT4 X)
 
int XLALModifyMultiAMCoeffsWeights (MultiNoiseWeights **multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes)
 
int XLALWeightMultiAMCoeffsShort (MultiAMCoeffs *multiAMcoef, const MultiNoiseWeights *multiWeights, REAL8 tShort, REAL8 tSFTOld, UINT4 numShortPerDet, MultiLIGOTimeGPSVector *multiTimes)
 (test function) used for weighting multi amplitude modulation coefficients More...
 
MultiAMCoeffsXLALComputeMultiAMCoeffsShort (const MultiDetectorStateSeries *multiDetStates, const MultiNoiseWeights *multiWeights, SkyPosition skypos, REAL8 tShort, REAL8 tSFTOld, UINT4 numShortPerDet, MultiLIGOTimeGPSVector *multiTimes)
 (test function) used for computing multi amplitude modulation weights More...
 
AMCoeffsXLALComputeAMCoeffsShort (const DetectorStateSeries *DetectorStates, SkyPosition skypos)
 (test function) used for computing amplitude modulation weights More...
 
UINT4 XLALCrossCorrNumShortPerDetector (const REAL8 resampTshort, const INT4 startTime, const INT4 endTime)
 Compute the number of tShort segments per detector. More...
 
REAL8TimeSeriesXLALCrossCorrGapFinderResamp (LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times)
 
REAL8TimeSeriesXLALCrossCorrGapFinderResampAlt (LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const SFTVector *_LAL_RESTRICT_ sfts)
 
int XLALEquipCrossCorrPairsWithScienceFlags (MultiResampSFTPairMultiIndexList *resampMultiPairs, MultiREAL8TimeSeries *scienceFlagVect)
 Demarcate pairs with flags about whether data exists in zero-padded timeseries. More...
 
int XLALCreateCrossCorrWorkspace (ResampCrossCorrWorkspace **wsOut, COMPLEX8 **ws1KFaX_kOut, COMPLEX8 **ws1KFbX_kOut, COMPLEX8 **ws2LFaX_kOut, COMPLEX8 **ws2LFbX_kOut, MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_aOut, MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_bOut, const PulsarDopplerParams binaryTemplateSpacings, FstatInput *resampFstatInput, const UINT4 numFreqBins, const REAL8 tCoh, const BOOLEAN treatWarningsAsErrors)
 ‍** (possible future function) does not work – would adjust timestamps of an SFT vector *‍/ More...
 
void XLALDestroySFTIndexList (SFTIndexList *sftIndices)
 Destroy a SFTIndexList structure. More...
 
void XLALDestroySFTPairIndexList (SFTPairIndexList *sftPairs)
 Destroy a SFTPairIndexList structure. More...
 
void XLALDestroyResampSFTIndexList (ResampSFTIndexList *sftResampList)
 
void XLALDestroyResampSFTMultiIndexList (ResampSFTMultiIndexList *sftResampMultiList)
 
void XLALDestroyResampSFTPairMultiIndexList (ResampSFTPairMultiIndexList *sftResampPairMultiList)
 
void XLALDestroyMultiResampSFTPairMultiIndexList (MultiResampSFTPairMultiIndexList *sftMultiPairsResamp)
 
void XLALDestroyMultiMatchList (MultiResampSFTMultiCountList *localMultiListOfLmatchingGivenMultiK)
 
void XLALDestroyResampCrossCorrWorkspace (void *workspace)
 

Go to the source code of this file.

Data Structures

struct  SFTIndex
 Index to refer to an SFT given a set of SFTs from several different detectors. More...
 
struct  SFTIndexList
 List of SFT indices. More...
 
struct  SFTPairIndex
 Index to refer to a pair of SFTs. More...
 
struct  SFTPairIndexList
 List of SFT pair indices. More...
 
struct  SFTCount
 Resampling Counter of matching SFTs for a given detector Y_K_X matching SFT K_X. More...
 
struct  ResampSFTMultiCountList
 INNER List of SFT indices. More...
 
struct  ResampSFTMultiCountListDet
 MIDDLE List of SFT indices. More...
 
struct  MultiResampSFTMultiCountList
 OUTER List of SFT indices. More...
 
struct  ResampSFTIndex
 Resampling: Index to refer to fields of an SFT given a specific index L_Y_K_X. More...
 
struct  ResampSFTIndexList
 Resampling: List of SFT indices L for a given detector Y_K_X: indexing method is nominally original vectors but may be affected by gaps. More...
 
struct  ResampSFTMultiIndexList
 Resampling: Multi List of indices of SFT L_Y, for a given sft K_X
More...
 
struct  ResampSFTPairMultiIndexList
 Resampling Multi List of SFT pair indices (L_Y_K), for a given detector X. More...
 
struct  MultiResampSFTPairMultiIndexList
 Resampling Multi List of all paired SFTs L_Y_K_X, top level (multiple pairs, multiple detectors) More...
 
struct  CrossCorrTimings_t
 
struct  ResampCrossCorrTimingInfo
 
struct  ResampCrossCorrWorkspace
 
struct  MultiUINT4Vector
 A collection of UINT4Vectors – one for each IFO
More...