Header-file for XLAL routines for v2 CW cross-correlation searches.
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 | XLALCreateSFTPairIndexListAccurateResamp (SFTPairIndexList **pairIndexList, UINT4VectorSequence **sftPairForTshortPair, SFTIndexList *indexList, MultiResampSFTPairMultiIndexList *resampPairs, const MultiLIGOTimeGPSVector *restrict multiTimes, const MultiLIGOTimeGPSVector *restrict resampMultiTimes) |
| Construct list of SFT pairs corresponding to a list of tShort pairs. 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 **resampMultiPairIndexList, const REAL8 maxLag, const BOOLEAN inclAutoCorr, const BOOLEAN inclSameDetector, const REAL8 Tsft, const MultiLIGOTimeGPSVector *restrict multiTimes) |
| With Tshort, and Resampling-modified: construct list of SFT pairs for inclusion in statistic. More... | |
| 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 | XLALCombineCrossCorrGammas (REAL8Vector **resampGamma, REAL8Vector *Gamma, UINT4VectorSequence *sftPairForTshortPair, REAL8 Tsft, REAL8 Tshort) |
| Combine G_alpha amplitudes for SFT pairs into amplitudes for Tshort pairs. 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 *restrict ccStatVector, REAL8Vector *restrict evSquaredVector, REAL8Vector *restrict numeEquivAve, REAL8Vector *restrict numeEquivCirc, const REAL8Vector *restrict resampCurlyGAmp, const MultiResampSFTPairMultiIndexList *restrict resampMultiPairs, const MultiNoiseWeights *restrict multiWeights, const PulsarDopplerParams *restrict binaryTemplateSpacings, const PulsarDopplerParams *restrict dopplerpos, const MultiCOMPLEX8TimeSeries *restrict multiTimeSeries_SRC_a, const MultiCOMPLEX8TimeSeries *restrict multiTimeSeries_SRC_b, ResampCrossCorrWorkspace *restrict ws, COMPLEX8 *restrict ws1KFaX_k, COMPLEX8 *restrict ws1KFbX_k, COMPLEX8 *restrict ws2LFaX_k, COMPLEX8 *restrict ws2LFbX_k) |
| Calculate multi-bin cross-correlation statistic using resampling. More... | |
| 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 *restrict G_alpha, const MultiResampSFTPairMultiIndexList *restrict resampMultiPairIndexList, const MultiLIGOTimeGPSVector *restrict timestamps, const MultiNoiseWeights *restrict multiWeights) |
| MODIFIED for Tshort: calculate metric diagonal components, also include the estimation of sensitivity E[rho]/(h_0)^2. More... | |
| LIGOTimeGPSVector * | XLALExtractTimestampsFromSFTsShort (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... | |
| LIGOTimeGPSVector * | XLALGenerateTshortTimestamps (const REAL8 tShort, const UINT4 numShortPerDet, const LIGOTimeGPS epoch) |
| Generate timestamps for one detector with tShort. More... | |
| LIGOTimeGPSVector * | XLALModifyTimestampsFromSFTsShort (REAL8TimeSeries **sciFlag, const LIGOTimeGPSVector *restrict Times, const REAL8 tShort, const UINT4 numShortPerDet) |
| Modify timestamps from one detector with tShort. More... | |
| MultiLIGOTimeGPSVector * | XLALGenerateMultiTshortTimestamps (const MultiLIGOTimeGPSVector *restrict multiTimes, const REAL8 tShort, const UINT4 numShortPerDet, const BOOLEAN alignTShorts) |
| Given a multi-SFT vector, return a MultiLIGOTimeGPSVector holding the modified SFT timestamps (production function using tShort) More... | |
| MultiLIGOTimeGPSVector * | XLALModifyMultiTimestampsFromSFTs (MultiREAL8TimeSeries **scienceFlagVect, const MultiLIGOTimeGPSVector *restrict multiTimes, const REAL8 tShort, const UINT4 numShortPerDet) |
| Given a multi-SFT vector, return a MultiLIGOTimeGPSVector holding the modified SFT timestamps (production function using tShort) More... | |
| MultiLIGOTimeGPSVector * | XLALExtractMultiTimestampsFromSFTsShort (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... | |
| DetectorStateSeries * | XLALGetDetectorStatesShort (const LIGOTimeGPSVector *timestamps, const LALDetector *detector, const EphemerisData *edat, REAL8 tOffset, REAL8 tShort, UINT4 numShortPerDet) |
| (test function) get detector states for tShort More... | |
| MultiDetectorStateSeries * | XLALGetMultiDetectorStatesShort (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 *restrict multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *restrict multiTimes, const UINT4 maxNumStepsOldIfGapless, const UINT4 X) |
| Modify amplitude weight coefficients for tShort. More... | |
| MultiNoiseWeights * | XLALModifyMultiWeights (const MultiNoiseWeights *restrict multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *restrict multiTimes) |
| Modify multiple detectors' amplitude weight coefficients for tShort. More... | |
| int | XLALModifyMultiAMCoeffsWeights (MultiNoiseWeights **multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *restrict multiTimes) |
| Modify multiple detectors' amplitude weight coefficients for tShort. More... | |
| 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... | |
| AMCoeffs * | XLALComputeAMCoeffsShort (const DetectorStateSeries *DetectorStates, SkyPosition skypos) |
| (test function) used for computing amplitude modulation weights More... | |
| MultiAMCoeffs * | XLALComputeMultiAMCoeffsShort (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... | |
| UINT4 | XLALCrossCorrNumShortPerDetector (const REAL8 resampTshort, const INT4 startTime, const INT4 endTime) |
| Compute the number of tShort segments per detector. More... | |
| REAL8TimeSeries * | XLALCrossCorrGapFinderResamp (LIGOTimeGPSVector *restrict timestamps, const LIGOTimeGPSVector *restrict Times) |
| Find gaps in the data given the SFTs. More... | |
| REAL8TimeSeries * | XLALCrossCorrGapFinderResampAlt (LIGOTimeGPSVector *restrict timestamps, const SFTVector *restrict sfts) |
| (test function) find gaps in the data given the SFTs More... | |
| 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) |
| static int | XLALApplyCrossCorrFreqShiftResamp (COMPLEX8 *restrict xOut, const COMPLEX8TimeSeries *restrict xIn, const PulsarDopplerParams *restrict doppler, const REAL8 freqShift, const UINT4 indexStartResamp, const UINT4 indexEndResamp, const UINT4 numSamplesIn, const UINT4 insertPoint) |
| Imported and modified from ComputeFstat_Resamp.c. More... | |
| static int | XLALComputeFaFb_CrossCorrResamp (ResampCrossCorrWorkspace *restrict ws, COMPLEX8 *restrict wsFaX_k, COMPLEX8 *restrict wsFbX_k, const MultiResampSFTPairMultiIndexList *restrict resampMultiPairs, const MultiCOMPLEX8TimeSeries *restrict multiTimeSeries_SRC_a, const MultiCOMPLEX8TimeSeries *restrict multiTimeSeries_SRC_b, const PulsarDopplerParams *restrict dopplerpos, const PulsarDopplerParams *restrict binaryTemplateSpacings, const REAL8 SRCsampPerTcoh, const UINT4 detX, const UINT4 sftK, const UINT4 detY, const BOOLEAN isL) |
| Compute the equivalent of Fa and Fb for CrossCorr's rho statistic. More... | |
| int | XLALCreateSFTPairIndexListAccurateResamp (SFTPairIndexList **pairIndexList, UINT4VectorSequence **sftPairForResampPair, SFTIndexList *indexList, MultiResampSFTPairMultiIndexList *resampPairs, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ resampMultiTimes) |
| int | XLALCreateSFTPairIndexListShortResamp (MultiResampSFTPairMultiIndexList **resampPairIndexList, const REAL8 maxLag, const BOOLEAN inclAutoCorr, const BOOLEAN inclSameDetector, const REAL8 Tsft, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes) |
| 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 | 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) |
| LIGOTimeGPSVector * | XLALModifyTimestampsFromSFTsShort (REAL8TimeSeries **sciFlag, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times, const REAL8 tShort, const UINT4 numShortPerDet) |
| MultiLIGOTimeGPSVector * | XLALGenerateMultiTshortTimestamps (const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet, const BOOLEAN alignTShorts) |
| MultiLIGOTimeGPSVector * | XLALModifyMultiTimestampsFromSFTs (MultiREAL8TimeSeries **scienceFlagVect, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet) |
| 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) |
| MultiNoiseWeights * | XLALModifyMultiWeights (const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes) |
| int | XLALModifyMultiAMCoeffsWeights (MultiNoiseWeights **multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes) |
| REAL8TimeSeries * | XLALCrossCorrGapFinderResamp (LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times) |
| REAL8TimeSeries * | XLALCrossCorrGapFinderResampAlt (LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const SFTVector *_LAL_RESTRICT_ sfts) |
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... | |
| 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.
| shiftedFreqs | Output list of shifted frequencies |
| lowestBins | Output list of bin indices |
| expSignalPhases | Output list of signal phases |
| sincList | Output list of sinc factors |
| numBins | Number of frequency bins to use |
| dopp | Doppler parameters for signal |
| sftIndices | List of indices for SFTs |
| inputSFTs | SFT data (needed for f0) |
| multiTimes | SSB or Binary times |
| badBins | List of bin indices contaminated by known lines |
| Tsft | SFT duration |
Definition at line 74 of file PulsarCrossCorr_v2.c.
| int XLALCreateSFTIndexListFromMultiSFTVect | ( | SFTIndexList ** | indexList, |
| MultiSFTVector * | sfts | ||
| ) |
Construct flat SFTIndexList out of a MultiSFTVector.
| indexList | Output: flat list of indices to locate SFTs |
| sfts | Input: set of per-detector SFT vectors |
Definition at line 197 of file PulsarCrossCorr_v2.c.
| int XLALCreateSFTPairIndexList | ( | SFTPairIndexList ** | pairIndexList, |
| SFTIndexList * | indexList, | ||
| MultiSFTVector * | sfts, | ||
| REAL8 | maxLag, | ||
| BOOLEAN | inclAutoCorr | ||
| ) |
Construct list of SFT pairs for inclusion in statistic.
| pairIndexList | Output: list of SFT pairs |
| indexList | Input: list of indices to locate SFTs |
| sfts | Input: set of per-detector SFT vectors |
| maxLag | Maximum allowed lag time |
| inclAutoCorr | Flag indicating whether a "pair" of an SFT with itself is allowed |
Definition at line 241 of file PulsarCrossCorr_v2.c.
| int XLALCreateSFTPairIndexListAccurateResamp | ( | SFTPairIndexList ** | pairIndexList, |
| UINT4VectorSequence ** | sftPairForTshortPair, | ||
| SFTIndexList * | indexList, | ||
| MultiResampSFTPairMultiIndexList * | resampPairs, | ||
| const MultiLIGOTimeGPSVector *restrict | multiTimes, | ||
| const MultiLIGOTimeGPSVector *restrict | resampMultiTimes | ||
| ) |
Construct list of SFT pairs corresponding to a list of tShort pairs.
| pairIndexList | Output: list of SFT pairs |
| sftPairForTshortPair | Output: list of SFT pairs corresponding to each Tshort pair |
| indexList | Input: list of indices to locate SFTs |
| resampPairs | Input: pairs of tShort segments for resampling |
| multiTimes | Input: timestamps containing SFT times for each detector |
| resampMultiTimes | Input: timestamps containing Tshort times for each detector |
Definition at line 313 of file PulsarCrossCorr_v2.c.
| 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.
| [out] | resampMultiPairIndexList | resamp list of SFT pairs |
| [out] | pairIndexList | list of SFT pairs |
| [in] | indexList | list of indices to locate SFTs |
| [in] | sfts | set of per-detector SFT vectors |
| [in] | maxLag | Maximum allowed lag time |
| [in] | inclAutoCorr | Flag indicating whether a "pair" of an SFT with itself is allowed |
| [in] | inclSameDetector | Flag indicating whether a "pair" of a detector with itself is allowed |
| [in] | Tsft | duration of a single SFT |
| [in] | Tshort | resampling Tshort |
Definition at line 436 of file PulsarCrossCorr_v2.c.
| int XLALCreateSFTPairIndexListShortResamp | ( | MultiResampSFTPairMultiIndexList ** | resampMultiPairIndexList, |
| const REAL8 | maxLag, | ||
| const BOOLEAN | inclAutoCorr, | ||
| const BOOLEAN | inclSameDetector, | ||
| const REAL8 | Tsft, | ||
| const MultiLIGOTimeGPSVector *restrict | multiTimes | ||
| ) |
With Tshort, and Resampling-modified: construct list of SFT pairs for inclusion in statistic.
| [out] | resampMultiPairIndexList | resamp list of SFT pairs |
| [in] | maxLag | Maximum allowed lag time |
| [in] | inclAutoCorr | Flag indicating whether a "pair" of an SFT with itself is allowed |
| [in] | inclSameDetector | Flag indicating whether a "pair" of a detector with itself is allowed |
| [in] | Tsft | duration of a single SFT |
| [in] | multiTimes | timestamps containing Tshort times for each detector |
Definition at line 690 of file PulsarCrossCorr_v2.c.
| int XLALTestResampPairIndexList | ( | MultiResampSFTPairMultiIndexList * | resampMultiPairIndexList | ) |
Check that the contents of a resampling multi-pair index list are sensible by inspection.
| [in] | resampMultiPairIndexList | resamp list of SFT pairs |
Definition at line 950 of file PulsarCrossCorr_v2.c.
| int XLALCalculateCrossCorrGammas | ( | REAL8Vector ** | Gamma_ave, |
| REAL8Vector ** | Gamma_circ, | ||
| SFTPairIndexList * | pairIndexList, | ||
| SFTIndexList * | indexList, | ||
| MultiAMCoeffs * | multiCoeffs | ||
| ) |
Construct vector of G_alpha amplitudes for each SFT pair.
| Gamma_ave | Output: vector of aa+bb values |
| Gamma_circ | Output: vector of ab-ba values |
| pairIndexList | Input: list of SFT pairs |
| indexList | Input: list of SFTs |
| multiCoeffs | Input: AM coefficients |
Definition at line 989 of file PulsarCrossCorr_v2.c.
| int XLALCalculateCrossCorrGammasResamp | ( | REAL8Vector ** | Gamma_ave, |
| REAL8Vector ** | Gamma_circ, | ||
| MultiResampSFTPairMultiIndexList * | resampMultiPairIndexList, | ||
| MultiAMCoeffs * | multiCoeffs | ||
| ) |
test function for RESAMPLING
Construct multi-dimensional array of G_L_Y_K_X amplitudes for each SFT pair
| Gamma_ave | Output: vector of aa+bb values |
| Gamma_circ | Output: vector of ab-ba values |
| resampMultiPairIndexList | Input: resamp list of SFT pairs |
| multiCoeffs | Input: AM coefficients |
Definition at line 1034 of file PulsarCrossCorr_v2.c.
| int XLALCalculateCrossCorrGammasResampShort | ( | REAL8Vector ** | Gamma_ave, |
| REAL8Vector ** | Gamma_circ, | ||
| MultiResampSFTPairMultiIndexList * | resampMultiPairIndexList, | ||
| MultiAMCoeffs * | multiCoeffs | ||
| ) |
test function for RESAMPLING with tShort
Construct multi-dimensional array of G_L_Y_K_X amplitudes for each SFT pair
| Gamma_ave | Output: vector of aa+bb values |
| Gamma_circ | Output: vector of ab-ba values |
| resampMultiPairIndexList | Input: resamp list of SFT pairs |
| multiCoeffs | Input: AM coefficients |
Definition at line 1093 of file PulsarCrossCorr_v2.c.
| int XLALCombineCrossCorrGammas | ( | REAL8Vector ** | resampGamma, |
| REAL8Vector * | Gamma, | ||
| UINT4VectorSequence * | sftPairForTshortPair, | ||
| REAL8 | Tsft, | ||
| REAL8 | Tshort | ||
| ) |
Combine G_alpha amplitudes for SFT pairs into amplitudes for Tshort pairs.
| resampGamma | Output: vector of combined values, one per Tshort pair | |
| Gamma | Input: vector of values, on per SFT pair | |
| sftPairForTshortPair | Input: list of SFT pair indices for each Tshort pair | |
| [in] | Tsft | duration of a single SFT |
| [in] | Tshort | resampling Tshort |
Definition at line 1179 of file PulsarCrossCorr_v2.c.
| 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.
multiWeights->data[detInd1]->data[sftNum1] * multiWeights->data[detInd2]->data[sftNum2]
| ccStat | Output: cross-correlation statistic rho |
| evSquared | Output: (E[rho]/h0^2)^2 |
| curlyGAmp | Input: Amplitude of curly G for each pair |
| expSignalPhases | Input: Phase of signal for each SFT |
| lowestBins | Input: Bin index to start with for each SFT |
| sincList | Input: input the sinc factors |
| sftPairs | Input: flat list of SFT pairs |
| sftIndices | Input: flat list of SFTs |
| inputSFTs | Input: SFT data |
| multiWeights | Input: nomalizeation factor S^-1 & weights for each SFT |
| numBins | Input Number of frequency bins to be taken into calc |
Definition at line 1221 of file PulsarCrossCorr_v2.c.
| int XLALCalculatePulsarCrossCorrStatisticResamp | ( | REAL8Vector *restrict | ccStatVector, |
| REAL8Vector *restrict | evSquaredVector, | ||
| REAL8Vector *restrict | numeEquivAve, | ||
| REAL8Vector *restrict | numeEquivCirc, | ||
| const REAL8Vector *restrict | resampCurlyGAmp, | ||
| const MultiResampSFTPairMultiIndexList *restrict | resampMultiPairs, | ||
| const MultiNoiseWeights *restrict | multiWeights, | ||
| const PulsarDopplerParams *restrict | binaryTemplateSpacings, | ||
| const PulsarDopplerParams *restrict | dopplerpos, | ||
| const MultiCOMPLEX8TimeSeries *restrict | multiTimeSeries_SRC_a, | ||
| const MultiCOMPLEX8TimeSeries *restrict | multiTimeSeries_SRC_b, | ||
| ResampCrossCorrWorkspace *restrict | ws, | ||
| COMPLEX8 *restrict | ws1KFaX_k, | ||
| COMPLEX8 *restrict | ws1KFbX_k, | ||
| COMPLEX8 *restrict | ws2LFaX_k, | ||
| COMPLEX8 *restrict | ws2LFbX_k | ||
| ) |
Calculate multi-bin cross-correlation statistic using resampling.
| [out] | ccStatVector | vector cross-correlation statistic rho |
| [out] | evSquaredVector | vector (E[rho]/h0^2) |
| [out] | numeEquivAve | vector for intermediate average statistic |
| [out] | numeEquivCirc | vector for intermediate circular statistic |
| [in] | resampCurlyGAmp | amplitude of curly G for each L_Y_K_X (alpha-indexed) |
| [in] | resampMultiPairs | resamp multi list of SFT pairs |
| [in] | multiWeights | normalization factor S^-1 & weights for each SFT |
| [in] | binaryTemplateSpacings | Set of spacings for search |
| [in] | dopplerpos | Doppler point to search |
| [in] | multiTimeSeries_SRC_a | resampled time series A |
| [in] | multiTimeSeries_SRC_b | resampled time series B |
| ws | [in/out] first workspace | |
| ws1KFaX_k | [in/out] holder for detector 1 Fa | |
| ws1KFbX_k | [in/out] holder for detector 1 Fb | |
| ws2LFaX_k | [in/out] holder for detector 2 Fa | |
| ws2LFbX_k | [in/out] holder for detector 2 Fb |
Definition at line 1338 of file PulsarCrossCorr_v2.c.
| 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
| phaseDerivs | Output: dPhi_K/dlambda_i; i is the "sequence" index, K is the "vector" index |
| dopplerPoint | Input: pulsar/binary orbit paramaters |
| edat | Input: Earth/Sun ephemeris |
| indexList | Input: list of SFT indices |
| multiTimes | Input: barycentered times of SFTs |
| coordSys | Input: coordinates with which to differentiate |
Definition at line 1433 of file PulsarCrossCorr_v2.c.
| int XLALCalculateCrossCorrPhaseDerivativesShort | ( | REAL8VectorSequence ** | resampPhaseDerivs, |
| const PulsarDopplerParams * | dopplerPoint, | ||
| const EphemerisData * | edat, | ||
| SFTIndexList * | indexList, | ||
| MultiResampSFTPairMultiIndexList * | resampMultiPairs, | ||
| MultiSSBtimes * | multiTimes, | ||
| const DopplerCoordinateSystem * | coordSys | ||
| ) |
(test function) MODIFIED for Tshort
calculate signal phase derivatives wrt Doppler coords, for each SFT
| [out] | resampPhaseDerivs | dPhi_K/dlambda_i; i is the "sequence" index, K is the "vector" index |
| [in] | dopplerPoint | pulsar/binary orbit paramaters |
| edat | [in Earth/Sun ephemeris | |
| [in] | indexList | list of SFT indices |
| [in] | resampMultiPairs | resamp list of SFT pairs |
| [in] | multiTimes | barycentered times of SFTs |
| [in] | coordSys | coordinates with which to differentiate |
Definition at line 1481 of file PulsarCrossCorr_v2.c.
| 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
This calculates the metric defined in (4.7) of Whelan et al 2015 and the parameter offset epsilon_i (not including the cosi-dependent prefactor) in defined in (4.8)
| g_ij | Output: parameter space metric |
| eps_i | Output: parameter offset vector from (4.8) of WSZP15 |
| sumGammaSq | Output: sum of (Gamma_ave)^2 for normalization and sensitivity |
| phaseDerivs | Input: dPhi_K/dlambda_i; i is the "sequence" index, K is the "vector" |
| pairIndexList | Input: list of SFT pairs |
| Gamma_ave | Input: vector of aa+bb values |
| Gamma_circ | Input: vector of ab-ba values |
| coordSys | Input: coordinate directions for metric |
Definition at line 1543 of file PulsarCrossCorr_v2.c.
| 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
calculate phase metric for CW cross-correlation search, as well as vector used for parameter offsets This calculates the metric defined in (4.7) of Whelan et al 2015 and the parameter offset epsilon_i (not including the cosi-dependent prefactor) in defined in (4.8)
| [out] | g_ij | parameter space metric |
| [out] | eps_i | parameter offset vector from (4.8) of WSZP15 |
| [out] | sumGammaSq | sum of (Gamma_ave)^2 for normalization and sensitivity |
| [in] | phaseDerivs | dPhi_K/dlambda_i; i is the "sequence" index, K is the "vector" |
| [in] | resampMultiPairs | resamp multi list of SFT pairs |
| [in] | Gamma_ave | : vector of aa+bb values |
| [in] | Gamma_circ | vector of ab-ba values |
| [in] | coordSys | coordinate directions for metric |
Definition at line 1634 of file PulsarCrossCorr_v2.c.
| 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 | ||
| ) |
| hSens | Output: sensitivity |
| g_ff | Output: Diagonal frequency metric element |
| g_aa | Output: Diagonal binary projected semimajor axis metric element |
| g_TT | Output: Diagonal reference time metric element |
| g_pp | Output: Diagonal orbital period metric element |
| weightedMuTAve | output: weighred T mean |
| DopplerParams | Input: pulsar/binary orbit paramaters |
| G_alpha | Input: vector of curlyGunshifted values |
| pairIndexList | Input: list of SFT pairs |
| indexList | Input: list of SFTs |
| sfts | Input: set of per-detector SFT vectors |
| multiWeights | Input: nomalizeation factor S^-1 & weights for each SFT |
Definition at line 1735 of file PulsarCrossCorr_v2.c.
| int XLALCalculateLMXBCrossCorrDiagMetricShort | ( | REAL8 * | hSens, |
| REAL8 * | g_ff, | ||
| REAL8 * | g_aa, | ||
| REAL8 * | g_TT, | ||
| REAL8 * | g_pp, | ||
| const PulsarDopplerParams | DopplerParams, | ||
| const REAL8Vector *restrict | G_alpha, | ||
| const MultiResampSFTPairMultiIndexList *restrict | resampMultiPairIndexList, | ||
| const MultiLIGOTimeGPSVector *restrict | timestamps, | ||
| const MultiNoiseWeights *restrict | multiWeights | ||
| ) |
MODIFIED for Tshort: calculate metric diagonal components, also include the estimation of sensitivity E[rho]/(h_0)^2.
| [out] | hSens | sensitivity |
| [out] | g_ff | Diagonal frequency metric element |
| [out] | g_aa | Diagonal binary projected semimajor axis metric element |
| [out] | g_TT | Diagonal reference time metric element |
| [out] | g_pp | Diagonal orbital period metric element |
| [in] | DopplerParams | pulsar/binary orbit paramaters |
| [in] | G_alpha | vector of curlyGunshifted values |
| [in] | resampMultiPairIndexList | resamp list of SFT pairs |
| [in] | timestamps | timestamps for resampling |
| [in] | multiWeights | normalization factor S^-1 & weights for each SFT |
Definition at line 1812 of file PulsarCrossCorr_v2.c.
| LIGOTimeGPSVector * XLALExtractTimestampsFromSFTsShort | ( | 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 */
Function to extract timestamps with tShort
| [out] | sciFlag | science or not flag |
| [in] | sfts | input SFT-vector |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
Definition at line 2024 of file PulsarCrossCorr_v2.c.
| LIGOTimeGPSVector * XLALGenerateTshortTimestamps | ( | const REAL8 | tShort, |
| const UINT4 | numShortPerDet, | ||
| const LIGOTimeGPS | epoch | ||
| ) |
Generate timestamps for one detector with tShort.
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
| [in] | epoch | start time for first Tshort |
Definition at line 2070 of file PulsarCrossCorr_v2.c.
| LIGOTimeGPSVector * XLALModifyTimestampsFromSFTsShort | ( | REAL8TimeSeries ** | sciFlag, |
| const LIGOTimeGPSVector *restrict | Times, | ||
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet | ||
| ) |
Modify timestamps from one detector with tShort.
| [out] | sciFlag | science or not flag |
| [in] | Times | input SFT-vector |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
Definition at line 2105 of file PulsarCrossCorr_v2.c.
| MultiLIGOTimeGPSVector * XLALGenerateMultiTshortTimestamps | ( | const MultiLIGOTimeGPSVector *restrict | multiTimes, |
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet, | ||
| const BOOLEAN | alignTShorts | ||
| ) |
Given a multi-SFT vector, return a MultiLIGOTimeGPSVector holding the modified SFT timestamps (production function using tShort)
| [in] | multiTimes | multiTimes, standard |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
| [in] | alignTShorts | align segments between detectors |
Definition at line 2154 of file PulsarCrossCorr_v2.c.
| MultiLIGOTimeGPSVector * XLALModifyMultiTimestampsFromSFTs | ( | MultiREAL8TimeSeries ** | scienceFlagVect, |
| const MultiLIGOTimeGPSVector *restrict | multiTimes, | ||
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet | ||
| ) |
Given a multi-SFT vector, return a MultiLIGOTimeGPSVector holding the modified SFT timestamps (production function using tShort)
| [out] | scienceFlagVect | science flag vector |
| [in] | multiTimes | multiTimes, standard |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
Definition at line 2218 of file PulsarCrossCorr_v2.c.
| MultiLIGOTimeGPSVector * XLALExtractMultiTimestampsFromSFTsShort | ( | 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)
| [out] | scienceFlagVect | science flag vector |
| [in] | multiSFTs | multiSFT vector |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
Definition at line 2282 of file PulsarCrossCorr_v2.c.
| int XLALFillDetectorTensorShort | ( | DetectorState * | detState, |
| const LALDetector * | detector | ||
| ) |
(test function) fill detector state with tShort, importing various slightly-modified LALPulsar functions for testing
| [in,out] | detState | : detector state: fill in detector-tensor |
| [in] | detector | : which detector |
Definition at line 2345 of file PulsarCrossCorr_v2.c.
| DetectorStateSeries * XLALGetDetectorStatesShort | ( | const LIGOTimeGPSVector * | timestamps, |
| const LALDetector * | detector, | ||
| const EphemerisData * | edat, | ||
| REAL8 | tOffset, | ||
| REAL8 | tShort, | ||
| UINT4 | numShortPerDet | ||
| ) |
(test function) get detector states for tShort
| timestamps | array of GPS timestamps t_i | |
| detector | detector info | |
| edat | ephemeris file data | |
| tOffset | compute detector states at timestamps SHIFTED by tOffset | |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort per detector |
Definition at line 2424 of file PulsarCrossCorr_v2.c.
| MultiDetectorStateSeries * XLALGetMultiDetectorStatesShort | ( | const MultiLIGOTimeGPSVector * | multiTS, |
| const MultiLALDetector * | multiIFO, | ||
| const EphemerisData * | edat, | ||
| REAL8 | tOffset, | ||
| REAL8 | tShort, | ||
| UINT4 | numShortPerDet | ||
| ) |
(test function) get multi detector states for tShort
| [in] | multiTS | multi-IFO timestamps |
| [in] | multiIFO | multi-IFO array holding detector info |
| [in] | edat | ephemeris data |
| [in] | tOffset | shift all timestamps by this amount |
| [in] | tShort | new time baseline, Tshort |
| [in] | numShortPerDet | number of Tshort segments per detector |
Definition at line 2531 of file PulsarCrossCorr_v2.c.
| int XLALModifyAMCoeffsWeights | ( | REAL8Vector ** | resampMultiWeightsX, |
| const MultiNoiseWeights *restrict | multiWeights, | ||
| const REAL8 | tShort, | ||
| const REAL8 | tSFTOld, | ||
| const UINT4 | numShortPerDet, | ||
| const MultiLIGOTimeGPSVector *restrict | multiTimes, | ||
| const UINT4 | maxNumStepsOldIfGapless, | ||
| const UINT4 | X | ||
| ) |
Modify amplitude weight coefficients for tShort.
| [out] | resampMultiWeightsX | new weights |
| [in] | multiWeights | old weights |
| [in] | tShort | new time baseline, Tshort |
| [in] | tSFTOld | old time baseline, tSFTOld |
| [in] | numShortPerDet | number of tShort segments per detector |
| [in] | multiTimes | multi-times vector to tell us when the SFTs were |
| [in] | maxNumStepsOldIfGapless | how many SFTs there would be without gaps |
| [in] | X | detector number |
Definition at line 2612 of file PulsarCrossCorr_v2.c.
| MultiNoiseWeights * XLALModifyMultiWeights | ( | const MultiNoiseWeights *restrict | multiWeights, |
| const REAL8 | tShort, | ||
| const REAL8 | tSFTOld, | ||
| const UINT4 | numShortPerDet, | ||
| const MultiLIGOTimeGPSVector *restrict | multiTimes | ||
| ) |
Modify multiple detectors' amplitude weight coefficients for tShort.
| [in] | multiWeights | old weights |
| [in] | tShort | new time baseline, Tshort |
| [in] | tSFTOld | old time baseline, tSFTOld |
| [in] | numShortPerDet | number of tShort segments per detector |
| [in] | multiTimes | multi-times vector to tell us when the SFTs were |
Definition at line 2699 of file PulsarCrossCorr_v2.c.
| int XLALModifyMultiAMCoeffsWeights | ( | MultiNoiseWeights ** | multiWeights, |
| const REAL8 | tShort, | ||
| const REAL8 | tSFTOld, | ||
| const UINT4 | numShortPerDet, | ||
| const MultiLIGOTimeGPSVector *restrict | multiTimes | ||
| ) |
Modify multiple detectors' amplitude weight coefficients for tShort.
| multiWeights | [in/out] old and new weights | |
| [in] | tShort | new time baseline, Tshort |
| [in] | tSFTOld | old time baseline, tSFTOld |
| [in] | numShortPerDet | number of tShort segments per detector |
| [in] | multiTimes | multi-times vector to tell us when the SFTs were |
Definition at line 2754 of file PulsarCrossCorr_v2.c.
| int XLALWeightMultiAMCoeffsShort | ( | MultiAMCoeffs * | multiAMcoef, |
| const MultiNoiseWeights * | multiWeights, | ||
| REAL8 | tShort, | ||
| REAL8 | tSFTOld, | ||
| UINT4 | numShortPerDet, | ||
| MultiLIGOTimeGPSVector * | multiTimes | ||
| ) |
(test function) used for weighting multi amplitude modulation coefficients
| multiAMcoef | [in/out] amplitude coefficients | |
| multiWeights | [in/out] weights | |
| [in] | tShort | new time baseline, Tshort |
| [in] | tSFTOld | old time baseline, tSFTOld |
| [in] | numShortPerDet | number of tShort segments per detector |
| [in] | multiTimes | multi-times vector to tell us when the SFTs were |
Definition at line 2801 of file PulsarCrossCorr_v2.c.
| AMCoeffs * XLALComputeAMCoeffsShort | ( | const DetectorStateSeries * | DetectorStates, |
| SkyPosition | skypos | ||
| ) |
(test function) used for computing amplitude modulation weights
| DetectorStates | timeseries of detector states |
| skypos | {alpha,delta} of the source |
Definition at line 2972 of file PulsarCrossCorr_v2.c.
| MultiAMCoeffs * XLALComputeMultiAMCoeffsShort | ( | 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
| [in] | multiDetStates | detector-states at timestamps t_i |
| [in] | multiWeights | noise-weights at timestamps t_i (can be NULL) |
| skypos | source sky-position [in equatorial coords!] | |
| [in] | tShort | new time baseline, Tshort |
| [in] | tSFTOld | old time baseline, tSFTOld |
| [in] | numShortPerDet | number of tShort segments per detector |
| [in] | multiTimes | multi-times vector to tell us when the SFTs were |
Definition at line 3046 of file PulsarCrossCorr_v2.c.
| UINT4 XLALCrossCorrNumShortPerDetector | ( | const REAL8 | resampTshort, |
| const INT4 | startTime, | ||
| const INT4 | endTime | ||
| ) |
Compute the number of tShort segments per detector.
| [in] | resampTshort | resampling Tshort |
| [in] | startTime | start time of observation |
| [in] | endTime | end time of observation |
Definition at line 3108 of file PulsarCrossCorr_v2.c.
| REAL8TimeSeries * XLALCrossCorrGapFinderResamp | ( | LIGOTimeGPSVector *restrict | timestamps, |
| const LIGOTimeGPSVector *restrict | Times | ||
| ) |
Find gaps in the data given the SFTs.
| [in] | timestamps | timestamps vector |
| [in] | Times | input SFT-vector |
Definition at line 3123 of file PulsarCrossCorr_v2.c.
| REAL8TimeSeries * XLALCrossCorrGapFinderResampAlt | ( | LIGOTimeGPSVector *restrict | timestamps, |
| const SFTVector *restrict | sfts | ||
| ) |
(test function) find gaps in the data given the SFTs
| [in] | timestamps | timestamps vector |
| [in] | sfts | input SFT-vector |
Definition at line 3179 of file PulsarCrossCorr_v2.c.
| int XLALEquipCrossCorrPairsWithScienceFlags | ( | MultiResampSFTPairMultiIndexList * | resampMultiPairs, |
| MultiREAL8TimeSeries * | scienceFlagVect | ||
| ) |
Demarcate pairs with flags about whether data exists in zero-padded timeseries.
| [in] | resampMultiPairs | resampling pairs |
| [in] | scienceFlagVect | science flags |
Definition at line 3234 of file PulsarCrossCorr_v2.c.
| 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 */
Generates a resampling workspace for CrossCorr
| [out] | wsOut | workspace for one cross-correlation |
| [out] | ws1KFaX_kOut | holder for detector 1 Fa |
| [out] | ws1KFbX_kOut | holder for detector 1 Fb |
| [out] | ws2LFaX_kOut | holder for detector 2 Fa |
| [out] | ws2LFbX_kOut | holder for detector 2 Fb |
| [out] | multiTimeSeries_SRC_aOut | resampling A time series |
| [out] | multiTimeSeries_SRC_bOut | resampling B time series |
| [in] | binaryTemplateSpacings | binary template spacings |
| [in] | resampFstatInput | resampling f-statistic input |
| [in] | numFreqBins | number of frequency bins |
| [in] | tCoh | Tcoh = 2 * max lag + resampling tShort |
| [in] | treatWarningsAsErrors | abort program if any warnings encountered |
Definition at line 3311 of file PulsarCrossCorr_v2.c.
| void XLALDestroySFTIndexList | ( | SFTIndexList * | sftIndices | ) |
Destroy a SFTIndexList structure.
Note, this is "NULL-robust" in the sense that it will not crash on NULL-entries anywhere in this struct, so it can be used for failure-cleanup even on incomplete structs
Definition at line 3446 of file PulsarCrossCorr_v2.c.
| void XLALDestroySFTPairIndexList | ( | SFTPairIndexList * | sftPairs | ) |
Destroy a SFTPairIndexList structure.
Note, this is "NULL-robust" in the sense that it will not crash on NULL-entries anywhere in this struct, so it can be used for failure-cleanup even on incomplete structs
Definition at line 3469 of file PulsarCrossCorr_v2.c.
| void XLALDestroyResampSFTIndexList | ( | ResampSFTIndexList * | sftResampList | ) |
Definition at line 3485 of file PulsarCrossCorr_v2.c.
| void XLALDestroyResampSFTMultiIndexList | ( | ResampSFTMultiIndexList * | sftResampMultiList | ) |
Definition at line 3500 of file PulsarCrossCorr_v2.c.
| void XLALDestroyResampSFTPairMultiIndexList | ( | ResampSFTPairMultiIndexList * | sftResampPairMultiList | ) |
Definition at line 3518 of file PulsarCrossCorr_v2.c.
| void XLALDestroyMultiResampSFTPairMultiIndexList | ( | MultiResampSFTPairMultiIndexList * | sftMultiPairsResamp | ) |
Definition at line 3537 of file PulsarCrossCorr_v2.c.
| void XLALDestroyMultiMatchList | ( | MultiResampSFTMultiCountList * | localMultiListOfLmatchingGivenMultiK | ) |
Definition at line 3563 of file PulsarCrossCorr_v2.c.
| void XLALDestroyResampCrossCorrWorkspace | ( | void * | workspace | ) |
Definition at line 3592 of file PulsarCrossCorr_v2.c.
|
static |
Imported and modified from ComputeFstat_Resamp.c.
| [out] | xOut | the spindown-corrected SRC-frame timeseries |
| [in] | xIn | the input SRC-frame timeseries |
| [in] | doppler | containing spindown parameters |
| [in] | freqShift | frequency-shift to apply, sign is "new - old" |
| [in] | indexStartResamp | index in resampling time series to start |
| [in] | indexEndResamp | index in resampling time series to end |
| [in] | numSamplesIn | number of samples in the output |
| [in] | insertPoint | number of sample indices offset to begin output |
Definition at line 3621 of file PulsarCrossCorr_v2.c.
|
static |
Compute the equivalent of Fa and Fb for CrossCorr's rho statistic.
| [out] | ws | contains modified Fa and Fb for cross-correlation |
| [out] | wsFaX_k | contains modified and normalized Fa for cross-correlation statistic |
| [out] | wsFbX_k | contains modified and normalized Fb for cross-correlation statistic |
| [in] | resampMultiPairs | resamp multi list of SFT pairs |
| [in] | multiTimeSeries_SRC_a | Resampled, heterodyned A(t)x(t) |
| [in] | multiTimeSeries_SRC_b | Resampled, heterodyned B(t)x(t) |
| [in] | dopplerpos | Doppler point to search |
| [in] | binaryTemplateSpacings | spacings for the search |
| [in] | SRCsampPerTcoh | floating point number of samples per tShort equivalent |
| [in] | detX | detector X index |
| [in] | sftK | tShort K index |
| [in] | detY | detector Y index |
| [in] | isL | Indicates that is not K but L, the second half of the cross-corr pair |
Definition at line 3697 of file PulsarCrossCorr_v2.c.
| int XLALCreateSFTPairIndexListAccurateResamp | ( | SFTPairIndexList ** | pairIndexList, |
| UINT4VectorSequence ** | sftPairForResampPair, | ||
| SFTIndexList * | indexList, | ||
| MultiResampSFTPairMultiIndexList * | resampPairs, | ||
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes, | ||
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | resampMultiTimes | ||
| ) |
| int XLALCreateSFTPairIndexListShortResamp | ( | MultiResampSFTPairMultiIndexList ** | resampPairIndexList, |
| const REAL8 | maxLag, | ||
| const BOOLEAN | inclAutoCorr, | ||
| const BOOLEAN | inclSameDetector, | ||
| const REAL8 | Tsft, | ||
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes | ||
| ) |
| 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 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 | ||
| ) |
| LIGOTimeGPSVector * XLALModifyTimestampsFromSFTsShort | ( | REAL8TimeSeries ** | sciFlag, |
| const LIGOTimeGPSVector *_LAL_RESTRICT_ | Times, | ||
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet | ||
| ) |
| MultiLIGOTimeGPSVector * XLALGenerateMultiTshortTimestamps | ( | const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes, |
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet, | ||
| const BOOLEAN | alignTShorts | ||
| ) |
| MultiLIGOTimeGPSVector * XLALModifyMultiTimestampsFromSFTs | ( | MultiREAL8TimeSeries ** | scienceFlagVect, |
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes, | ||
| const REAL8 | tShort, | ||
| const UINT4 | numShortPerDet | ||
| ) |
| 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 | ||
| ) |
| MultiNoiseWeights * XLALModifyMultiWeights | ( | const MultiNoiseWeights *_LAL_RESTRICT_ | multiWeights, |
| const REAL8 | tShort, | ||
| const REAL8 | tSFTOld, | ||
| const UINT4 | numShortPerDet, | ||
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes | ||
| ) |
| int XLALModifyMultiAMCoeffsWeights | ( | MultiNoiseWeights ** | multiWeights, |
| const REAL8 | tShort, | ||
| const REAL8 | tSFTOld, | ||
| const UINT4 | numShortPerDet, | ||
| const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ | multiTimes | ||
| ) |
| REAL8TimeSeries * XLALCrossCorrGapFinderResamp | ( | LIGOTimeGPSVector *_LAL_RESTRICT_ | timestamps, |
| const LIGOTimeGPSVector *_LAL_RESTRICT_ | Times | ||
| ) |
| REAL8TimeSeries * XLALCrossCorrGapFinderResampAlt | ( | LIGOTimeGPSVector *_LAL_RESTRICT_ | timestamps, |
| const SFTVector *_LAL_RESTRICT_ | sfts | ||
| ) |