Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
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, 2024 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
26extern "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"
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 */
86typedef 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 */
92typedef struct tagSFTIndexList {
93#ifdef SWIG /* SWIG interface directives */
94 SWIGLAL( ARRAY_1D( SFTIndexList, SFTIndex, data, UINT4, length ) );
95#endif /* SWIG */
96 UINT4 length; /**< number of SFTs */
97 SFTIndex *data; /**< array of SFT indices */
99
100/** Index to refer to a pair of SFTs */
101typedef struct tagSFTPairIndex {
102#if 0
103 SFTIndex sftInd1; /**< index of 1st SFT in pair */
104 SFTIndex sftInd2; /**< index of 2nd SFT in pair */
105#endif
106 UINT4 sftNum[2]; /**< ordinal numbers of first and second SFTs */
108
109/** List of SFT pair indices */
110typedef struct tagSFTPairIndexList {
111#ifdef SWIG /* SWIG interface directives */
112 SWIGLAL( ARRAY_1D( SFTPairIndexList, SFTPairIndex, data, UINT4, length ) );
113#endif /* SWIG */
114 UINT4 length; /**< number of SFT Pairs */
115 SFTPairIndex *data; /**< array of SFT Pair indices */
117
118
119/** Resampling Counter of matching SFTs for a given detector Y_K_X matching SFT K_X */
120typedef struct tagSFTCount {
121 UINT4 detInd; /**< original vector index of detector Y */
122 UINT4 sftCount; /**< number of matching SFTs */
123} SFTCount;
124
125/* Special multi-types for resampling */
126/* These are MULTI-structs */
127
128
129/** INNER List of SFT indices */
130typedef struct tagResampSFTMultiCountList {
131#ifdef SWIG /* SWIG interface directives */
132 SWIGLAL( ARRAY_1D( ResampSFTMultiCountList, SFTCount, data, UINT4, length ) );
133#endif /* SWIG */
134 UINT4 length; /**< number of detectors Y */
135 UINT4 sftInd; /**< original vector index of sft K */
136 SFTCount *data; /**< arrays of count of SFTs L_Y, at given detector Y, that match each given SFT K_X */
138
139/** MIDDLE List of SFT indices */
140typedef struct tagResampSFTMultiCountListDet {
141#ifdef SWIG /* SWIG interface directives */
142 SWIGLAL( ARRAY_1D( ResampSFTMultiCountListDet, ResampSFTMultiCountList, data, UINT4, length ) );
143#endif /* SWIG */
144 UINT4 length; /**< number of SFTs K_X */
145 UINT4 detInd; /**< original vector index of detector X */
146 ResampSFTMultiCountList *data; /**< array of SFT L_Y indices for given SFT K_X at detector X */
148
149/** OUTER List of SFT indices */
150typedef struct tagMultiResampSFTMultiCountList {
151#ifdef SWIG /* SWIG interface directives */
152 SWIGLAL( ARRAY_1D( MultiResampSFTMultiCountList, ResampSFTMultiCountListDet, data, UINT4, length ) );
153#endif /* SWIG */
154 UINT4 length; /**< number of detectors X */
155 ResampSFTMultiCountListDet *data; /**< list per-detectors X */
157
158
159/** Resampling: Index to refer to fields of an SFT given a specific index L_Y_K_X */
160typedef struct tagResampSFTIndex {
161 UINT4 detInd; /**< index of detector in list */
162 UINT4 sftInd; /**< index of SFT in list for this detector L */
163 UINT4 flatInd; /**< index in the flat SFT list */
164 UINT4 pairInd; /**< index of SFT among overall pairs */
165 REAL4 sciFlag; /**< science-or-not value: 1 is science mode, 0 is in a gap */
167
168/** 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 */
169typedef struct tagResampSFTIndexList {
170#ifdef SWIG /* SWIG interface directives */
171 SWIGLAL( ARRAY_1D( ResampSFTIndexList, ResampSFTIndex, data, UINT4, length ) );
172#endif /* SWIG */
173 UINT4 length; /**< number of SFTs */
174 UINT4 detInd; /**< original vector index of detector Y */
175 ResampSFTIndex *data; /**< array of SFT indices */
177
178/** Resampling: Multi List of indices of SFT L_Y, for a given sft K_X */
179typedef struct tagResampSFTMultiIndexList {
180#ifdef SWIG /* SWIG interface directives */
181 SWIGLAL( ARRAY_1D( ResampSFTMultiIndexList, ResampSFTIndexList, data, UINT4, length ) );
182#endif /* SWIG */
183 UINT4 length; /**< number of detectors Y */
184 UINT4 sftInd; /**< original vector index of sft K */
185 UINT4 flatInd; /**< index in the flat SFT list */
186 REAL4 sciFlag; /**< science-or-not value: 1 is science mode, 0 is in a gap */
187 ResampSFTIndexList *data; /**< array, per detector, of lists of SFT indices */
189
190/** Resampling Multi List of SFT pair indices (L_Y_K), for a given detector X */
191typedef struct tagResampSFTPairMultiIndexList {
192#ifdef SWIG /* SWIG interface directives */
193 SWIGLAL( ARRAY_1D( ResampSFTPairMultiIndexList, ResampSFTMultiIndexList, data, UINT4, length ) );
194#endif /* SWIG */
195 UINT4 length; /**< number of K SFTs at detector X */
196 UINT4 detInd; /**< original vector index of detector X */
197 ResampSFTMultiIndexList *data; /**< array of SFT L_Y indices for given SFT K_X at detector X */
199
200
201/** Resampling Multi List of all paired SFTs L_Y_K_X, top level (multiple pairs, multiple detectors) */
202typedef struct tagMultiResampSFTPairMultiIndexList {
203#ifdef SWIG /* SWIG interface directives */
205#endif /* SWIG */
206 UINT4 allPairCount; /**< count of all pairs */
207 UINT4 oldPairCount; /**< count of sft pairs, old-style */
208 UINT4 sftTotalCount; /**< count of all sfts */
209 REAL8 maxLag; /**< Maximum allowed lag time */
210 BOOLEAN inclAutoCorr; /**< Do we include auto-correlation? */
211 BOOLEAN inclSameDetector; /**< Do we include same detectors? */
212 REAL8 Tsft; /**< Duration of one SFT */
213 REAL8 Tshort; /**< Duration of resamp Tshort */
214 SFTIndexList *indexList; /**< Make an overall flat index list */
215 SFTPairIndexList *pairIndexList; /**< Make a classic pair index list */
216 UINT4 length; /**< number of detectors X */
217 ResampSFTPairMultiIndexList *data; /**< list per-detector X */
219
220// ----- local types ----------
221
222typedef struct tagCrossCorrTimings_t {
223 REAL8 Total; //!< total time spent in XLALComputeFstatResamp()
224 REAL8 Bary; //!< time spent in barycentric resampling
225 REAL8 Spin; //!< time spent in spindown+frequency correction
226 REAL8 FFT; //!< time spent in FFT
227 REAL8 Norm; //!< time spent normalizing the final Fa,Fb
228 REAL8 Fab2F; //!< time to compute Fstat from {Fa,Fb}
229 REAL8 Mem; //!< time to realloc and memset-0 arrays
230 REAL8 SumFabX; //!< time to sum_X Fab^X
231 REAL8 F1Buf; //!< Resampling timing 'constant': Fstat time per template per detector for a 'buffered' case (same skypos, same numFreqBins)
232 REAL8 F1NoBuf; //!< Resampling timing 'constant': Fstat time per template per detector for an 'unbuffered' usage (different skypos and numFreqBins)
234
235typedef struct tagResampCrossCorrTimingInfo {
236 // NOTE: all times refer to a single-detector timing case
237 BOOLEAN collectTiming; //!< turn on/off the collection of F-stat-method-specific timing-data (stored in workspace)
238
239 UINT4 numFreqBins; //!< number of frequency bins to compute F-stat for
240 UINT4 numSFTs; //!< total number of SFTs used
241 UINT4 numDetectors; //!< number of detectors
242 UINT4 numSamplesFFT0; //!< 'original' length of barycentered timeseries to be FFT'ed
243 UINT4 numSamplesFFT; //!< actual length of FFT, potentially rounded up to power-of-2 for efficiency
246
247// ----- workspace ----------
248
249typedef struct tagResampCrossCorrWorkspace {
250 // intermediate quantities to interpolate and operate on SRC-frame timeseries
251 COMPLEX8Vector *TStmp1_SRC; //!< can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
252 COMPLEX8Vector *TStmp2_SRC; //!< can hold a single-detector SRC-frame spindown-corrected timeseries [without zero-padding]
253 REAL8Vector *SRCtimes_DET; //!< holds uniformly-spaced SRC-frame timesteps translated into detector frame [for interpolation]
254
255 // input padded timeseries ts(t) and output Fab(f) of length 'numSamplesFFT' and corresponding fftw plan
256 UINT4 numSamplesFFT; //!< allocated number of zero-padded SRC-frame time samples (related to dFreq)
257 UINT4 decimateFFT; //!< output every n-th frequency bin, with n>1 iff (dFreq > 1/Tspan), and was internally decreased by n
258 fftwf_plan fftplan; //!< buffer FFT plan for given numSamplesOut length
259 COMPLEX8 *TS_FFT; //!< zero-padded, spindown-corr SRC-frame TS
260 COMPLEX8 *FabX_Raw; //!< raw full-band FFT result Fa,Fb
261
262 // arrays of size numFreqBinsOut over frequency bins f_k:
263 UINT4 numFreqBinsOut; //!< number of output frequency bins {f_k}
264 COMPLEX8 *FaX_k; //!< properly normalized F_a^X(f_k) over output bins
265 COMPLEX8 *FbX_k; //!< properly normalized F_b^X(f_k) over output bins
266 COMPLEX8 *Fa_k; //!< properly normalized F_a(f_k) over output bins
267 COMPLEX8 *Fb_k; //!< properly normalized F_b(f_k) over output bins
268 UINT4 numFreqBinsAlloc; //!< internal: keep track of allocated length of frequency-arrays
269
270 ResampCrossCorrTimingInfo *timingInfo; //!< pointer to storage for collecting timing data (which lives in ResampMethodData)
272
273/* end Resampling multi-types */
274
275/** A collection of UINT4Vectors -- one for each IFO */
276/* Probably belongs in SFTUtils.h */
277typedef struct tagMultiUINT4Vector {
278#ifdef SWIG /* SWIG interface directives */
279 SWIGLAL( ARRAY_1D( MultiUINT4Vector, UINT4Vector *, data, UINT4, length ) );
280#endif /* SWIG */
281 UINT4 length; /**< number of ifos */
282 UINT4Vector **data; /**< unit4vector for each ifo */
284
285/*
286* Functions Declarations (i.e., prototypes).
287*/
288
290(
291 REAL8Vector *shiftedFreqs,
292 UINT4Vector *lowestBins,
293 COMPLEX8Vector *expSignalPhases,
294 REAL8VectorSequence *sincList,
297 SFTIndexList *sfts,
298 MultiSFTVector *inputSFTs,
299 MultiSSBtimes *multiTimes,
300 MultiUINT4Vector *badBins,
301 REAL8 Tsft
302)
303;
304
306(
307 SFTIndexList **indexList,
308 MultiSFTVector *sfts
309)
310;
311
313(
314 SFTPairIndexList **pairIndexList,
315 SFTIndexList *indexList,
316 MultiSFTVector *sfts,
317 REAL8 maxLag,
318 BOOLEAN inclAutoCorr
319)
320;
321
323(
324 SFTPairIndexList **pairIndexList,
325 UINT4VectorSequence **sftPairForResampPair,
326 SFTIndexList *indexList,
328 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
329 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ resampMultiTimes
330)
331;
332
334(
335 MultiResampSFTPairMultiIndexList **resampPairIndexList,
336 SFTPairIndexList **pairIndexList,
337 SFTIndexList *indexList,
338 MultiSFTVector *sfts,
339 REAL8 maxLag,
340 BOOLEAN inclAutoCorr,
341 BOOLEAN inclSameDetector,
342 REAL8 Tsft,
343 REAL8 Tshort
344)
345;
346
348(
349 MultiResampSFTPairMultiIndexList **resampPairIndexList,
350 const REAL8 maxLag,
351 const BOOLEAN inclAutoCorr,
352 const BOOLEAN inclSameDetector,
353 const REAL8 Tsft,
354 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes
355)
356;
357
359(
360 MultiResampSFTPairMultiIndexList *resampMultiPairIndexList
361);
362
364(
365 REAL8Vector **Gamma_ave,
366 REAL8Vector **Gamma_circ,
367 SFTPairIndexList *pairIndexList,
368 SFTIndexList *indexList,
369 MultiAMCoeffs *multiCoeffs
370)
371;
372
374(
375 REAL8Vector **Gamma_ave,
376 REAL8Vector **Gamma_circ,
377 MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
378 MultiAMCoeffs *multiCoeffs
379)
380;
381
383(
384 REAL8Vector **Gamma_ave,
385 REAL8Vector **Gamma_circ,
386 MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
387 MultiAMCoeffs *multiCoeffs
388)
389;
390
392(
393 REAL8Vector **resampGamma,
394 REAL8Vector *Gamma,
395 UINT4VectorSequence *sftPairForTshortPair,
396 REAL8 Tsft,
397 REAL8 Tshort
398)
399;
400
402(
403 REAL8 *ccStat,
404 REAL8 *evSquared,
405 REAL8Vector *curlyGAmp,
406 COMPLEX8Vector *expSignalPhases,
407 UINT4Vector *lowestBins,
408 REAL8VectorSequence *sincList,
409 SFTPairIndexList *sftPairs,
410 SFTIndexList *sftIndices,
411 MultiSFTVector *inputSFTs,
412 MultiNoiseWeights *multiWeights,
414)
415;
416
418(
419 REAL8Vector *_LAL_RESTRICT_ ccStatVector,
420 REAL8Vector *_LAL_RESTRICT_ evSquaredVector,
421 REAL8Vector *_LAL_RESTRICT_ numeEquivAve,
422 REAL8Vector *_LAL_RESTRICT_ numeEquivCirc,
423 const REAL8Vector *_LAL_RESTRICT_ resampCurlyGAmp,
424 const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampPairs,
425 const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights,
426 const PulsarDopplerParams *_LAL_RESTRICT_ binaryTemplateSpacings,
427 const PulsarDopplerParams *_LAL_RESTRICT_ dopplerpos,
428 const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_a,
429 const MultiCOMPLEX8TimeSeries *_LAL_RESTRICT_ multiTimeSeries_SRC_b,
430 ResampCrossCorrWorkspace *_LAL_RESTRICT_ ws,
431 COMPLEX8 *_LAL_RESTRICT_ ws1KFaX_k,
432 COMPLEX8 *_LAL_RESTRICT_ ws1KFbX_k,
433 COMPLEX8 *_LAL_RESTRICT_ ws2LFaX_k,
434 COMPLEX8 *_LAL_RESTRICT_ ws2LFbX_k
435)
436;
437
439(
440 REAL8VectorSequence **phaseDerivs,
441 const PulsarDopplerParams *dopplerPoint,
442 const EphemerisData *edat,
443 SFTIndexList *indexList,
444 MultiSSBtimes *multiTimes,
445 const DopplerCoordinateSystem *coordSys
446)
447;
448
450(
451 REAL8VectorSequence **resampPhaseDerivs,
452 const PulsarDopplerParams *dopplerPoint,
453 const EphemerisData *edat,
454 SFTIndexList *indexList,
455 MultiResampSFTPairMultiIndexList *resampMultiPairIndexList,
456 MultiSSBtimes *multiTimes,
457 const DopplerCoordinateSystem *coordSys
458)
459;
460
462(
463 gsl_matrix **g_ij,
464 gsl_vector **eps_i,
465 REAL8 *sumGammaSq,
466 const REAL8VectorSequence *phaseDerivs,
467 const SFTPairIndexList *pairIndexList,
468 const REAL8Vector *Gamma_ave,
469 const REAL8Vector *Gamma_circ,
470 const DopplerCoordinateSystem *coordSys
471);
472
474(
475 gsl_matrix **g_ij,
476 gsl_vector **eps_i,
477 REAL8 *sumGammaSq,
478 const REAL8VectorSequence *phaseDerivs,
479 const MultiResampSFTPairMultiIndexList *resampMultiPairs,
480 const REAL8Vector *Gamma_ave,
481 const REAL8Vector *Gamma_circ,
482 const DopplerCoordinateSystem *coordSys
483);
484
486(
487 REAL8 *hSens,
488 REAL8 *g_ff,
489 REAL8 *g_aa,
490 REAL8 *g_TT,
491 REAL8 *g_pp,
492 REAL8 *weightedMuTAve,
493 PulsarDopplerParams DopplerParams,
494 REAL8Vector *G_alpha,
495 SFTPairIndexList *pairIndexList,
496 SFTIndexList *indexList,
497 MultiSFTVector *sfts,
498 MultiNoiseWeights *multiWeights
499)
500;
501
503(
504 REAL8 *hSens,
505 REAL8 *g_ff,
506 REAL8 *g_aa,
507 REAL8 *g_TT,
508 REAL8 *g_pp,
509 const PulsarDopplerParams DopplerParams,
510 const REAL8Vector *_LAL_RESTRICT_ G_alpha,
511 const MultiResampSFTPairMultiIndexList *_LAL_RESTRICT_ resampMultiPairIndexList,
512 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
513 const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights
514)
515;
516
517/* (possible future function) */
518//int
519//XLALBesselCrossCorrOrbitalSpaceStep(
520// COMPLEX8 * xTildeOut,
521// COMPLEX8 * xTildeIn,
522// PulsarDopplerParams * dopplerpos,
523// PulsarDopplerParams * binaryTemplateSpacings,
524// INT4 aStep,
525// INT4 pStep,
526// INT4 tStep,
527// INT4 nFFT,
528// UINT4 nTermsMax
529//);
530
533 const REAL8 tShort,
534 const UINT4 numShortPerDet,
535 const LIGOTimeGPS epoch
536);
537
538/* DEPRECATED function; use XLALGenerateTshortTimestamps() instead */
539
542 REAL8TimeSeries **sciFlag,
543 const LIGOTimeGPSVector *_LAL_RESTRICT_ Times,
544 const REAL8 tShort,
545 const UINT4 numShortPerDet
546);
547
550 REAL8TimeSeries **sciFlag,
551 const SFTVector *sfts,
552 REAL8 tShort,
553 UINT4 numShortPerDet
554);
555
556
559 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
560 const REAL8 tShort,
561 const UINT4 numShortPerDet,
562 const BOOLEAN alignTShorts
563);
564
565/* DEPRECATED function; use XLALGenerateMultiTshortTimestamps() instead */
568 MultiREAL8TimeSeries **scienceFlagVect,
569 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
570 const REAL8 tShort,
571 const UINT4 numShortPerDet
572);
573
576 MultiREAL8TimeSeries **scienceFlagVect,
578 REAL8 tShort,
579 UINT4 numShortPerDet
580);
581
582int
584 DetectorState *detState,
585 const LALDetector *detector
586);
587
591 const LALDetector *detector,
592 const EphemerisData *edat,
593 REAL8 tOffset,
594 REAL8 tShort,
595 UINT4 numShortPerDet
596);
597
601 const MultiLALDetector *multiIFO,
602 const EphemerisData *edat,
603 REAL8 tOffset,
604 REAL8 tShort,
605 UINT4 numShortPerDet
606);
607
608
609int
611 REAL8Vector **resampMultiWeightsX,
612 const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights,
613 const REAL8 tShort,
614 const REAL8 tSFTOld,
615 const UINT4 numShortPerDet,
616 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes,
617 const UINT4 maxNumStepsOldIfGapless,
618 const UINT4 X
619);
620
623 const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights,
624 const REAL8 tShort,
625 const REAL8 tSFTOld,
626 const UINT4 numShortPerDet,
627 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes
628);
629
630/* DEPRECATED function; use XLALModifyMultiWeights() instead */
631#ifdef SWIG // SWIG interface directives
632SWIGLAL( INOUT_STRUCTS( MultiNoiseWeights **, multiWeights ) );
633#endif
634int
636 MultiNoiseWeights **multiWeights,
637 const REAL8 tShort,
638 const REAL8 tSFTOld,
639 const UINT4 numShortPerDet,
640 const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes
641);
642
643int
645 MultiAMCoeffs *multiAMcoef,
646 const MultiNoiseWeights *multiWeights,
647 REAL8 tShort,
648 REAL8 tSFTOld,
649 UINT4 numShortPerDet,
650 MultiLIGOTimeGPSVector *multiTimes
651);
652
655 const MultiDetectorStateSeries *multiDetStates,
656 const MultiNoiseWeights *multiWeights,
657 SkyPosition skypos,
658 REAL8 tShort,
659 REAL8 tSFTOld,
660 UINT4 numShortPerDet,
661 MultiLIGOTimeGPSVector *multiTimes
662);
663
666 const DetectorStateSeries *DetectorStates,
667 SkyPosition skypos
668);
669
670
671UINT4
673 const REAL8 resampTshort,
674 const INT4 startTime,
675 const INT4 endTime
676);
677
680 LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
681 const LIGOTimeGPSVector *_LAL_RESTRICT_ Times
682);
683
686 LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps,
687 const SFTVector *_LAL_RESTRICT_ sfts
688);
689
690int
692 MultiResampSFTPairMultiIndexList *resampMultiPairs,
693 MultiREAL8TimeSeries *scienceFlagVect
694);
695
696//(possible future function)
697//MultiSFTVector
698//*XLALModifyCrossCorrTimestampsIntoSFTVector(
699// const MultiLIGOTimeGPSVector *multiTimes
700//);
701
702int
705 COMPLEX8 **ws1KFaX_kOut,
706 COMPLEX8 **ws1KFbX_kOut,
707 COMPLEX8 **ws2LFaX_kOut,
708 COMPLEX8 **ws2LFbX_kOut,
709 MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_a,
710 MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_b,
711 const PulsarDopplerParams binaryTemplateSpacings,
712 FstatInput *resampFstatInput,
713 const UINT4 numFreqBins,
714 const REAL8 tCoh,
715 const BOOLEAN treatWarningsAsErrors
716);
717
718/** @} */
719
720void XLALDestroySFTIndexList( SFTIndexList *sftIndices );
721
723
725
727
729
731
732void XLALDestroyMultiMatchList( MultiResampSFTMultiCountList *localMultiListOfLmatchingGivenMultiK );
733
734void XLALDestroyResampCrossCorrWorkspace( void *workspace );
735
736#ifdef __cplusplus
737} /* Close C++ protection */
738#endif
739
740
741#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
MultiLIGOTimeGPSVector * XLALGenerateMultiTshortTimestamps(const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet, const BOOLEAN alignTShorts)
REAL8TimeSeries * XLALCrossCorrGapFinderResampAlt(LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const SFTVector *_LAL_RESTRICT_ sfts)
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
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.
MultiLIGOTimeGPSVector * XLALModifyMultiTimestampsFromSFTs(MultiREAL8TimeSeries **scienceFlagVect, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const REAL8 tShort, const UINT4 numShortPerDet)
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.
LIGOTimeGPSVector * XLALGenerateTshortTimestamps(const REAL8 tShort, const UINT4 numShortPerDet, const LIGOTimeGPS epoch)
Generate timestamps for one detector with tShort.
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 XLALCombineCrossCorrGammas(REAL8Vector **resampGamma, REAL8Vector *Gamma, UINT4VectorSequence *sftPairForTshortPair, REAL8 Tsft, REAL8 Tshort)
Combine G_alpha amplitudes for SFT pairs into amplitudes for Tshort pairs.
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)
MultiNoiseWeights * XLALModifyMultiWeights(const MultiNoiseWeights *_LAL_RESTRICT_ multiWeights, const REAL8 tShort, const REAL8 tSFTOld, const UINT4 numShortPerDet, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes)
int XLALCreateSFTPairIndexListAccurateResamp(SFTPairIndexList **pairIndexList, UINT4VectorSequence **sftPairForResampPair, SFTIndexList *indexList, MultiResampSFTPairMultiIndexList *resampPairs, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ multiTimes, const MultiLIGOTimeGPSVector *_LAL_RESTRICT_ resampMultiTimes)
LIGOTimeGPSVector * XLALModifyTimestampsFromSFTsShort(REAL8TimeSeries **sciFlag, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times, const REAL8 tShort, const UINT4 numShortPerDet)
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
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
REAL8TimeSeries * XLALCrossCorrGapFinderResamp(LIGOTimeGPSVector *_LAL_RESTRICT_ timestamps, const LIGOTimeGPSVector *_LAL_RESTRICT_ Times)
float data[BLOCKSIZE]
Definition: hwinject.c:360
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