LALPulsar  6.1.0.1-bede9b2
Header ComputeFstat.h

Detailed Description

The \( \mathcal{F} \) -statistic.

Authors
Badri Krishnan, Bernd Machenschalk, Chris Messenger, David Keitel, Holger Pletsch, John T. Whelan, Karl Wette, Pinkesh Patel, Reinhard Prix, Xavier Siemens

This module provides a API for computing the \( \mathcal{F} \) -statistic [9] , using various different methods. All data required to compute the \( \mathcal{F} \) -statistic are contained in the opaque structure FstatInput, which is shared by all methods. A function XLALCreateFstatInput() is provided for creating an FstatInput structure configured for the particular method. The FstatInput structure is passed to the function XLALComputeFstat(), which computes the \( \mathcal{F} \) -statistic using the chosen method, and fills a FstatResults structure with the results.

Note
The \( \mathcal{F} \) -statistic method codes are partly descended from earlier implementations found in:
  • LALDemod.[ch] by Jolien Creighton, Maria Alessandra Papa, Reinhard Prix, Steve Berukoff, Xavier Siemens, Bruce Allen
  • ComputeSky.[ch] by Jolien Creighton, Reinhard Prix, Steve Berukoff
  • LALComputeAM.[ch] by Jolien Creighton, Maria Alessandra Papa, Reinhard Prix, Steve Berukoff, Xavier Siemens
  • ComputeFStatistic_resamp.c by Pinkesh Patel, Xavier Siemens, Reinhard Prix, Iraj Gholami, Yousuke Itoh, Maria Alessandra Papa

Modules

 Module ComputeFstat_Demod.c
 Implements the Demod Dirichlet kernel-based demodulation algorithm for computing the \( \mathcal{F} \) -statistic [36] .
 
 Module ComputeFstat_Resamp_CUDA.cu
 Implements a CUDA version [7] of the Resamp FFT-based resampling algorithm for computing the \(\mathcal{F}\)-statistic [9] .
 
 Module ComputeFstat_Resamp_Generic.c
 Implements a generic version [22] of the Resamp FFT-based resampling algorithm for computing the \( \mathcal{F} \) -statistic [9] .
 

Prototypes

REAL8 XLALFstatMaximumSFTLength (const REAL8 maxFreq, const REAL8 binaryMaxAsini, const REAL8 binaryMinPeriod, const REAL8 mu_SFT)
 Compute the maximum SFT length that can safely be used as input to XLALComputeFstat(), given the desired range limits in search parameters. More...
 
int XLALFstatCheckSFTLengthMismatch (const REAL8 Tsft, const REAL8 maxFreq, const REAL8 binaryMaxAsini, const REAL8 binaryMinPeriod, const REAL8 allowedMismatch)
 Check that the SFT length \( T_{\textrm{SFT}} \) does not exceed the result of XLALFstatMaximumSFTLength(), in which case a large unexpected mismatch would be produced. More...
 
FstatInputVectorXLALCreateFstatInputVector (const UINT4 length)
 Create a FstatInputVector of the given length, for example for setting up F-stat searches over several segments. More...
 
void XLALDestroyFstatInputVector (FstatInputVector *inputs)
 Free all memory associated with a FstatInputVector structure. More...
 
FstatAtomVectorXLALCreateFstatAtomVector (const UINT4 length)
 Create a FstatAtomVector of the given length. More...
 
void XLALDestroyFstatAtomVector (FstatAtomVector *atoms)
 Free all memory associated with a FstatAtomVector structure. More...
 
MultiFstatAtomVectorXLALCreateMultiFstatAtomVector (const UINT4 length)
 Create a MultiFstatAtomVector of the given length. More...
 
void XLALDestroyMultiFstatAtomVector (MultiFstatAtomVector *multiAtoms)
 Free all memory associated with a MultiFstatAtomVector structure. More...
 
FstatInput * XLALCreateFstatInput (const SFTCatalog *SFTcatalog, const REAL8 minCoverFreq, const REAL8 maxCoverFreq, const REAL8 dFreq, const EphemerisData *ephemerides, const FstatOptionalArgs *optionalArgs)
 Create a fully-setup FstatInput structure for computing the \( \mathcal{F} \) -statistic using XLALComputeFstat(). More...
 
