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
ComputeFstat_Demod.c
Go to the documentation of this file.
1//
2// Copyright (C) 2012--2015 Karl Wette
3// Copyright (C) 2005--2007, 2009, 2010, 2012, 2014 Reinhard Prix
4// Copyright (C) 2007--2010, 2012 Bernd Machenschalk
5// Copyright (C) 2007 Chris Messenger
6// Copyright (C) 2006 John T. Whelan, Badri Krishnan
7//
8// This program is free software; you can redistribute it and/or modify
9// it under the terms of the GNU General Public License as published by
10// the Free Software Foundation; either version 2 of the License, or
11// (at your option) any later version.
12//
13// This program is distributed in the hope that it will be useful,
14// but WITHOUT ANY WARRANTY; without even the implied warranty of
15// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16// GNU General Public License for more details.
17//
18// You should have received a copy of the GNU General Public License
19// along with with program; see the file COPYING. If not, write to the
20// Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
21// MA 02110-1301 USA
22//
23
24#include <stdlib.h>
25#include <stdio.h>
26#include <math.h>
27
29
30#include <lal/Factorial.h>
31#include <lal/LogPrintf.h>
32#include <lal/SinCosLUT.h>
33
34///
35/// \defgroup ComputeFstat_Demod_c Module ComputeFstat_Demod.c
36/// \ingroup ComputeFstat_h
37/// \brief Implements the \a Demod Dirichlet kernel-based demodulation algorithm
38/// for computing the \f$ \mathcal{F} \f$ -statistic \cite Williams1999 .
39///
40
41/// @{
42
43// ========== Demod internals ==========
44
45// ----- local types ----------
46
47// ---------- BEGIN: Demod-specific timing model data ----------
48typedef struct tagFstatTimingDemod {
49 REAL4 Nsft; // (average) number of SFTs per detector
50 REAL4 tau0_coreLD; // 'fundamental' timing coefficient: fully-buffered time to compute F-stat for one SFT
51 REAL4 tau0_bufferLD; // 'fundamental' timing coefficient: time to re-compute buffered quantities for one SFT
53
54static const char FstatTimingDemodHelp[] =
55 "%%%% ----- Demod-specific F-statistic timing model -----\n"
56 "%%%% Nsft: (average) number of SFTs per detector\n"
57 "%%%% tau0_coreLD: timing coefficient for core Demod F-stat time\n"
58 "%%%% tau0_bufferLD: timing coefficient for computation of buffered quantities\n"
59 "%%%%\n"
60 "%%%% Demod F-statistic timing model:\n"
61 "%%%% tauF_core = Nsft * tau0_coreLD\n"
62 "%%%% tauF_buffer = Nsft * tau0_bufferLD / NFbin\n"
63 "%%%%"
64 "";
65// ---------- END: Demod-specific timing model data ----------
66
67
68typedef struct {
69 int ( *computefafb_func )( // XLALComputeFaFb_...() function for the selected Demod hotloop
70 COMPLEX8 *, COMPLEX8 *, FstatAtomVector **, const SFTVector *, const PulsarSpins, const SSBtimes *, const AMCoeffs *, const UINT4 Dterms
71 );
72 UINT4 Dterms; // Number of terms to keep in Dirichlet kernel
73 MultiSFTVector *multiSFTs; // Input multi-detector SFTs
74 REAL8 prevAlpha, prevDelta; // buffering: previous skyposition computed
75 LIGOTimeGPS prevRefTime; // buffering: keep track of previous refTime for SSBtimes buffering
76 MultiSSBtimes *prevMultiSSBtimes; // buffering: previous multiSSB times, unique to skypos + SFTs
77 MultiAMCoeffs *prevMultiAMcoef; // buffering: previous AM-coeffs, unique to skypos + SFTs
78
79 // ----- timing -----
80 BOOLEAN collectTiming; // flag whether or not to collect timing information
81 FstatTimingGeneric timingGeneric; // measured generic F-statistic timing values
82 FstatTimingDemod timingDemod; // storage for measured Demod-specific timing model
83
85
86// ----- local prototypes ----------
87
88int XLALSetupFstatDemod( void **method_data, FstatCommon *common, FstatMethodFuncs *funcs, MultiSFTVector *multiSFTs, const FstatOptionalArgs *optArgs );
89
90int XLALComputeFaFb_Generic( COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts,
91 const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms );
92
93int XLALComputeFaFb_OptC( COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts,
94 const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms );
95
96#ifdef HAVE_ALTIVEC
97int XLALComputeFaFb_Altivec( COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts,
98 const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms );
99#endif
100
101#ifdef HAVE_SSE_COMPILER
102int XLALComputeFaFb_SSE( COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts,
103 const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms );
104#endif
105
106int XLALGetFstatTiming_Demod( const void *method_data, FstatTimingGeneric *timingGeneric, FstatTimingModel *timingModel );
107void *XLALFstatInputTimeslice_Demod( const void *method_data, const UINT4 iStart[PULSAR_MAX_DETECTORS], const UINT4 iEnd[PULSAR_MAX_DETECTORS] );
108void XLALDestroyFstatInputTimeslice_Demod( void *method_data );
109
110// ----- local function definitions ----------
111static int
113 const FstatCommon *common,
114 void *method_data
115 )
116{
117 // Check input
118 XLAL_CHECK( Fstats != NULL, XLAL_EFAULT );
119 XLAL_CHECK( common != NULL, XLAL_EFAULT );
120 XLAL_CHECK( method_data != NULL, XLAL_EFAULT );
121
122 DemodMethodData *demod = ( DemodMethodData * ) method_data;
123 // get internal timing info
124 BOOLEAN collectTiming = demod->collectTiming;
125 BOOLEAN BufferRecomputed = 0; // keep track if buffer was recomputed in this call or not
126 REAL8 Tau_buffer = 0;
127 REAL8 tic = 0, toc = 0;
128
129 // Get which F-statistic quantities to compute
130 const FstatQuantities whatToCompute = Fstats->whatWasComputed;
131
132 // handy shortcuts
133 BOOLEAN returnAtoms = ( whatToCompute & FSTATQ_ATOMS_PER_DET );
134 PulsarDopplerParams thisPoint = Fstats->doppler;
135 const REAL8 fStart = thisPoint.fkdot[0];
136 const MultiSFTVector *multiSFTs = demod->multiSFTs;
137 const MultiNoiseWeights *multiWeights = common->multiNoiseWeights;
138 const MultiDetectorStateSeries *multiDetStates = common->multiDetectorStates;
139
140 UINT4 numDetectors = multiSFTs->length;
141 XLAL_CHECK( multiDetStates->length == numDetectors, XLAL_EINVAL );
142 XLAL_CHECK( multiWeights == NULL || ( multiWeights->length == numDetectors ), XLAL_EINVAL );
143
144 // initialize timing info struct
145 if ( collectTiming ) {
146 tic = XLALGetCPUTime();
147 }
148
149 MultiSSBtimes *multiSSB = NULL;
150 MultiAMCoeffs *multiAMcoef = NULL;
151 // ----- check if we have buffered SSB+AMcoef for current sky-position
152 if ( ( demod->prevAlpha == thisPoint.Alpha ) && ( demod->prevDelta == thisPoint.Delta ) &&
153 ( demod->prevMultiSSBtimes != NULL ) && ( XLALGPSDiff( &demod->prevRefTime, &thisPoint.refTime ) == 0 ) && // have SSB times for same reftime?
154 ( demod->prevMultiAMcoef != NULL )
155 ) {
156 // if yes ==> reuse
157 multiSSB = demod->prevMultiSSBtimes;
158 multiAMcoef = demod->prevMultiAMcoef;
159 BufferRecomputed = 0;
160 } else {
161 // if not, compute SSB + AMcoef for this skyposition
162 BufferRecomputed = 1;
164
165 SkyPosition skypos;
167 skypos.longitude = thisPoint.Alpha;
168 skypos.latitude = thisPoint.Delta;
169 XLAL_CHECK( ( multiSSB = XLALGetMultiSSBtimes( multiDetStates, skypos, thisPoint.refTime, common->SSBprec ) ) != NULL, XLAL_EFUNC );
170 XLAL_CHECK( ( multiAMcoef = XLALComputeMultiAMCoeffs( multiDetStates, multiWeights, skypos ) ) != NULL, XLAL_EFUNC );
171
172 // store these for possible later re-use in buffer
174 demod->prevMultiSSBtimes = multiSSB;
175 demod->prevRefTime = thisPoint.refTime;
177 demod->prevMultiAMcoef = multiAMcoef;
178 demod->prevAlpha = thisPoint.Alpha;
179 demod->prevDelta = thisPoint.Delta;
180 } // if could not reuse previously buffered quantites
181
182 if ( collectTiming ) {
183 toc = XLALGetCPUTime();
184 Tau_buffer = ( toc - tic );
185 }
186
187 MultiSSBtimes *multiBinary = NULL;
188 MultiSSBtimes *multiSSBTotal = NULL;
189 // handle binary-orbital timing corrections, if applicable
190 if ( thisPoint.asini > 0 ) {
191 // compute binary time corrections to the SSB time delays and SSB time derivitive
192 XLAL_CHECK( XLALAddMultiBinaryTimes( &multiBinary, multiSSB, &thisPoint ) == XLAL_SUCCESS, XLAL_EFUNC );
193 multiSSBTotal = multiBinary;
194 } else {
195 multiSSBTotal = multiSSB;
196 }
197
198 // ----- compute final Fstatistic-value -----
199 REAL4 Ad = multiAMcoef->Mmunu.Ad;
200 REAL4 Bd = multiAMcoef->Mmunu.Bd;
201 REAL4 Cd = multiAMcoef->Mmunu.Cd;
202 REAL4 Ed = multiAMcoef->Mmunu.Ed;;
203 REAL4 Dd_inv = 1.0 / multiAMcoef->Mmunu.Dd;
204
205 // ---------- Compute F-stat for each frequency bin ----------
206 for ( UINT4 k = 0; k < Fstats->numFreqBins; k++ ) {
207 // Set frequency to search at
208 thisPoint.fkdot[0] = fStart + k * Fstats->dFreq;
209
210 COMPLEX8 Fa = 0; // complex amplitude Fa
211 COMPLEX8 Fb = 0; // complex amplitude Fb
212 MultiFstatAtomVector *multiFstatAtoms = NULL; // per-IFO, per-SFT arrays of F-stat 'atoms', ie quantities required to compute F-stat
213
214 // prepare return of 'FstatAtoms' if requested
215 if ( returnAtoms ) {
216 XLAL_CHECK( ( multiFstatAtoms = XLALMalloc( sizeof( *multiFstatAtoms ) ) ) != NULL, XLAL_ENOMEM );
217 multiFstatAtoms->length = numDetectors;
218 XLAL_CHECK( ( multiFstatAtoms->data = XLALMalloc( numDetectors * sizeof( *multiFstatAtoms->data ) ) ) != NULL, XLAL_ENOMEM );
219 } // if returnAtoms
220
221 // loop over detectors and compute all detector-specific quantities
222 for ( UINT4 X = 0; X < numDetectors; X ++ ) {
223 COMPLEX8 FaX, FbX;
224 FstatAtomVector *FstatAtoms = NULL;
225 FstatAtomVector **FstatAtoms_p = returnAtoms ? ( &FstatAtoms ) : NULL;
226
227 // call XLALComputeFaFb_...() function for the user-requested hotloop variant
228 XLAL_CHECK( ( demod->computefafb_func )( &FaX, &FbX, FstatAtoms_p, multiSFTs->data[X], thisPoint.fkdot,
229 multiSSBTotal->data[X], multiAMcoef->data[X], demod->Dterms ) == XLAL_SUCCESS, XLAL_EFUNC );
230
231 if ( returnAtoms ) {
232 multiFstatAtoms->data[X] = FstatAtoms; // copy pointer to IFO-specific Fstat-atoms 'contents'
233 }
234
235 XLAL_CHECK( isfinite( creal( FaX ) ) && isfinite( cimag( FaX ) ) && isfinite( creal( FbX ) ) && isfinite( cimag( FbX ) ), XLAL_EFPOVRFLW );
236
237 if ( whatToCompute & FSTATQ_FAFB_PER_DET ) {
238 Fstats->FaPerDet[X][k] = FaX;
239 Fstats->FbPerDet[X][k] = FbX;
240 }
241
242 // compute single-IFO F-stats, if requested
243 if ( whatToCompute & FSTATQ_2F_PER_DET ) {
244 REAL4 AdX = multiAMcoef->data[X]->A;
245 REAL4 BdX = multiAMcoef->data[X]->B;
246 REAL4 CdX = multiAMcoef->data[X]->C;
247 REAL4 EdX = 0;
248 REAL4 DdX_inv = 1.0 / multiAMcoef->data[X]->D;
249
250 // compute final single-IFO F-stat
251 Fstats->twoFPerDet[X][k] = compute_fstat_from_fa_fb( FaX, FbX, AdX, BdX, CdX, EdX, DdX_inv );
252
253 } // if FSTATQ_2F_PER_DET
254
255 /* Fa = sum_X Fa_X */
256 Fa += FaX;
257
258 /* Fb = sum_X Fb_X */
259 Fb += FbX;
260
261 } // for X < numDetectors
262
263 if ( whatToCompute & FSTATQ_2F ) {
264 Fstats->twoF[k] = compute_fstat_from_fa_fb( Fa, Fb, Ad, Bd, Cd, Ed, Dd_inv );
265 }
266 if ( Fstats->whatWasComputed & FSTATQ_2F_CUDA ) {
267 XLAL_ERROR( XLAL_EINVAL, "Not implemented for FSTATQ_2F_CUDA" );
268 }
269 if ( Fstats->whatWasComputed & FSTATQ_FAFB_CUDA ) {
270 XLAL_ERROR( XLAL_EINVAL, "Not implemented for FSTATQ_FAFB_CUDA" );
271 }
272
273 // Return multi-detector Fa & Fb
274 if ( whatToCompute & FSTATQ_FAFB ) {
275 Fstats->Fa[k] = Fa;
276 Fstats->Fb[k] = Fb;
277 }
278
279 // Return F-atoms per detector
280 if ( whatToCompute & FSTATQ_ATOMS_PER_DET ) {
282 Fstats->multiFatoms[k] = multiFstatAtoms;
283 }
284
285 } // for k < Fstats->numFreqBins
286
287 // this needs to be free'ed, as it's currently not buffered
288 XLALDestroyMultiSSBtimes( multiBinary );
289
290 // Return amplitude modulation coefficients
291 Fstats->Mmunu = demod->prevMultiAMcoef->Mmunu;
292
293 // return per-detector antenna-pattern matrices
294 for ( UINT4 X = 0; X < numDetectors; X ++ ) {
295 Fstats->MmunuX[X].Ad = multiAMcoef->data[X]->A;
296 Fstats->MmunuX[X].Bd = multiAMcoef->data[X]->B;
297 Fstats->MmunuX[X].Cd = multiAMcoef->data[X]->C;
298 Fstats->MmunuX[X].Dd = multiAMcoef->data[X]->D;
299 Fstats->MmunuX[X].Ed = 0;
300 }
301
302 if ( collectTiming ) {
303 toc = XLALGetCPUTime();
304 REAL8 Tau_total = ( toc - tic );
305
306 FstatTimingGeneric *tiGen = &( demod->timingGeneric ); // generic part
307 FstatTimingDemod *tiLD = &( demod->timingDemod ); // Demod-method specific part
308
309 XLAL_CHECK( numDetectors == tiGen->Ndet, XLAL_EINVAL, "Inconsistent number of detectors between XLALCreateSetup() [%d] and XLALComputeFstat() [%d]\n", tiGen->Ndet, numDetectors );
310
311 // rescale all relevant timings to per-detector
312 Tau_total /= numDetectors;
313 Tau_buffer /= numDetectors;
314
315 // compute generic F-stat timing model contributions
316 UINT4 NFbin = Fstats->numFreqBins;
317 REAL8 tauF_eff = Tau_total / NFbin;
318 REAL8 tauF_core = ( Tau_total - Tau_buffer ) / NFbin;
319
320 // compute Demod timing model coefficients
321 REAL8 NsftPerDet = tiLD->Nsft;
322 REAL8 tau0_coreLD = tauF_core / NsftPerDet;
323
324 // update the averaged timing-model quantities
325 tiGen->NCalls ++; // keep track of number of Fstat-calls for timing
326#define updateAvgF(q) tiGen->q = ((tiGen->q *(tiGen->NCalls-1) + q)/(tiGen->NCalls))
327 updateAvgF( tauF_eff );
328 updateAvgF( tauF_core );
329 // we also average NFbin, which can be different betnween different calls to XLALComputeFstat() (contrary to Ndet)
330 updateAvgF( NFbin );
331
332#define updateAvgLD(q) tiLD->q = ((tiLD->q *(tiGen->NCalls-1) + q)/(tiGen->NCalls))
333 updateAvgLD( tau0_coreLD );
334
335 // buffer-quantities only updated if buffer was actually recomputed
336 if ( BufferRecomputed ) {
337 REAL8 tau0_bufferLD = Tau_buffer / NsftPerDet;
338 REAL8 tauF_buffer = Tau_buffer / NFbin;
339
340 updateAvgF( tauF_buffer );
341 updateAvgLD( tau0_bufferLD );
342 } // if BufferRecomputed
343
344 } // if collect timing
345
346 return XLAL_SUCCESS;
347
348} // XLALComputeFstatDemod()
349
350
351static void
352XLALDestroyDemodMethodData( void *method_data )
353{
354
355 DemodMethodData *demod = ( DemodMethodData * ) method_data;
356
360 XLALFree( demod );
361
362} // XLALDestroyDemodMethodData()
363
364int
365XLALSetupFstatDemod( void **method_data,
366 FstatCommon *common,
367 FstatMethodFuncs *funcs,
369 const FstatOptionalArgs *optArgs
370 )
371{
372 // Check input
373 XLAL_CHECK( method_data != NULL, XLAL_EFAULT );
374 XLAL_CHECK( common != NULL, XLAL_EFAULT );
375 XLAL_CHECK( funcs != NULL, XLAL_EFAULT );
376 XLAL_CHECK( multiSFTs != NULL, XLAL_EFAULT );
377 XLAL_CHECK( optArgs != NULL, XLAL_EFAULT );
378
379 // Allocate method data
380 DemodMethodData *demod = *method_data = XLALCalloc( 1, sizeof( *demod ) );
381 XLAL_CHECK( demod != NULL, XLAL_ENOMEM );
382
383 // Set method function pointers
386 funcs->workspace_destroy_func = NULL;
387
388 // Save pointer to SFTs
389 demod->multiSFTs = multiSFTs;
390
391 // Save Dterms
392 demod->Dterms = optArgs->Dterms;
393
394 // turn on timing collection if requested
395 demod->collectTiming = optArgs->collectTiming;
396
397 // initialize struct for collecting timing data, store invariant 'meta' quantities about this setup
398 if ( demod->collectTiming ) {
400 XLAL_INIT_MEM( demod->timingDemod );
402 UINT4 numSFTs = 0;
403 for ( UINT4 X = 0; X < numDetectors; X ++ ) {
404 numSFTs += demod->multiSFTs->data[X]->length;
405 }
407 demod->timingDemod.Nsft = 1.0 * numSFTs / numDetectors; // average number of sfts *per detector*
408 } // if collectTiming
409
410 // Select XLALComputeFaFb_...() function for the user-requested hotloop variant
411 switch ( optArgs->FstatMethod ) {
414 break;
417 break;
418#ifdef HAVE_ALTIVEC
420 demod->computefafb_func = XLALComputeFaFb_Altivec;
421 break;
422#endif
423#ifdef HAVE_SSE_COMPILER
425 demod->computefafb_func = XLALComputeFaFb_SSE;
426 break;
427#endif
428 default:
429 XLAL_ERROR( XLAL_EINVAL, "Invalid Demod hotloop optArgs->FstatMethod='%d'", optArgs->FstatMethod );
430 break;
431 }
432
433 return XLAL_SUCCESS;
434
435} // XLALSetupFstatDemod()
436
437
438int
439XLALGetFstatTiming_Demod( const void *method_data, FstatTimingGeneric *timingGeneric, FstatTimingModel *timingModel )
440{
441 XLAL_CHECK( method_data != NULL, XLAL_EINVAL );
442 XLAL_CHECK( timingGeneric != NULL, XLAL_EINVAL );
443 XLAL_CHECK( timingModel != NULL, XLAL_EINVAL );
444
445 const DemodMethodData *demod = ( const DemodMethodData * ) method_data;
446 XLAL_CHECK( demod != NULL, XLAL_EINVAL );
447
448 ( *timingGeneric ) = demod->timingGeneric; // struct-copy generic timing measurements
449
450 const FstatTimingDemod *tiLD = &( demod->timingDemod );
451
452 // return method-specific timing model values
453 XLAL_INIT_MEM( ( *timingModel ) );
454
455 UINT4 i = 0;
456 timingModel->names[i] = "Nsft";
457 timingModel->values[i] = tiLD->Nsft;
458
459 i++;
460 timingModel->names[i] = "tau0_coreLD";
461 timingModel->values[i] = tiLD->tau0_coreLD;
462
463 i++;
464 timingModel->names[i] = "tau0_bufferLD";
465 timingModel->values[i] = tiLD->tau0_bufferLD;
466
467 timingModel->numVariables = i + 1;
468 timingModel->help = FstatTimingDemodHelp;
469
470 return XLAL_SUCCESS;
471} // XLALGetFstatTiming_Demod()
472
473// Generates an FstatInput Timeslice based on minStartGPS and maxStartGPS according to XLALCWGPSinRange().
474void *
475XLALFstatInputTimeslice_Demod( const void *method_data,
476 const UINT4 iStart[PULSAR_MAX_DETECTORS],
477 const UINT4 iEnd[PULSAR_MAX_DETECTORS]
478 )
479{
480 XLAL_CHECK_NULL( method_data != NULL, XLAL_EINVAL );
481 XLAL_CHECK_NULL( iStart != NULL, XLAL_EINVAL );
482 XLAL_CHECK_NULL( iEnd != NULL, XLAL_EINVAL );
483
484 const DemodMethodData *demod_input = ( const DemodMethodData * )method_data;
485
486 // allocate memory and copy the input method_data struct
487 DemodMethodData *demod_slice;
488 XLAL_CHECK_NULL( ( demod_slice = XLALCalloc( 1, sizeof( *demod_input ) ) ) != NULL, XLAL_ENOMEM );
489 memcpy( demod_slice, demod_input, sizeof( *demod_input ) );
490
491 UINT4 numIFOs = demod_input->multiSFTs->length;
492
493 // prepare MultiSFTs struct
495 XLAL_CHECK_NULL( ( multiSFTs = XLALCalloc( 1, sizeof( *multiSFTs ) ) ) != NULL, XLAL_ENOMEM );
496 XLAL_CHECK_NULL( ( multiSFTs->data = XLALCalloc( numIFOs, sizeof( *multiSFTs->data ) ) ) != NULL, XLAL_ENOMEM );
497 multiSFTs->length = numIFOs;
498
499 // loop over all detectors
500 for ( UINT4 X = 0; X < numIFOs; X ++ ) {
501 XLAL_CHECK_NULL( ( multiSFTs->data[X] = XLALCalloc( 1, sizeof( *multiSFTs->data[X] ) ) ) != NULL, XLAL_EFUNC );
502
503 if ( iStart[X] > iEnd[X] ) {
504 continue; // empty slice
505 }
506 UINT4 sliceLength = iEnd[X] - iStart[X] + 1; // guaranteed >= 1
507 multiSFTs->data[X]->length = sliceLength;
508 multiSFTs->data[X]->data = &( demod_input->multiSFTs->data[X]->data[iStart[X]] );
509
510 } // for loop over all detectors
511
512 demod_slice->multiSFTs = multiSFTs;
513
514 // empty all buffering quantities
515 demod_slice->prevAlpha = 0;
516 demod_slice->prevDelta = 0;
517 XLAL_INIT_MEM( demod_slice->prevRefTime );
518 demod_slice->prevMultiSSBtimes = NULL;
519 demod_slice->prevMultiAMcoef = NULL;
520
521 // reset timing counters
522 XLAL_INIT_MEM( demod_slice->timingGeneric );
523 XLAL_INIT_MEM( demod_slice->timingDemod );
524
525 return demod_slice;
526
527} // XLALFstatInputTimeslice_Demod()
528
529///
530/// Free all memory not needed by the orginal FstatInput structure
531///
532void
534{
535 if ( !method_data ) {
536 return;
537 }
538
539 DemodMethodData *demod = ( DemodMethodData * ) method_data;
540
543
544 for ( UINT4 X = 0; X < demod->multiSFTs->length; X ++ ) {
545 XLALFree( demod->multiSFTs->data[X] );
546 }
547
548 XLALFree( demod->multiSFTs->data );
549 XLALFree( demod->multiSFTs );
550 XLALFree( demod );
551
552 return;
553
554} // XLALDestroyFstatInputTimeslice_Demod()
555
556/// @}
#define updateAvgF(q)
#define updateAvgLD(q)
static REAL4 compute_fstat_from_fa_fb(COMPLEX8 Fa, COMPLEX8 Fb, REAL4 A, REAL4 B, REAL4 C, REAL4 E, REAL4 Dinv)
int k
int XLALComputeFaFb_OptC(COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts, const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms)
static void XLALDestroyDemodMethodData(void *method_data)
void * XLALFstatInputTimeslice_Demod(const void *method_data, const UINT4 iStart[PULSAR_MAX_DETECTORS], const UINT4 iEnd[PULSAR_MAX_DETECTORS])
int XLALComputeFaFb_Generic(COMPLEX8 *Fa, COMPLEX8 *Fb, FstatAtomVector **FstatAtoms, const SFTVector *sfts, const PulsarSpins fkdot, const SSBtimes *tSSB, const AMCoeffs *amcoe, const UINT4 Dterms)
static const char FstatTimingDemodHelp[]
int XLALSetupFstatDemod(void **method_data, FstatCommon *common, FstatMethodFuncs *funcs, MultiSFTVector *multiSFTs, const FstatOptionalArgs *optArgs)
int XLALGetFstatTiming_Demod(const void *method_data, FstatTimingGeneric *timingGeneric, FstatTimingModel *timingModel)
static int XLALComputeFstatDemod(FstatResults *Fstats, const FstatCommon *common, void *method_data)
void XLALDestroyFstatInputTimeslice_Demod(void *method_data)
Free all memory not needed by the orginal FstatInput structure.
FstatQuantities
Bit-field of -statistic quantities which can be computed by XLALComputeFstat().
Definition: ComputeFstat.h:94
void XLALDestroyMultiFstatAtomVector(MultiFstatAtomVector *multiAtoms)
Free all memory associated with a MultiFstatAtomVector structure.
Definition: ComputeFstat.c:338
@ FMETHOD_DEMOD_SSE
Demod: SSE hotloop with precalc divisors, uses fixed
Definition: ComputeFstat.h:120
@ FMETHOD_DEMOD_OPTC
Demod: gptimized C hotloop using Akos' algorithm, only works for
Definition: ComputeFstat.h:118
@ FMETHOD_DEMOD_GENERIC
Demod: generic C hotloop, works for any number of Dirichlet kernel terms
Definition: ComputeFstat.h:117
@ FMETHOD_DEMOD_ALTIVEC
Demod: Altivec hotloop variant, uses fixed
Definition: ComputeFstat.h:119
@ FSTATQ_2F
Compute multi-detector .
Definition: ComputeFstat.h:96
@ FSTATQ_FAFB_PER_DET
Compute and for each detector.
Definition: ComputeFstat.h:99
@ FSTATQ_FAFB
Compute multi-detector and .
Definition: ComputeFstat.h:97
@ FSTATQ_2F_CUDA
Compute multi-detector , store results on CUDA GPU (CUDA implementation of Resamp only).
Definition: ComputeFstat.h:101
@ FSTATQ_ATOMS_PER_DET
Compute per-SFT -statistic atoms for each detector (Demod only).
Definition: ComputeFstat.h:100
@ FSTATQ_2F_PER_DET
Compute for each detector.
Definition: ComputeFstat.h:98
@ FSTATQ_FAFB_CUDA
Compute multi-detector and , store results on CUDA GPU (CUDA implementation of Resamp only).
Definition: ComputeFstat.h:102
void XLALDestroyMultiAMCoeffs(MultiAMCoeffs *multiAMcoef)
Destroy a MultiAMCoeffs structure.
Definition: LALComputeAM.c:469
MultiAMCoeffs * XLALComputeMultiAMCoeffs(const MultiDetectorStateSeries *multiDetStates, const MultiNoiseWeights *multiWeights, SkyPosition skypos)
Multi-IFO version of XLALComputeAMCoeffs().
Definition: LALComputeAM.c:379
unsigned char BOOLEAN
#define XLAL_INIT_MEM(x)
double REAL8
uint32_t UINT4
float complex COMPLEX8
float REAL4
void * XLALMalloc(size_t n)
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
REAL8 XLALGetCPUTime(void)
#define PULSAR_MAX_DETECTORS
maximal number of detectors we can handle (for static arrays of detector quantities)
REAL8 PulsarSpins[PULSAR_MAX_SPINS]
Typedef for fixed-size array holding GW frequency and derivatives fk = d^k Freq/dt^k|(tau_ref)
void XLALDestroyMultiSFTVector(MultiSFTVector *multvect)
Destroy a multi SFT-vector.
Definition: SFTtypes.c:424
MultiSSBtimes * XLALGetMultiSSBtimes(const MultiDetectorStateSeries *multiDetStates, SkyPosition skypos, LIGOTimeGPS refTime, SSBprecision precision)
Multi-IFO version of XLALGetSSBtimes().
Definition: SSBtimes.c:657
void XLALDestroyMultiSSBtimes(MultiSSBtimes *multiSSB)
Destroy a MultiSSBtimes structure.
Definition: SSBtimes.c:845
int XLALAddMultiBinaryTimes(MultiSSBtimes **multiSSBOut, const MultiSSBtimes *multiSSBIn, const PulsarDopplerParams *Doppler)
Multi-IFO version of XLALAddBinaryTimes().
Definition: SSBtimes.c:413
COORDINATESYSTEM_EQUATORIAL
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
#define XLAL_CHECK_NULL(assertion,...)
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EINVAL
XLAL_EFPOVRFLW
REAL8 XLALGPSDiff(const LIGOTimeGPS *t1, const LIGOTimeGPS *t0)
fStart
size_t numDetectors
This structure contains the per-SFT (weighted) antenna-pattern functions , with the SFT-index,...
Definition: LALComputeAM.h:63
REAL4 B
summed antenna-pattern matrix coefficient:
Definition: LALComputeAM.h:67
REAL4 A
summed antenna-pattern matrix coefficient:
Definition: LALComputeAM.h:66
REAL4 C
summed antenna-pattern matrix coefficient:
Definition: LALComputeAM.h:68
REAL4 D
determinant
Definition: LALComputeAM.h:69
REAL4 Dd
determinant factor , such that
Definition: LALComputeAM.h:132
A vector of COMPLEX8FrequencySeries.
COMPLEX8FrequencySeries * data
Pointer to the data array.
UINT4 length
Number of elements in array.
MultiSSBtimes * prevMultiSSBtimes
LIGOTimeGPS prevRefTime
FstatTimingGeneric timingGeneric
FstatTimingDemod timingDemod
MultiSFTVector * multiSFTs
int(* computefafb_func)(COMPLEX8 *, COMPLEX8 *, FstatAtomVector **, const SFTVector *, const PulsarSpins, const SSBtimes *, const AMCoeffs *, const UINT4 Dterms)
MultiAMCoeffs * prevMultiAMcoef
A vector of -statistic 'atoms', i.e.
Definition: ComputeFstat.h:175
MultiDetectorStateSeries * multiDetectorStates
SSBprecision SSBprec
MultiNoiseWeights * multiNoiseWeights
int(* compute_func)(FstatResults *, const FstatCommon *, void *)
void(* workspace_destroy_func)(void *)
void(* method_data_destroy_func)(void *)
Struct of optional 'advanced level' and (potentially method-specific) arguments to be passed to the ...
Definition: ComputeFstat.h:137
BOOLEAN collectTiming
a flag to turn on/off the collection of F-stat-method-specific timing-data
Definition: ComputeFstat.h:147
UINT4 Dterms
Number of Dirichlet kernel terms, used by some Demod methods; see FstatMethodType.
Definition: ComputeFstat.h:140
FstatMethodType FstatMethod
Method to use for computing the -statistic.
Definition: ComputeFstat.h:142
XLALComputeFstat() computed results structure.
Definition: ComputeFstat.h:202
AntennaPatternMatrix MmunuX[PULSAR_MAX_DETECTORS]
Per detector antenna pattern matrix , used in computing .
Definition: ComputeFstat.h:231
COMPLEX8 * Fb
Definition: ComputeFstat.h:260
AntennaPatternMatrix Mmunu
Antenna pattern matrix , used in computing .
Definition: ComputeFstat.h:228
REAL4 * twoFPerDet[PULSAR_MAX_DETECTORS]
If whatWasComputed & FSTATQ_2F_PER_DET is true, the values computed at numFreqBins frequencies space...
Definition: ComputeFstat.h:277
FstatQuantities whatWasComputed
Bit-field of which -statistic quantities were computed.
Definition: ComputeFstat.h:234
REAL8 dFreq
Spacing in frequency between each computed -statistic.
Definition: ComputeFstat.h:214
COMPLEX8 * FaPerDet[PULSAR_MAX_DETECTORS]
If whatWasComputed & FSTATQ_PARTS_PER_DET is true, the and values computed at numFreqBins frequenci...
Definition: ComputeFstat.h:287
UINT4 numFreqBins
Number of frequencies at which the were computed.
Definition: ComputeFstat.h:217
REAL4 * twoF
If whatWasComputed & FSTATQ_2F is true, the multi-detector values computed at numFreqBins frequencie...
Definition: ComputeFstat.h:242
COMPLEX8 * Fa
If whatWasComputed & FSTATQ_PARTS is true, the multi-detector and computed at numFreqBins frequenci...
Definition: ComputeFstat.h:259
COMPLEX8 * FbPerDet[PULSAR_MAX_DETECTORS]
Definition: ComputeFstat.h:288
PulsarDopplerParams doppler
Doppler parameters, including the starting frequency, at which the were computed.
Definition: ComputeFstat.h:206
MultiFstatAtomVector ** multiFatoms
If whatWasComputed & FSTATQ_ATOMS_PER_DET is true, the per-SFT -statistic multi-atoms computed at nu...
Definition: ComputeFstat.h:297
Generic F-stat timing coefficients (times in seconds) [see https://dcc.ligo.org/LIGO-T1600531-v4 for ...
Definition: ComputeFstat.h:316
Struct to carry the -statistic method-specific timing model in terms of a list of variable names and...
Definition: ComputeFstat.h:335
REAL4 values[TIMING_MODEL_MAX_VARS]
Definition: ComputeFstat.h:338
const char * help
Definition: ComputeFstat.h:339
const char * names[TIMING_MODEL_MAX_VARS]
Definition: ComputeFstat.h:337
Multi-IFO container for antenna-pattern coefficients and atenna-pattern matrix .
Definition: LALComputeAM.h:137
AMCoeffs ** data
noise-weighted AM-coeffs , and
Definition: LALComputeAM.h:142
AntennaPatternMatrix Mmunu
antenna-pattern matrix
Definition: LALComputeAM.h:143
Multi-IFO time-series of DetectorStates.
UINT4 length
number of detectors
A multi-detector vector of FstatAtomVector.
Definition: ComputeFstat.h:187
FstatAtomVector ** data
Array of FstatAtomVector pointers, one for each detector X.
Definition: ComputeFstat.h:192
UINT4 length
Number of detectors.
Definition: ComputeFstat.h:191
One noise-weight (number) per SFT (therefore indexed over IFOs and SFTs.
Definition: PSDutils.h:71
UINT4 length
number of detectors
Definition: PSDutils.h:75
A collection of SFT vectors – one for each IFO in a multi-IFO search.
Definition: SFTfileIO.h:179
UINT4 length
number of ifos
Definition: SFTfileIO.h:183
SFTVector ** data
sftvector for each ifo
Definition: SFTfileIO.h:184
Multi-IFO container for SSB timings.
Definition: SSBtimes.h:67
SSBtimes ** data
array of SSBtimes (pointers)
Definition: SSBtimes.h:72
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
PulsarSpins fkdot
Intrinsic spins: [Freq, f1dot, f2dot, ... ] where fkdot = d^kFreq/dt^k.
REAL8 Delta
Sky position: DEC (latitude) in equatorial coords and radians.
LIGOTimeGPS refTime
Reference time of pulsar parameters (in SSB!)
REAL8 Alpha
Sky position: RA (longitude) in equatorial coords and radians.
REAL8 asini
Binary: projected, normalized orbital semi-major axis (s).
Simple container for two REAL8-vectors, namely the SSB-timings DeltaT_alpha and Tdot_alpha,...
Definition: SSBtimes.h:60
REAL8 longitude
REAL8 latitude
CoordinateSystem system