LALPulsar  6.1.0.1-fe68b98
PulsarCrossCorr_v2.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2012, 2013 John Whelan, Shane Larson and Badri Krishnan
3  * Copyright (C) 2013, 2014 Badri Krishnan, John Whelan, Yuanhao Zhang
4  * Copyright (C) 2016, 2017 Grant David Meadors
5  * Copyright (C) 2023 John Whelan
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with with program; see the file COPYING. If not, write to the
19  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
20  * MA 02110-1301 USA
21  */
22 #ifndef _PULSARCROSSCORRV2_H
23 #define _PULSARCROSSCORRV2_H
24 
25 #ifdef __cplusplus
26 extern "C" {
27 #endif
28 
29 /**
30  * \defgroup PulsarCrossCorr_v2_h Header PulsarCrossCorr_v2.h
31  * \ingroup lalpulsar_crosscorr
32  * \author John Whelan, Yuanhao Zhang, Shane Larson, Badri Krishnan
33  * \date 2012, 2013, 2014
34  * \brief Header-file for XLAL routines for v2 CW cross-correlation searches
35  *
36  */
37 /** @{ */
38 
39 #include <sys/types.h>
40 #include <sys/stat.h>
41 #include <fcntl.h>
42 #include <string.h>
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <math.h>
46 #if HAVE_GLOB_H
47 #include <glob.h>
48 #endif
49 #include <time.h>
50 #include <errno.h>
51 #include <lal/AVFactories.h>
52 #include <lal/Date.h>
53 #include <lal/DetectorSite.h>
54 #include <lal/LALDatatypes.h>
55 #include <lal/LALHough.h>
56 #include <lal/RngMedBias.h>
57 #include <lal/LALRunningMedian.h>
58 #include <lal/Velocity.h>
59 #include <lal/Statistics.h>
60 #include <lal/ComputeFstat.h>
61 #include <lal/LALConstants.h>
62 #include <lal/SFTfileIO.h>
63 #include <lal/NormalizeSFTRngMed.h>
64 #include <lal/LALInitBarycenter.h>
65 #include <lal/SFTClean.h>
66 #include <gsl/gsl_cdf.h>
67 #include <gsl/gsl_sf_trig.h>
68 #include <lal/FrequencySeries.h>
69 #include <lal/Sequence.h>
70 #include <lal/SinCosLUT.h>
71 #include <lal/LogPrintf.h>
72 #include <lal/UniversalDopplerMetric.h>
73 #include <lal/ExtrapolatePulsarSpins.h>
74 #include <lal/TimeSeries.h>
75 #include <lal/Units.h>
76 #include <lal/FFTWMutex.h>
77 #include <fftw3.h>
78 #include "ComputeFstat.h"
79 #include "GeneratePulsarSignal.h"
80 
81 /* ******************************************************************
82  * Structure, enum, union, etc., typdefs.
83  */
84 
85 /** Index to refer to an SFT given a set of SFTs from several different detectors */
86 typedef struct tagSFTIndex {
87  UINT4 detInd; /**< index of detector in list */
88  UINT4 sftInd; /**< index of SFT in list for this detector */
89 } SFTIndex;
90 
91 /** List of SFT indices */
92 typedef struct tagSFTIndexList {
93  UINT4 length; /**< number of SFTs */
94  SFTIndex *data; /**< array of SFT indices */
95 } SFTIndexList;
96 
97 /** Index to refer to a pair of SFTs */
98 typedef struct tagSFTPairIndex {
99 #if 0
100  SFTIndex sftInd1; /**< index of 1st SFT in pair */
101  SFTIndex sftInd2; /**< index of 2nd SFT in pair */
102 #endif
103  UINT4 sftNum[2]; /**< ordinal numbers of first and second SFTs */
104 } SFTPairIndex;
105 
106 /** List of SFT pair indices */
107 typedef struct tagSFTPairIndexList {
108  UINT4 length; /**< number of SFT Pairs */
109  SFTPairIndex *data; /**< array of SFT Pair indices */
111 
112 
113 /** Resampling Counter of matching SFTs for a given detector Y_K_X matching SFT K_X */
114 typedef struct tagSFTCount {
115  UINT4 detInd; /**< original vector index of detector Y */
116  UINT4 sftCount; /**< number of matching SFTs */
117 } SFTCount;
118 
119 /* Special multi-types for resampling */
120 /* These are MULTI-structs */
121 
122 
123 /** INNER List of SFT indices */
124 typedef struct tagResampSFTMultiCountList {
125  UINT4 length; /**< number of detectors Y */
126  UINT4 sftInd; /**< original vector index of sft K */
127  SFTCount *data; /**< arrays of count of SFTs L_Y, at given detector Y, that match each given SFT K_X */
129 
130 /** MIDDLE List of SFT indices */
131 typedef struct tagResampSFTMultiCountListDet {
132  UINT4 length; /**< number of SFTs K_X */
133  UINT4 detInd; /**< original vector index of detector X */
134  ResampSFTMultiCountList *data; /**< array of SFT L_Y indices for given SFT K_X at detector X */
136 
137 /** OUTER List of SFT indices */
138 typedef struct tagMultiResampSFTMultiCountList {
139  UINT4 length; /**< number of detectors X */
140  ResampSFTMultiCountListDet *data; /**< list per-detectors X */
142 
143 
144 /** Resampling: Index to refer to fields of an SFT given a specific index L_Y_K_X */
145 typedef struct tagResampSFTIndex {
146  UINT4 detInd; /**< index of detector in list */
147  UINT4 sftInd; /**< index of SFT in list for this detector L */
148  UINT4 flatInd; /**< index in the flat SFT list */
149  UINT4 pairInd; /**< index of SFT among overall pairs */
150  REAL4 sciFlag; /**< science-or-not value: 1 is science mode, 0 is in a gap */
152 
153 /** 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 */
154 typedef struct tagResampSFTIndexList {
155  UINT4 length; /**< number of SFTs */
156  UINT4 detInd; /**< original vector index of detector Y */
157  ResampSFTIndex *data; /**< array of SFT indices */
159 
160 /** Resampling: Multi List of indices of SFT L_Y, for a given sft K_X */
161 typedef struct tagResampSFTMultiIndexList {
162  UINT4 length; /**< number of detectors Y */
163  UINT4 sftInd; /**< original vector index of sft K */
164  UINT4 flatInd; /**< index in the flat SFT list */
165  REAL4 sciFlag; /**< science-or-not value: 1 is science mode, 0 is in a gap */
166  ResampSFTIndexList *data; /**< array, per detector, of lists of SFT indices */
168 
169 /** Resampling Multi List of SFT pair indices (L_Y_K), for a given detector X */
170 typedef struct tagResampSFTPairMultiIndexList {
171  UINT4 length; /**< number of K SFTs at detector X */
172  UINT4 detInd; /**< original vector index of detector X */
173  ResampSFTMultiIndexList *data; /**< array of SFT L_Y indices for given SFT K_X at detector X */
175 
176 
177 /** Resampling Multi List of all paired SFTs L_Y_K_X, top level (multiple pairs, multiple detectors) */
178 typedef struct tagMultiResampSFTPairIndexMultiList {
179  UINT4 allPairCount; /**< count of all pairs */
180  UINT4 oldPairCount; /**< count of sft pairs, old-style */
181  UINT4 sftTotalCount; /**< count of all sfts */
182  REAL8 maxLag; /**< Maximum allowed lag time */
183  BOOLEAN inclAutoCorr; /**< Do we include auto-correlation? */
184  BOOLEAN inclSameDetector; /**< Do we include same detectors? */
185  REAL8 Tsft; /**< Duration of one SFT */
186  REAL8 Tshort; /**< Duration of resamp Tshort */
187  SFTIndexList *indexList; /**< Make an overall flat index list */
188  SFTPairIndexList *pairIndexList; /**< Make a classic pair index list */
189  UINT4 length; /**< number of detectors X */
190  ResampSFTPairMultiIndexList *data; /**< list per-detector X */
192 
193 // ----- local types ----------
194 
195 typedef struct tagCrossCorrTimings_t {
196  REAL8 Total; //!< total time spent in XLALComputeFstatResamp()
197  REAL8 Bary; //!< time spent in barycentric resampling
198  REAL8 Spin; //!< time spent in spindown+frequency correction
199  REAL8 FFT; //!< time spent in FFT
200  REAL8 Norm; //!< time spent normalizing the final Fa,Fb
201  REAL8 Fab2F; //!< time to compute Fstat from {Fa,Fb}
202  REAL8 Mem; //!< time to realloc and memset-0 arrays
203  REAL8 SumFabX; //!< time to sum_X Fab^X
204  REAL8 F1Buf; //!< Resampling timing 'constant': Fstat time per template per detector for a 'buffered' case (same skypos, same numFreqBins)
205  REAL8 F1NoBuf; //!< Resampling timing 'constant': Fstat time per template per detector for an 'unbuffered' usage (different skypos and numFreqBins)
207 
208 typedef struct tagResampCrossCorrTimingInfo {
209  // NOTE: all times refer to a single-detector timing case
210  BOOLEAN collectTiming; //!< turn on/off the collection of F-stat-method-specific timing-data (stored in workspace)
211 
212  UINT4 numFreqBins; //!< number of frequency bins to compute F-stat for
213  UINT4 numSFTs; //!< total number of SFTs used
214  UINT4 numDetectors; //!< number of detectors
215  UINT4 numSamplesFFT0; //!< 'original' length of barycentered timeseries to be FFT'ed
216  UINT4 numSamplesFFT; //!< actual length of FFT, potentially rounded up to power-of-2 for efficiency
219 
220 // ----- workspace ----------
221 
222 typedef struct tagResampCrossCorrWorkspace {
223  // intermediate quantities to interpolate and operate on SRC-frame timeseries
224  COMPLEX8Vector *TStmp1_SRC; //!< can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
225  COMPLEX8Vector *TStmp2_SRC; //!< can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
226  REAL8Vector *SRCtimes_DET; //!< holds uniformly-spaced SRC-frame timesteps translated into detector frame [for interpolation]
227 
228  // input padded timeseries ts(t) and output Fab(f) of length 'numSamplesFFT' and corresponding fftw plan
229  UINT4 numSamplesFFT; //!< allocated number of zero-padded SRC-frame time samples (related to dFreq)
230  UINT4 decimateFFT; //!< output every n-th frequency bin, with n>1 iff (dFreq > 1/Tspan), and was internally decreased by n
231  fftwf_plan fftplan; //!< buffer FFT plan for given numSamplesOut length
232  COMPLEX8 *TS_FFT; //!< zero-padded, spindown-corr SRC-frame TS
233  COMPLEX8 *FabX_Raw; //!< raw full-band FFT result Fa,Fb
234 
235  // arrays of size numFreqBinsOut over frequency bins f_k:
236  UINT4 numFreqBinsOut; //!< number of output frequency bins {f_k}
237  COMPLEX8 *FaX_k; //!< properly normalized F_a^X(f_k) over output bins
238  COMPLEX8 *FbX_k; //!< properly normalized F_b^X(f_k) over output bins
239  COMPLEX8 *Fa_k; //!< properly normalized F_a(f_k) over output bins
240  COMPLEX8 *Fb_k; //!< properly normalized F_b(f_k) over output bins
241  UINT4 numFreqBinsAlloc; //!< internal: keep track of allocated length of frequency-arrays
242 
243  ResampCrossCorrTimingInfo *timingInfo; //!< pointer to storage for collecting timing data (which lives in ResampMethodData)
245 
246 /* end Resampling multi-types */
247 
248 /** A collection of UINT4Vectors -- one for each IFO */
249 /* Probably belongs in SFTUtils.h */
250 typedef struct tagMultiUINT4Vector {
251 #ifdef SWIG /* SWIG interface directives */
252  SWIGLAL( ARRAY_1D( MultiUINT4Vector, UINT4Vector *, data, UINT4, length ) );
253 #endif /* SWIG */
254  UINT4 length; /**< number of ifos */
255  UINT4Vector **data; /**< unit4vector for each ifo */
257 
258 /*
259 * Functions Declarations (i.e., prototypes).
260 */
261 
263 (
264  REAL8Vector *shiftedFreqs,
265  UINT4Vector *lowestBins,
266  COMPLEX8Vector *expSignalPhases,
267  REAL8VectorSequence *sincList,
268  UINT4 numBins,
269  PulsarDopplerParams *dopp,
270  SFTIndexList *sfts,
271  MultiSFTVector *inputSFTs,
272  MultiSSBtimes *multiTimes,
273  MultiUINT4Vector *badBins,
274  REAL8 Tsft
275 )
276 ;
277 
279 (
280  SFTIndexList **indexList,
281  MultiSFTVector *sfts
282 )
283 ;
284 
286 (
287  SFTPairIndexList **pairIndexList,
288  SFTIndexList *indexList,
289  MultiSFTVector *sfts,
290  REAL8 maxLag,
291  BOOLEAN inclAutoCorr
292 )
293 ;
294 
296 (
297  MultiResampSFTPairMultiIndexList **resampPairIndexList,
298  SFTPairIndexList **pairIndexList,
299  SFTIndexList *indexList,
300  MultiSFTVector *sfts,
301  REAL8 maxLag,
302  BOOLEAN inclAutoCorr,
303  BOOLEAN inclSameDetector,
304  REAL8 Tsft,
305  REAL8 Tshort
306 )
307 ;
308 
310 (
311  MultiResampSFTPairMultiIndexList **resampPairIndexList,
312  const REAL8 maxLag,
313  const BOOLEAN inclAutoCorr,
314  const BOOLEAN inclSameDetector,
315  const REAL8 Tsft,
316  const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes
317 )
318 ;
319 
321 (
322  MultiResampSFTPairMultiIndexList *resampMultiPairIndexList
323 );
324 
326 (
327  REAL8Vector **Gamma_ave,
328  REAL8Vector **Gamma_circ,
329  SFTPairIndexList *pairIndexList,
330  SFTIndexList *indexList,
331  MultiAMCoeffs *multiCoeffs
332 )
333 ;
334 
336 (
337  REAL8Vector **Gamma_ave,
338  REAL8Vector **Gamma_circ,
339  MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
340  MultiAMCoeffs *multiCoeffs
341 )
342 ;
343 
345 (
346  REAL8Vector **Gamma_ave,
347  REAL8Vector **Gamma_circ,
348  MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
349  MultiAMCoeffs *multiCoeffs
350 )
351 ;
352 
354 (
355  REAL8 *ccStat,
356  REAL8 *evSquared,
357  REAL8Vector *curlyGAmp,
358  COMPLEX8Vector *expSignalPhases,
359  UINT4Vector *lowestBins,
360  REAL8VectorSequence *sincList,
361  SFTPairIndexList *sftPairs,
362  SFTIndexList *sftIndices,
363  MultiSFTVector *inputSFTs,
364  MultiNoiseWeights *multiWeights,
365  UINT4 numBins
366 )
367 ;
368 
370 (
371  REAL8Vector *_LAL_RESTRICT_ ccStatVector,
372  REAL8Vector *_LAL_RESTRICT_ evSquaredVector,
373  REAL8Vector *_LAL_RESTRICT_ numeEquivAve,
374  REAL8Vector *_LAL_RESTRICT_ numeEquivCirc,
375  const REAL8Vector *_LAL_RESTRICT_ resampCurlyGAmp,
376  const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampPairs,
377  const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights,
378  const PulsarDopplerParams *_LAL_RESTRICT_ binaryTemplateSpacings,
379  const PulsarDopplerParams *_LAL_RESTRICT_ dopplerpos,
380  const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_a,
381  const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_b,
382  ResampCrossCorrWorkspace *_LAL_RESTRICT_ ws,
383  COMPLEX8 *_LAL_RESTRICT_ ws1KFaX_k,
384  COMPLEX8 *_LAL_RESTRICT_ ws1KFbX_k,
385  COMPLEX8 *_LAL_RESTRICT_ ws2LFaX_k,
386  COMPLEX8 *_LAL_RESTRICT_ ws2LFbX_k
387 )
388 ;
389 
391 (
392  REAL8VectorSequence **phaseDerivs,
393  const PulsarDopplerParams *dopplerPoint,
394  const EphemerisData *edat,
395  SFTIndexList *indexList,
396  MultiSSBtimes *multiTimes,
397  const DopplerCoordinateSystem *coordSys
398 )
399 ;
400 
402 (
403  REAL8VectorSequence **resampPhaseDerivs,
404  const PulsarDopplerParams *dopplerPoint,
405  const EphemerisData *edat,
406  SFTIndexList *indexList,
407  MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
408  MultiSSBtimes *multiTimes,
409  const DopplerCoordinateSystem *coordSys
410 )
411 ;
412 
414 (
415  gsl_matrix **g_ij,
416  gsl_vector **eps_i,
417  REAL8 *sumGammaSq,
418  const REAL8VectorSequence *phaseDerivs,
419  const SFTPairIndexList *pairIndexList,
420  const REAL8Vector *Gamma_ave,
421  const REAL8Vector *Gamma_circ,
422  const DopplerCoordinateSystem *coordSys
423 );
424 
426 (
427  gsl_matrix **g_ij,
428  gsl_vector **eps_i,
429  REAL8 *sumGammaSq,
430  const REAL8VectorSequence *phaseDerivs,
431  const MultiResampSFTPairMultiIndexList *resampMultiPairs,
432  const REAL8Vector *Gamma_ave,
433  const REAL8Vector *Gamma_circ,
434  const DopplerCoordinateSystem *coordSys
435 );
436 
438 (
439  REAL8 *hSens,
440  REAL8 *g_ff,
441  REAL8 *g_aa,
442  REAL8 *g_TT,
443  REAL8 *g_pp,
444  REAL8 *weightedMuTAve,
445  PulsarDopplerParams DopplerParams,
446  REAL8Vector *G_alpha,
447  SFTPairIndexList *pairIndexList,
448  SFTIndexList *indexList,
449  MultiSFTVector *sfts,
450  MultiNoiseWeights *multiWeights
451 )
452 ;
453 
455 (
456  REAL8 *hSens,
457  REAL8 *g_ff,
458  REAL8 *g_aa,
459  REAL8 *g_TT,
460  REAL8 *g_pp,
461  const PulsarDopplerParams DopplerParams,
462  const REAL8Vector *_LAL_RESTRICT_ G_alpha,
463  const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampMultiPairIndexList,
464  const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
465  const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights
466 )
467 ;
468 
469 /* (possible future function) */
470 //int
471 //XLALBesselCrossCorrOrbitalSpaceStep(
472 // COMPLEX8 * xTildeOut,
473 // COMPLEX8 * xTildeIn,
474 // PulsarDopplerParams * dopplerpos,
475 // PulsarDopplerParams * binaryTemplateSpacings,
476 // INT4 aStep,
477 // INT4 pStep,
478 // INT4 tStep,
479 // INT4 nFFT,
480 // UINT4 nTermsMax
481 //);
482 
485  REAL8TimeSeries **sciFlag,
486  const LIGOTimeGPSVector *_LAL_RESTRICT_ Times,
487  const REAL8 tShort,
488  const UINT4 numShortPerDet
489 );
490 
493  REAL8TimeSeries **sciFlag,
494  const SFTVector *sfts,
495  REAL8 tShort,
496  UINT4 numShortPerDet
497 );
498 
501  MultiREAL8TimeSeries **scienceFlagVect,
502  const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
503  const REAL8 tShort,
504  const UINT4 numShortPerDet
505 );
506 
509  MultiREAL8TimeSeries **scienceFlagVect,
510  const MultiSFTVector *multiSFTs,
511  REAL8 tShort,
512  UINT4 numShortPerDet
513 );
514 
515 int
517  DetectorState *detState,
518  const LALDetector *detector
519 );
520 
524  const LALDetector *detector,
525  const EphemerisData *edat,
526  REAL8 tOffset,
527  REAL8 tShort,
528  UINT4 numShortPerDet
529 );
530 
534  const MultiLALDetector *multiIFO,
535  const EphemerisData *edat,
536  REAL8 tOffset,
537  REAL8 tShort,
538  UINT4 numShortPerDet
539 );
540 
541 
542 int
544  REAL8Vector **resampMultiWeightsX,
545  const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights,
546  const REAL8 tShort,
547  const REAL8 tSFTOld,
548  const UINT4 numShortPerDet,
549  const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
550  const UINT4 maxNumStepsOldIfGapless,
551  const UINT4 X
552 );
553 
554 #ifdef SWIG // SWIG interface directives
555 SWIGLAL( INOUT_STRUCTS( MultiNoiseWeights **, multiWeights ) );
556 #endif
557 int
559  MultiNoiseWeights **multiWeights,
560  const REAL8 tShort,
561  const REAL8 tSFTOld,
562  const UINT4 numShortPerDet,
563  const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes
564 );
565 
566 
567 int
569  MultiAMCoeffs *multiAMcoef,
570  const MultiNoiseWeights *multiWeights,
571  REAL8 tShort,
572  REAL8 tSFTOld,
573  UINT4 numShortPerDet,
574  MultiLIGOTimeGPSVector *multiTimes
575 );
576 
579  const MultiDetectorStateSeries *multiDetStates,
580  const MultiNoiseWeights *multiWeights,
581  SkyPosition skypos,
582  REAL8 tShort,
583  REAL8 tSFTOld,
584  UINT4 numShortPerDet,
585  MultiLIGOTimeGPSVector *multiTimes
586 );
587 
588 AMCoeffs
590  const DetectorStateSeries *DetectorStates,
591  SkyPosition skypos
592 );
593 
594 
595 UINT4
597  const REAL8 resampTshort,
598  const INT4 startTime,
599  const INT4 endTime
600 );
601 
604  LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
605  const LIGOTimeGPSVector *_LAL_RESTRICT_ Times
606 );
607 
610  LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
611  const SFTVector *_LAL_RESTRICT_ sfts
612 );
613 
614 int
616  MultiResampSFTPairMultiIndexList *resampMultiPairs,
617  MultiREAL8TimeSeries *scienceFlagVect
618 );
619 
620 //(possible future function)
621 //MultiSFTVector
622 //*XLALModifyCrossCorrTimestampsIntoSFTVector(
623 // const MultiLIGOTimeGPSVector *multiTimes
624 //);
625 
626 int
629  COMPLEX8 **ws1KFaX_kOut,
630  COMPLEX8 **ws1KFbX_kOut,
631  COMPLEX8 **ws2LFaX_kOut,
632  COMPLEX8 **ws2LFbX_kOut,
633  MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_a,
634  MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_b,
635  const PulsarDopplerParams binaryTemplateSpacings,
636  FstatInput *resampFstatInput,
637  const UINT4 numFreqBins,
638  const REAL8 tCoh,
639  const BOOLEAN treatWarningsAsErrors
640 );
641 
642 /** @} */
643 
644 void XLALDestroySFTIndexList( SFTIndexList *sftIndices );
645 
647 
649 
651 
653 
655 
656 void XLALDestroyMultiMatchList( MultiResampSFTMultiCountList *localMultiListOfLmatchingGivenMultiK );
657 
658 void XLALDestroyResampCrossCorrWorkspace( void *workspace );
659 
660 #ifdef __cplusplus
661 } /* Close C++ protection */
662 #endif
663 
664 
665 #endif /* Close double-include protection _PULSARCROSSCORR_H */
LIGOTimeGPSVector * timestamps
unsigned char BOOLEAN
double REAL8
uint32_t UINT4
float complex COMPLEX8
int32_t INT4
float REAL4
int XLALEquipCrossCorrPairsWithScienceFlags(MultiResampSFTPairMultiIndexList *resampMultiPairs, MultiREAL8TimeSeries *scienceFlagVect)
Demarcate pairs with flags about whether data exists in zero-padded timeseries.
void XLALDestroyResampSFTIndexList(ResampSFTIndexList *sftResampList)
int XLALModifyMultiAMCoeffsWeights(MultiNoiseWeights **multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes)
void XLALDestroySFTPairIndexList(SFTPairIndexList *sftPairs)
Destroy a SFTPairIndexList structure.
void XLALDestroyResampCrossCorrWorkspace(void *workspace)
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 *‍/
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
int XLALCreateSFTPairIndexList(SFTPairIndexList **pairIndexList, SFTIndexList *indexList, MultiSFTVector *sfts, REAL8 maxLag, BOOLEAN inclAutoCorr)
Construct list of SFT pairs for inclusion in statistic.
int XLALCreateSFTPairIndexListShortResamp(MultiResampSFTPairMultiIndexList **resampPairIndexList, const REAL8 maxLag, const BOOLEAN inclAutoCorr, const BOOLEAN inclSameDetector, const REAL8 Tsft, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes)
AMCoeffs * XLALComputeAMCoeffsShort(const DetectorStateSeries *DetectorStates, SkyPosition skypos)
(test function) used for computing amplitude modulation weights
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 XLALCalculateCrossCorrGammasResampShort(REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, MultiResampSFTPairMultiIndexList *resampMultiPairIndexList, MultiAMCoeffs *multiCoeffs)
test function for RESAMPLING with tShort
MultiLIGOTimeGPSVector * XLALModifyMultiTimestampsFromSFTs(MultiREAL8TimeSeries **scienceFlagVect, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet)
int XLALFillDetectorTensorShort(DetectorState *detState, const LALDetector *detector)
(test function) fill detector state with tShort, importing various slightly-modified LALPulsar functi...
int XLALWeightMultiAMCoeffsShort(MultiAMCoeffs *multiAMcoef, const MultiNoiseWeights *multiWeights, REAL8 tShort, REAL8 tSFTOld, UINT4 numShortPerDet, MultiLIGOTimeGPSVector *multiTimes)
(test function) used for weighting multi amplitude modulation coefficients
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....
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)
int XLALCalculateCrossCorrGammas(REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, SFTPairIndexList *pairIndexList, SFTIndexList *indexList, MultiAMCoeffs *multiCoeffs)
Construct vector of G_alpha amplitudes for each SFT pair.
REAL8TimeSeries * XLALCrossCorrGapFinderResamp(LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times)
int XLALTestResampPairIndexList(MultiResampSFTPairMultiIndexList *resampMultiPairIndexList)
Check that the contents of a resampling multi-pair index list are sensible by inspection.
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.
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
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
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.
void XLALDestroySFTIndexList(SFTIndexList *sftIndices)
Destroy a SFTIndexList structure.
void XLALDestroyMultiResampSFTPairMultiIndexList(MultiResampSFTPairMultiIndexList *sftMultiPairsResamp)
void XLALDestroyResampSFTMultiIndexList(ResampSFTMultiIndexList *sftResampMultiList)
int XLALCreateSFTIndexListFromMultiSFTVect(SFTIndexList **indexList, MultiSFTVector *sfts)
Construct flat SFTIndexList out of a MultiSFTVector.
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.
int XLALCalculateCrossCorrGammasResamp(REAL8Vector **Gamma_ave, REAL8Vector **Gamma_circ, MultiResampSFTPairMultiIndexList *resampMultiPairIndexList, MultiAMCoeffs *multiCoeffs)
test function for RESAMPLING
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
void XLALDestroyMultiMatchList(MultiResampSFTMultiCountList *localMultiListOfLmatchingGivenMultiK)
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 u...
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)
UINT4 XLALCrossCorrNumShortPerDetector(const REAL8 resampTshort, const INT4 startTime, const INT4 endTime)
Compute the number of tShort segments per detector.
DetectorStateSeries * XLALGetDetectorStatesShort(const LIGOTimeGPSVector *timestamps, const LALDetector *detector, const EphemerisData *edat, REAL8 tOffset, REAL8 tShort, UINT4 numShortPerDet)
(test function) get detector states for tShort
REAL8TimeSeries * XLALCrossCorrGapFinderResampAlt(LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const SFTVector *_LAL_RESTRICT_ sfts)
LIGOTimeGPSVector * XLALModifyTimestampsFromSFTsShort(REAL8TimeSeries **sciFlag, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times, const REAL8 tShort, const UINT4 numShortPerDet)
void XLALDestroyResampSFTPairMultiIndexList(ResampSFTPairMultiIndexList *sftResampPairMultiList)
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
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 XLALCalculateCrossCorrPhaseDerivativesShort(REAL8VectorSequence **resampPhaseDerivs, const PulsarDopplerParams *dopplerPoint, const EphemerisData *edat, SFTIndexList *indexList, MultiResampSFTPairMultiIndexList *resampMultiPairs, MultiSSBtimes *multiTimes, const DopplerCoordinateSystem *coordSys)
(test function) MODIFIED for Tshort
This structure contains the per-SFT (weighted) antenna-pattern functions , with the SFT-index,...
Definition: LALComputeAM.h:63
A vector of COMPLEX8FrequencySeries.
REAL8 F1NoBuf
Resampling timing 'constant': Fstat time per template per detector for an 'unbuffered' usage (differe...
REAL8 Norm
time spent normalizing the final Fa,Fb
REAL8 F1Buf
Resampling timing 'constant': Fstat time per template per detector for a 'buffered' case (same skypos...
REAL8 FFT
time spent in FFT
REAL8 Bary
time spent in barycentric resampling
REAL8 Total
total time spent in XLALComputeFstatResamp()
REAL8 Spin
time spent in spindown+frequency correction
REAL8 Mem
time to realloc and memset-0 arrays
REAL8 Fab2F
time to compute Fstat from {Fa,Fb}
REAL8 SumFabX
time to sum_X Fab^X
State-info about position, velocity and LMST of a detector together with corresponding EarthState.
Timeseries of DetectorState's, representing the detector-info at different timestamps.
type describing a Doppler coordinate system: lists the number of dimensions and the symbolic names of...
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
A vector of 'timestamps' of type LIGOTimeGPS.
Definition: SFTfileIO.h:188
Multi-IFO container for antenna-pattern coefficients and atenna-pattern matrix .
Definition: LALComputeAM.h:137
Multi-IFO container for COMPLEX8 resampled timeseries.
Definition: LFTandTSutils.h:52
Multi-IFO time-series of DetectorStates.
array of detectors definitions 'LALDetector'
A collection of (multi-IFO) LIGOTimeGPSVector time-stamps vectors.
Definition: SFTfileIO.h:198
One noise-weight (number) per SFT (therefore indexed over IFOs and SFTs.
Definition: PSDutils.h:71
A collection of (multi-IFO) REAL8 time-series.
OUTER List of SFT indices.
ResampSFTMultiCountListDet * data
list per-detectors X
UINT4 length
number of detectors X
Resampling Multi List of all paired SFTs L_Y_K_X, top level (multiple pairs, multiple detectors)
UINT4 oldPairCount
count of sft pairs, old-style
REAL8 Tshort
Duration of resamp Tshort.
REAL8 maxLag
Maximum allowed lag time.
SFTPairIndexList * pairIndexList
Make a classic pair index list.
UINT4 sftTotalCount
count of all sfts
BOOLEAN inclSameDetector
Do we include same detectors?
ResampSFTPairMultiIndexList * data
list per-detector X
UINT4 length
number of detectors X
REAL8 Tsft
Duration of one SFT.
SFTIndexList * indexList
Make an overall flat index list.
UINT4 allPairCount
count of all pairs
BOOLEAN inclAutoCorr
Do we include auto-correlation?
A collection of SFT vectors – one for each IFO in a multi-IFO search.
Definition: SFTfileIO.h:179
Multi-IFO container for SSB timings.
Definition: SSBtimes.h:67
A collection of UINT4Vectors – one for each IFO
UINT4Vector ** data
unit4vector for each ifo
UINT4 length
number of ifos
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
UINT4 numSamplesFFT
actual length of FFT, potentially rounded up to power-of-2 for efficiency
UINT4 numFreqBins
number of frequency bins to compute F-stat for
BOOLEAN collectTiming
turn on/off the collection of F-stat-method-specific timing-data (stored in workspace)
UINT4 numSFTs
total number of SFTs used
UINT4 numSamplesFFT0
'original' length of barycentered timeseries to be FFT'ed
UINT4 numDetectors
number of detectors
COMPLEX8 * FabX_Raw
raw full-band FFT result Fa,Fb
COMPLEX8 * Fa_k
properly normalized F_a(f_k) over output bins
COMPLEX8Vector * TStmp2_SRC
can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
COMPLEX8 * FaX_k
properly normalized F_a^X(f_k) over output bins
REAL8Vector * SRCtimes_DET
holds uniformly-spaced SRC-frame timesteps translated into detector frame [for interpolation]
UINT4 numFreqBinsOut
number of output frequency bins {f_k}
ResampCrossCorrTimingInfo * timingInfo
pointer to storage for collecting timing data (which lives in ResampMethodData)
fftwf_plan fftplan
buffer FFT plan for given numSamplesOut length
COMPLEX8 * TS_FFT
zero-padded, spindown-corr SRC-frame TS
COMPLEX8 * Fb_k
properly normalized F_b(f_k) over output bins
UINT4 numSamplesFFT
allocated number of zero-padded SRC-frame time samples (related to dFreq)
COMPLEX8 * FbX_k
properly normalized F_b^X(f_k) over output bins
UINT4 numFreqBinsAlloc
internal: keep track of allocated length of frequency-arrays
UINT4 decimateFFT
output every n-th frequency bin, with n>1 iff (dFreq > 1/Tspan), and was internally decreased by n
COMPLEX8Vector * TStmp1_SRC
can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
Resampling: Index to refer to fields of an SFT given a specific index L_Y_K_X.
REAL4 sciFlag
science-or-not value: 1 is science mode, 0 is in a gap
UINT4 detInd
index of detector in list
UINT4 pairInd
index of SFT among overall pairs
UINT4 sftInd
index of SFT in list for this detector L
UINT4 flatInd
index in the flat SFT list
Resampling: List of SFT indices L for a given detector Y_K_X: indexing method is nominally original v...
UINT4 length
number of SFTs
ResampSFTIndex * data
array of SFT indices
UINT4 detInd
original vector index of detector Y
MIDDLE List of SFT indices.
UINT4 length
number of SFTs K_X
UINT4 detInd
original vector index of detector X
ResampSFTMultiCountList * data
array of SFT L_Y indices for given SFT K_X at detector X
INNER List of SFT indices.
UINT4 sftInd
original vector index of sft K
UINT4 length
number of detectors Y
SFTCount * data
arrays of count of SFTs L_Y, at given detector Y, that match each given SFT K_X
Resampling: Multi List of indices of SFT L_Y, for a given sft K_X
ResampSFTIndexList * data
array, per detector, of lists of SFT indices
UINT4 sftInd
original vector index of sft K
UINT4 length
number of detectors Y
UINT4 flatInd
index in the flat SFT list
REAL4 sciFlag
science-or-not value: 1 is science mode, 0 is in a gap
Resampling Multi List of SFT pair indices (L_Y_K), for a given detector X.
UINT4 length
number of K SFTs at detector X
ResampSFTMultiIndexList * data
array of SFT L_Y indices for given SFT K_X at detector X
UINT4 detInd
original vector index of detector X
Resampling Counter of matching SFTs for a given detector Y_K_X matching SFT K_X.
UINT4 sftCount
number of matching SFTs
UINT4 detInd
original vector index of detector Y
Index to refer to an SFT given a set of SFTs from several different detectors.
UINT4 sftInd
index of SFT in list for this detector
UINT4 detInd
index of detector in list
List of SFT indices.
UINT4 length
number of SFTs
SFTIndex * data
array of SFT indices
Index to refer to a pair of SFTs.
List of SFT pair indices.
UINT4 length
number of SFT Pairs
SFTPairIndex * data
array of SFT Pair indices