int XLALGetFstatInputSFTBand (const FstatInput *input, REAL8 *minFreqFull, REAL8 *maxFreqFull)
 Returns the frequency band loaded from input SFTs. More...
 
const CHARXLALGetFstatInputMethodName (const FstatInput *input)
 Returns the human-readable name of the \( \mathcal{F} \) -statistic method being used by a FstatInput structure. More...
 
const MultiLALDetectorXLALGetFstatInputDetectors (const FstatInput *input)
 Returns the detector information stored in a FstatInput structure. More...
 
const MultiLIGOTimeGPSVectorXLALGetFstatInputTimestamps (const FstatInput *input)
 Returns the SFT timestamps stored in a FstatInput structure. More...
 
MultiNoiseWeightsXLALGetFstatInputNoiseWeights (const FstatInput *input)
 Returns the multi-detector noise weights stored in a FstatInput structure. More...
 
const MultiDetectorStateSeriesXLALGetFstatInputDetectorStates (const FstatInput *input)
 Returns the multi-detector state series stored in a FstatInput structure. More...
 
int XLALComputeFstat (FstatResults **Fstats, FstatInput *input, const PulsarDopplerParams *doppler, const UINT4 numFreqBins, const FstatQuantities whatToCompute)
 Compute the \( \mathcal{F} \) -statistic over a band of frequencies. More...
 
void XLALDestroyFstatInput (FstatInput *input)
 Free all memory associated with a FstatInput structure. More...
 
void XLALDestroyFstatResults (FstatResults *Fstats)
 Free all memory associated with a FstatResults structure. More...
 
REAL4 XLALComputeFstatFromFaFb (COMPLEX8 Fa, COMPLEX8 Fb, REAL4 A, REAL4 B, REAL4 C, REAL4 E, REAL4 Dinv)
 Compute the \( \mathcal{F} \) -statistic from the complex \( F_a \) and \( F_b \) components and the antenna pattern matrix. More...
 
REAL4 XLALComputeFstatFromAtoms (const MultiFstatAtomVector *multiFstatAtoms, const INT4 X)
 Compute single-or multi-IFO Fstat '2F' from multi-IFO 'atoms'. More...
 
static int XLALSelectBestFstatMethod (FstatMethodType *method)
 If user asks for a 'best' FstatMethodType, find and select it. More...
 
int XLALFstatMethodIsAvailable (FstatMethodType method)
 Return true if given FstatMethodType corresponds to a valid and available Fstat method, false otherwise. More...
 
const CHARXLALFstatMethodName (FstatMethodType method)
 Return pointer to a static string giving the name of the FstatMethodType method. More...
 
const UserChoicesXLALFstatMethodChoices (void)
 Return pointer to a static array of all (available) FstatMethodType choices. More...
 
int XLALGetFstatTiming (const FstatInput *input, FstatTimingGeneric *timingGeneric, FstatTimingModel *timingModel)
 Return measured values and details about generic F-statistic timing and method-specific timing model,. More...
 
int XLALAppendFstatTiming2File (const FstatInput *input, FILE *fp, BOOLEAN printHeader)
 
int XLALExtractResampledTimeseries (MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_a, MultiCOMPLEX8TimeSeries **multiTimeSeries_SRC_b, FstatInput *input)
 Extracts the resampled timeseries from a given FstatInput structure (must have been initialized previously by XLALCreateFstatInput() and a call to XLALComputeFstat()). More...
 
int XLALFstatInputTimeslice (FstatInput **slice, const FstatInput *input, const LIGOTimeGPS *minStartGPS, const LIGOTimeGPS *maxStartGPS)
 Create and return an FstatInput 'timeslice' for given input FstatInput object [must be using LALDemod Fstat method!] which covers only the time interval ['minStartGPS','maxStartGPS'), using conventions of XLALCWGPSinRange(). More...
 
static void XLALDestroyFstatInputTimeslice_common (FstatCommon *common)
 

Data Structures

struct  FstatInputVector
 A vector of XLALComputeFstat() input data structures, for e.g. More...
 
struct  FstatOptionalArgs
 Struct of optional 'advanced level' and (potentially method-specific) arguments to be passed to the \( \mathcal{F} \) -statistic setup function XLALCreateFstatInput(). More...
 
struct  FstatAtom
 An \( \mathcal{F} \) -statistic 'atom', i.e. More...
 
struct  FstatAtomVector
 A vector of \( \mathcal{F} \) -statistic 'atoms', i.e. More...
 
struct  MultiFstatAtomVector
 A multi-detector vector of FstatAtomVector. More...
 
struct  FstatResults
 XLALComputeFstat() computed results structure. More...
 
struct  FstatTimingGeneric
 Generic F-stat timing coefficients (times in seconds) [see https://dcc.ligo.org/LIGO-T1600531-v4 for details] tauF_eff = tauF_core + b * tauF_buffer where 'b' denotes the fraction of times per call to XLALComputeFstat() the buffer needed to be recomputed ==> b = NBufferMisses / NCalls
". More...
 
struct  FstatTimingModel
 Struct to carry the \( \mathcal{F} \) -statistic method-specific timing model in terms of a list of variable names and corresponding REAL4 values, including a help-string documenting the timing model variables. More...
 

Enumerations

enum  FstatQuantities {
  FSTATQ_NONE = 0x00 , FSTATQ_2F = 0x01 , FSTATQ_FAFB = 0x02 , FSTATQ_2F_PER_DET = 0x04 ,
  FSTATQ_FAFB_PER_DET = 0x08 , FSTATQ_ATOMS_PER_DET = 0x10 , FSTATQ_2F_CUDA = 0x20 , FSTATQ_LAST = 0x80
}
 Bit-field of \( \mathcal{F} \) -statistic quantities which can be computed by XLALComputeFstat(). More...
 
enum  FstatMethodType {
  FMETHOD_DEMOD_GENERIC , FMETHOD_DEMOD_OPTC , FMETHOD_DEMOD_ALTIVEC , FMETHOD_DEMOD_SSE ,
  FMETHOD_DEMOD_BEST , FMETHOD_RESAMP_GENERIC , FMETHOD_RESAMP_CUDA , FMETHOD_RESAMP_BEST
}
 Different methods available to compute the \( \mathcal{F} \) -statistic, falling into two broad classes: More...
 

Macros

#define DEFAULT_MAX_MISMATCH_FROM_SFT_LENGTH   0.05
 default maximum allowed F-stat mismatch from SFTs being too long, to be used in XLALFstatCheckSFTLengthMismatch() More...
 
#define TIMING_MODEL_MAX_VARS   10
 

Variables

const FstatOptionalArgs FstatOptionalArgsDefaults
 Global initializer for setting FstatOptionalArgs to default values. More...
 

Function Documentation

◆ XLALFstatMaximumSFTLength()

REAL8 XLALFstatMaximumSFTLength ( const REAL8  maxFreq,
const REAL8  binaryMaxAsini,
const REAL8  binaryMinPeriod,
const REAL8  mu_SFT 
)

Compute the maximum SFT length that can safely be used as input to XLALComputeFstat(), given the desired range limits in search parameters.

The \( \mathcal{F} \) -statistic algorithms implemented in this module assume that the input SFTs supplied to XLALCreateFstatInput() are of such a length that, within the time span of each SFT, the spectrum of any CW signal being searched for will be mostly contained within one SFT bin:

  • The Demod algorithm makes explicit use of this assumption in order to sum up only a few bins in the Dirichlet kernel, thereby speeding up the computation.
  • While the Resamp algorithm is in principle independent of the SFT length (as the data are converted from SFTs back into a time series), in practise it makes use of this assumption for convenience to linearly approximate the phase over the time span of a single SFT.

In order for this assumption to be valid, the SFT bins must be of a certain minimum size, and hence the time spans of the SFTs must be of a certain maximum length. This maximum length is dependent upon the parameter space being searched: searches for isolated CW sources rarely encounter this restriction by using standard 1800-second SFTs, but searches for CW sources in binary systems typically must use shorter SFTs.

An expression giving the maximum allowed SFT length is given by [16] , eq. (C2):

\[ T_{\textrm{SFT-max}} = \sqrt{ \frac{ 6 \sqrt{ 5 \mu_{\textrm{SFT}} } }{ \pi (a \sin\iota / c)_{\textrm{max}} f_{\textrm{max}} \Omega_{\textrm{max}}^2 } } \,, \]

This function computes this expression with a user-given acceptable maximum mismatch \( \mu_{\textrm{SFT}} \) , and \( (a \sin\iota / c)_{\textrm{max}} \) and \( \Omega_{\textrm{max}} \) set to either the binary motion of the source, as specified by binaryMaxAsini and binaryMinPeriod, or the sidereal motion of the Earth, i.e. with \( (a \sin\iota / c) \sim 0.02~\textrm{s} \) and \( \Omega \sim 2\pi / 1~\textrm{day} \) .

Parameters
[in]maxFreqMaximum signal frequency
[in]binaryMaxAsiniMaximum projected semi-major axis a*sini/c (= 0 for isolated sources)
[in]binaryMinPeriodMinimum orbital period (s)
[in]mu_SFTMaximum allowed fractional mismatch

Definition at line 163 of file ComputeFstat.c.

◆ XLALFstatCheckSFTLengthMismatch()

int XLALFstatCheckSFTLengthMismatch ( const REAL8  Tsft,
const REAL8  maxFreq,
const REAL8  binaryMaxAsini,
const REAL8  binaryMinPeriod,
const REAL8  allowedMismatch 
)

Check that the SFT length \( T_{\textrm{SFT}} \) does not exceed the result of XLALFstatMaximumSFTLength(), in which case a large unexpected mismatch would be produced.

See documentation of that function for parameter dependence. If allowedMismatch==0, the default value is taken from DEFAULT_MAX_MISMATCH_FROM_SFT_LENGTH .

Parameters
[in]Tsftactual SFT length
[in]maxFreqMaximum signal frequency
[in]binaryMaxAsiniMaximum projected semi-major axis a*sini/c (= 0 for isolated sources)
[in]binaryMinPeriodMinimum orbital period (s)
[in]allowedMismatchMaximum allowed fractional mismatch

Definition at line 203 of file ComputeFstat.c.

◆ XLALCreateFstatInputVector()

FstatInputVector * XLALCreateFstatInputVector ( const UINT4  length)

Create a FstatInputVector of the given length, for example for setting up F-stat searches over several segments.

Parameters
[in]lengthLength of the FstatInputVector.

Definition at line 231 of file ComputeFstat.c.

◆ XLALDestroyFstatInputVector()

void XLALDestroyFstatInputVector ( FstatInputVector inputs)

Free all memory associated with a FstatInputVector structure.

Parameters
[in]inputsFstatInputVector structure to be freed.

Definition at line 252 of file ComputeFstat.c.

◆ XLALCreateFstatAtomVector()

FstatAtomVector * XLALCreateFstatAtomVector ( const UINT4  length)

Create a FstatAtomVector of the given length.

Parameters
[in]lengthLength of the FstatAtomVector.

Definition at line 276 of file ComputeFstat.c.

◆ XLALDestroyFstatAtomVector()

void XLALDestroyFstatAtomVector ( FstatAtomVector atoms)

Free all memory associated with a FstatAtomVector structure.

Parameters
[in]atomsFstatAtomVector structure to be freed.

Definition at line 297 of file ComputeFstat.c.

◆ XLALCreateMultiFstatAtomVector()

MultiFstatAtomVector * XLALCreateMultiFstatAtomVector ( const UINT4  length)

Create a MultiFstatAtomVector of the given length.

Parameters
[in]lengthLength of the MultiFstatAtomVector.

Definition at line 317 of file ComputeFstat.c.

◆ XLALDestroyMultiFstatAtomVector()

void XLALDestroyMultiFstatAtomVector ( MultiFstatAtomVector multiAtoms)

Free all memory associated with a MultiFstatAtomVector structure.

Parameters
[in]multiAtomsMultiFstatAtomVector structure to be freed.

Definition at line 338 of file ComputeFstat.c.

◆ XLALCreateFstatInput()

FstatInput * XLALCreateFstatInput ( const SFTCatalog SFTcatalog,
const REAL8  minCoverFreq,
const REAL8  maxCoverFreq,
const REAL8  dFreq,
const EphemerisData ephemerides,
const FstatOptionalArgs optionalArgs 
)

Create a fully-setup FstatInput structure for computing the \( \mathcal{F} \) -statistic using XLALComputeFstat().

Parameters
[in]SFTcatalogCatalog of SFTs to either load from files, or generate in memory. The locator field of each SFTDescriptor must be !=NULL for SFT loading, and ==NULL for SFT generation. *
[in]minCoverFreqMinimum instantaneous frequency which will be covered over the SFT time span.
[in]maxCoverFreqMaximum instantaneous frequency which will be covered over the SFT time span.
[in]dFreqRequested spacing of \( \mathcal{F} \) -statistic frequency bins. May be zero only for single-frequency searches.
[in]ephemeridesEphemerides for the time-span of the SFTs.
[in]optionalArgsOptional 'advanced-level' and method-specific extra arguments; NULL: use defaults from FstatOptionalArgsDefaults.

Definition at line 359 of file ComputeFstat.c.

◆ XLALGetFstatInputSFTBand()

int XLALGetFstatInputSFTBand ( const FstatInput *  input,
REAL8 minFreqFull,
REAL8 maxFreqFull 
)

Returns the frequency band loaded from input SFTs.

Parameters
[in]inputFstatInput structure.
[out]minFreqFullMinimum frequency loaded from input SFTs
[out]maxFreqFullMaximum frequency loaded from input SFTs

Definition at line 650 of file ComputeFstat.c.

◆ XLALGetFstatInputMethodName()

const CHAR * XLALGetFstatInputMethodName ( const FstatInput *  input)

Returns the human-readable name of the \( \mathcal{F} \) -statistic method being used by a FstatInput structure.

Parameters
[in]inputFstatInput structure.

Definition at line 672 of file ComputeFstat.c.

◆ XLALGetFstatInputDetectors()

const MultiLALDetector * XLALGetFstatInputDetectors ( const FstatInput *  input)

Returns the detector information stored in a FstatInput structure.

Parameters
[in]inputFstatInput structure.

Definition at line 687 of file ComputeFstat.c.

◆ XLALGetFstatInputTimestamps()

const MultiLIGOTimeGPSVector * XLALGetFstatInputTimestamps ( const FstatInput *  input)

Returns the SFT timestamps stored in a FstatInput structure.

Parameters
[in]inputFstatInput structure.

Definition at line 701 of file ComputeFstat.c.

◆ XLALGetFstatInputNoiseWeights()

MultiNoiseWeights * XLALGetFstatInputNoiseWeights ( const FstatInput *  input)

Returns the multi-detector noise weights stored in a FstatInput structure.

Note: the returned object is allocated here and needs to be free'ed by caller.

Parameters
[in]inputFstatInput structure.

Definition at line 716 of file ComputeFstat.c.

◆ XLALGetFstatInputDetectorStates()

const MultiDetectorStateSeries * XLALGetFstatInputDetectorStates ( const FstatInput *  input)

Returns the multi-detector state series stored in a FstatInput structure.

Parameters
[in]inputFstatInput structure.

Definition at line 746 of file ComputeFstat.c.

◆ XLALComputeFstat()

int XLALComputeFstat ( FstatResults **  Fstats,
FstatInput *  input,
const PulsarDopplerParams doppler,
const UINT4  numFreqBins,
const FstatQuantities  whatToCompute 
)

Compute the \( \mathcal{F} \) -statistic over a band of frequencies.

Parameters
Fstats[in/out] Address of a pointer to a FstatResults results structure; if NULL, allocate here.
[in]inputInput data structure created by one of the setup functions.
[in]dopplerDoppler parameters, including starting frequency, at which to compute \( 2\mathcal{F} \)
[in]numFreqBinsNumber of frequencies at which the \( 2\mathcal{F} \) are to be computed. Must be 1 if XLALCreateFstatInput() was passed zero dFreq.
[in]whatToComputeBit-field of which \( \mathcal{F} \) -statistic quantities to compute.

Definition at line 760 of file ComputeFstat.c.

◆ XLALDestroyFstatInput()

void XLALDestroyFstatInput ( FstatInput *  input)

Free all memory associated with a FstatInput structure.

Parameters
[in]inputFstatInput structure to be freed.

Definition at line 885 of file ComputeFstat.c.

◆ XLALDestroyFstatResults()

void XLALDestroyFstatResults ( FstatResults Fstats)

Free all memory associated with a FstatResults structure.

Parameters
[in]FstatsFstatResults structure to be freed.

Definition at line 934 of file ComputeFstat.c.

◆ XLALComputeFstatFromFaFb()

REAL4 XLALComputeFstatFromFaFb ( COMPLEX8  Fa,
COMPLEX8  Fb,
REAL4  A,
REAL4  B,
REAL4  C,
REAL4  E,
REAL4  Dinv 
)

Compute the \( \mathcal{F} \) -statistic from the complex \( F_a \) and \( F_b \) components and the antenna pattern matrix.

Note for developers: This is just a wrapper to the internal function compute_fstat_from_fa_fb() so that it can be re-used externally. Inside this module, better use that function directly!

Parameters
[in]FaF-stat component \( F_a \)
[in]FbF-stat component \( F_b \)
[in]Aantenna pattern matrix component
[in]Bantenna pattern matrix component
[in]Cantenna pattern matrix component
[in]Eantenna pattern matrix component
[in]Dinvinverse determinant

Definition at line 977 of file ComputeFstat.c.

◆ XLALComputeFstatFromAtoms()

REAL4 XLALComputeFstatFromAtoms ( const MultiFstatAtomVector multiFstatAtoms,
const INT4  X 
)

Compute single-or multi-IFO Fstat '2F' from multi-IFO 'atoms'.

Parameters
[in]multiFstatAtomsMulti-detector atoms
[in]XDetector number, give -1 for multi-Fstat

Definition at line 994 of file ComputeFstat.c.

◆ XLALSelectBestFstatMethod()

static int XLALSelectBestFstatMethod ( FstatMethodType method)
static

If user asks for a 'best' FstatMethodType, find and select it.

Definition at line 1052 of file ComputeFstat.c.

◆ XLALFstatMethodIsAvailable()

int XLALFstatMethodIsAvailable ( FstatMethodType  method)

Return true if given FstatMethodType corresponds to a valid and available Fstat method, false otherwise.

Definition at line 1085 of file ComputeFstat.c.

◆ XLALFstatMethodName()

const CHAR * XLALFstatMethodName ( FstatMethodType  method)

Return pointer to a static string giving the name of the FstatMethodType method.

Definition at line 1135 of file ComputeFstat.c.

◆ XLALFstatMethodChoices()

const UserChoices * XLALFstatMethodChoices ( void  )

Return pointer to a static array of all (available) FstatMethodType choices.

This data is used by the UserInput module to parse a user enumeration.

Definition at line 1147 of file ComputeFstat.c.

◆ XLALGetFstatTiming()

int XLALGetFstatTiming ( const FstatInput *  input,
FstatTimingGeneric timingGeneric,
FstatTimingModel timingModel 
)

Return measured values and details about generic F-statistic timing and method-specific timing model,.

See also https://dcc.ligo.org/LIGO-T1600531-v4 for details about the timing model.

Definition at line 1190 of file ComputeFstat.c.

◆ XLALAppendFstatTiming2File()

int XLALAppendFstatTiming2File ( const FstatInput *  input,
FILE *  fp,
BOOLEAN  printHeader 
)

Definition at line 1225 of file ComputeFstat.c.

◆ XLALExtractResampledTimeseries()

int XLALExtractResampledTimeseries ( MultiCOMPLEX8TimeSeries **  multiTimeSeries_SRC_a,
MultiCOMPLEX8TimeSeries **  multiTimeSeries_SRC_b,
FstatInput *  input 
)

Extracts the resampled timeseries from a given FstatInput structure (must have been initialized previously by XLALCreateFstatInput() and a call to XLALComputeFstat()).

Note
This function only returns a pointer to the time-series that are still part of the FstatInput structure, and will by managed by F-statistic API calls. Therefore do not free the resampled timeseries with XLALDestroyMultiCOMPLEX8TimeSeries(), only free the complete Fstat input structure with XLALDestroyFstatInputVector() at the end once its no longer needed.
Parameters
[out]multiTimeSeries_SRC_amulti-detector SRC-frame timeseries, multiplied by AM function a(t).
[out]multiTimeSeries_SRC_bmulti-detector SRC-frame timeseries, multiplied by AM function b(t).
[in]inputFstatInput structure.

Definition at line 1270 of file ComputeFstat.c.

◆ XLALFstatInputTimeslice()

int XLALFstatInputTimeslice ( FstatInput **  slice,
const FstatInput *  input,
const LIGOTimeGPS minStartGPS,
const LIGOTimeGPS maxStartGPS 
)

Create and return an FstatInput 'timeslice' for given input FstatInput object [must be using LALDemod Fstat method!] which covers only the time interval ['minStartGPS','maxStartGPS'), using conventions of XLALCWGPSinRange().

The returned FstatInput structure is equivalent to one that would have been produced for the given time-interval in the first place, except that it uses "containers" and references the data in the original FstatInput object. A special flag is set in the FstatInput object to notify the destructor to only free the data-containers.

Note: if a timeslice is empty for any detector, the corresponding timeseries 'container' will be allocated and will have length==0 and data==NULL.

Parameters
[out]sliceAddress of a pointer to a FstatInput structure
[in]inputInput data structure
[in]minStartGPSMinimum starting GPS time
[in]maxStartGPSMaximum starting GPS time

Definition at line 1308 of file ComputeFstat.c.

◆ XLALDestroyFstatInputTimeslice_common()

static void XLALDestroyFstatInputTimeslice_common ( FstatCommon common)
static

Definition at line 1439 of file ComputeFstat.c.

Enumeration Type Documentation

◆ FstatQuantities

Bit-field of \( \mathcal{F} \) -statistic quantities which can be computed by XLALComputeFstat().

Not all options are supported by all \( \mathcal{F} \) -statistic methods.

Enumerator
FSTATQ_NONE 

Do not compute \( \mathcal{F} \) -statistic, still compute buffered quantities.

FSTATQ_2F 

Compute multi-detector \( 2\mathcal{F} \) .

FSTATQ_FAFB 

Compute multi-detector \( F_a \) and \( F_b \) .

FSTATQ_2F_PER_DET 

Compute \( 2\mathcal{F} \) for each detector.

FSTATQ_FAFB_PER_DET 

Compute \( F_a \) and \( F_b \) for each detector.

FSTATQ_ATOMS_PER_DET 

Compute per-SFT \( \mathcal{F} \) -statistic atoms for each detector (Demod only).

FSTATQ_2F_CUDA 

Compute multi-detector \( 2\mathcal{F} \) , store results on CUDA GPU (CUDA implementation of Resamp only).

FSTATQ_LAST 

Definition at line 94 of file ComputeFstat.h.

◆ FstatMethodType

Different methods available to compute the \( \mathcal{F} \) -statistic, falling into two broad classes:

  • Demod: Dirichlet kernel-based demodulation [36] .
  • Resamp: FFT-based resampling [9] [22] [7] .
Enumerator
FMETHOD_DEMOD_GENERIC 

Demod: generic C hotloop, works for any number of Dirichlet kernel terms \( \text{Dterms} \)

FMETHOD_DEMOD_OPTC 

Demod: gptimized C hotloop using Akos' algorithm, only works for \( \text{Dterms} \lesssim 20 \)

FMETHOD_DEMOD_ALTIVEC 

Demod: Altivec hotloop variant, uses fixed \( \text{Dterms} = 8 \)

FMETHOD_DEMOD_SSE 

Demod: SSE hotloop with precalc divisors, uses fixed \( \text{Dterms} = 8 \)

FMETHOD_DEMOD_BEST 

Demod: best guess of the fastest available hotloop

FMETHOD_RESAMP_GENERIC 

Resamp: generic implementation [22]

FMETHOD_RESAMP_CUDA 

Resamp: CUDA resampling [7]

FMETHOD_RESAMP_BEST 

Resamp: best guess of the fastest available implementation

Definition at line 110 of file ComputeFstat.h.

Macro Definition Documentation

◆ DEFAULT_MAX_MISMATCH_FROM_SFT_LENGTH

#define DEFAULT_MAX_MISMATCH_FROM_SFT_LENGTH   0.05

default maximum allowed F-stat mismatch from SFTs being too long, to be used in XLALFstatCheckSFTLengthMismatch()

  • this value allows a search using 1800-second SFTs for isolated CW signals with frequencies <~ 2054 Hz, while preventing a search using 1800-second SFTs for a Sco X-1-like binary CW signal (which would require <~ 200-second SFTs)

Definition at line 70 of file ComputeFstat.h.

◆ TIMING_MODEL_MAX_VARS

#define TIMING_MODEL_MAX_VARS   10

Definition at line 318 of file ComputeFstat.h.

Variable Documentation

◆ FstatOptionalArgsDefaults

const FstatOptionalArgs FstatOptionalArgsDefaults
extern

Global initializer for setting FstatOptionalArgs to default values.

Definition at line 90 of file ComputeFstat.c.