LALPulsar 7.1.2.1-bf6a62b
ValidateHoughMulti.c
Go to the documentation of this file.
1/*
2 * Copyright (C) 2005 Badri Krishnan, Alicia Sintes
3 *
4 * This program is free software; you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation; either version 2 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with with program; see the file COPYING. If not, write to the
16 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17 * MA 02110-1301 USA
18 */
19
20/**
21 * \file
22 * \ingroup lalpulsar_bin_Hough
23 * \author Badri Krishnan, Alicia Sintes
24 * \brief Driver code for performing Hough transform search on non-demodulated
25 * data using SFTs from possible multiple IFOs
26 *
27 * History: Created by Sintes and Krishnan July 15, 2006
28 */
29
30
31#include "DriveHoughColor.h"
32#include "MCInjectHoughMulti.h"
33
34/* globals, constants and defaults */
35
36/* boolean global variables for controlling output */
38
39/* #define EARTHEPHEMERIS "./earth05-09.dat" */
40/* #define SUNEPHEMERIS "./sun05-09.dat" */
41
42#define EARTHEPHEMERIS "/home/badkri/lscsoft/share/lal/earth05-09.dat"
43#define SUNEPHEMERIS "/home/badkri/lscsoft/share/lal/sun05-09.dat"
44
45#define MAXFILENAMELENGTH 512 /* maximum # of characters of a filename */
46
47#define DIROUT "./outMulti" /* output directory */
48#define BASENAMEOUT "HM" /* prefix file output */
49
50#define THRESHOLD 1.6 /* thresold for peak selection, with respect to the
51 the averaged power in the search band */
52#define FALSEALARM 1.0e-9 /* Hough false alarm for candidate selection */
53#define SKYFILE "./skypatchfile"
54#define F0 310.0 /* frequency to build the LUT and start search */
55#define FBAND 0.05 /* search frequency band */
56#define NFSIZE 21 /* n-freq. span of the cylinder, to account for spin-down search */
57#define BLOCKSRNGMED 101 /* Running median window size */
59#define TRUE (1==1)
60#define FALSE (1==0)
61
62/* local function prototype */
63
64
65/******************************************/
67int main( int argc, char *argv[] )
68{
69
70 /* LALStatus pointer */
71 static LALStatus status;
72
73 /* time and velocity */
74 static LIGOTimeGPSVector timeV;
75 static REAL8Cart3CoorVector velV;
76 static REAL8Vector timeDiffV;
77 LIGOTimeGPS firstTimeStamp;
78
79 /* standard pulsar sft types */
80 MultiSFTVector *inputSFTs = NULL;
81 UINT4 binsSFT;
82 UINT4 sftFminBin;
83 UINT4 numsft;
84
85 INT4 k;
86
87 /* information about all the ifos */
88 MultiDetectorStateSeries *mdetStates = NULL;
89 UINT4 numifo;
90
91 /* vector of weights */
92 REAL8Vector weightsV;
93
94 /* ephemeris */
95 EphemerisData *edat = NULL;
96
97 static UCHARPeakGram pg1;
98 static HoughTemplate pulsarTemplate;
99 static REAL8Vector foft;
100
101 /* miscellaneous */
102 UINT4 mObsCoh;
103 REAL8 timeBase, deltaF;
104 REAL8 numberCount;
105
106 /* sft constraint variables */
107 LIGOTimeGPS startTimeGPS, endTimeGPS;
108 LIGOTimeGPSVector *inputTimeStampsVector = NULL;
109
110
111 REAL8 alphaPeak, meanN, sigmaN;
112
113 /* user input variables */
114 BOOLEAN uvar_weighAM, uvar_weighNoise;
115 INT4 uvar_blocksRngMed, uvar_nfSizeCylinder, uvar_maxBinsClean;
117 REAL8 uvar_fStart, uvar_peakThreshold, uvar_fSearchBand;
118 REAL8 uvar_Alpha, uvar_Delta, uvar_Freq, uvar_fdot;
119 REAL8 uvar_AlphaWeight, uvar_DeltaWeight;
120 CHAR *uvar_earthEphemeris = NULL;
121 CHAR *uvar_sunEphemeris = NULL;
122 CHAR *uvar_sftDir = NULL;
123 CHAR *uvar_timeStampsFile = NULL;
124 CHAR *uvar_outfile = NULL;
125 LALStringVector *uvar_linefiles = NULL;
126
127
128 /* Set up the default parameters */
129
130 /* LAL error-handler */
132
133 uvar_weighAM = TRUE;
134 uvar_weighNoise = TRUE;
136 uvar_nfSizeCylinder = NFSIZE;
137 uvar_fStart = F0;
138 uvar_fSearchBand = FBAND;
139 uvar_peakThreshold = THRESHOLD;
140 uvar_maxBinsClean = 100;
141 uvar_startTime = 0;
143
144 uvar_Alpha = 1.0;
145 uvar_Delta = 1.0;
146 uvar_Freq = 310.0;
147 uvar_fdot = 0.0;
148
149 uvar_AlphaWeight = uvar_Alpha;
150 uvar_DeltaWeight = uvar_Delta;
151
152 uvar_outfile = ( CHAR * )LALCalloc( MAXFILENAMELENGTH, sizeof( CHAR ) );
153 strcpy( uvar_outfile, "./tempout" );
154
155 uvar_earthEphemeris = ( CHAR * )LALCalloc( MAXFILENAMELENGTH, sizeof( CHAR ) );
156 strcpy( uvar_earthEphemeris, EARTHEPHEMERIS );
157
158 uvar_sunEphemeris = ( CHAR * )LALCalloc( MAXFILENAMELENGTH, sizeof( CHAR ) );
159 strcpy( uvar_sunEphemeris, SUNEPHEMERIS );
160
161
162 /* register user input variables */
163 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_fStart, "fStart", REAL8, 'f', OPTIONAL, "Start search frequency" ) == XLAL_SUCCESS, XLAL_EFUNC );
164 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_fSearchBand, "fSearchBand", REAL8, 'b', OPTIONAL, "Search frequency band" ) == XLAL_SUCCESS, XLAL_EFUNC );
165 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_startTime, "startTime", REAL8, 0, OPTIONAL, "GPS start time of observation (SFT timestamps must be >= this)" ) == XLAL_SUCCESS, XLAL_EFUNC );
166 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_endTime, "endTime", REAL8, 0, OPTIONAL, "GPS end time of observation (SFT timestamps must be < this)" ) == XLAL_SUCCESS, XLAL_EFUNC );
167 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_timeStampsFile, "timeStampsFile", STRING, 0, OPTIONAL, "Input time-stamps file" ) == XLAL_SUCCESS, XLAL_EFUNC );
168 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_peakThreshold, "peakThreshold", REAL8, 0, OPTIONAL, "Peak selection threshold" ) == XLAL_SUCCESS, XLAL_EFUNC );
169 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_weighAM, "weighAM", BOOLEAN, 0, OPTIONAL, "Use amplitude modulation weights" ) == XLAL_SUCCESS, XLAL_EFUNC );
170 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_weighNoise, "weighNoise", BOOLEAN, 0, OPTIONAL, "Use SFT noise weights" ) == XLAL_SUCCESS, XLAL_EFUNC );
171 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_earthEphemeris, "earthEphemeris", STRING, 'E', OPTIONAL, "Earth Ephemeris file" ) == XLAL_SUCCESS, XLAL_EFUNC );
172 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_sunEphemeris, "sunEphemeris", STRING, 'S', OPTIONAL, "Sun Ephemeris file" ) == XLAL_SUCCESS, XLAL_EFUNC );
173 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_sftDir, "sftDir", STRING, 'D', REQUIRED, "SFT filename pattern. Possibilities are:\n"
174 " - '<SFT file>;<SFT file>;...', where <SFT file> may contain wildcards\n - 'list:<file containing list of SFT files>'" ) == XLAL_SUCCESS, XLAL_EFUNC );
175 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_linefiles, "linefiles", STRINGVector, 0, OPTIONAL, "Comma separated List of linefiles (filenames must contain IFO name)" ) == XLAL_SUCCESS, XLAL_EFUNC );
176
177 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_Alpha, "Alpha", REAL8, 0, OPTIONAL, "Sky location (longitude)" ) == XLAL_SUCCESS, XLAL_EFUNC );
178 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_Delta, "Delta", REAL8, 0, OPTIONAL, "Sky location (latitude)" ) == XLAL_SUCCESS, XLAL_EFUNC );
179 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_Freq, "Freq", REAL8, 0, OPTIONAL, "Template frequency" ) == XLAL_SUCCESS, XLAL_EFUNC );
180 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_fdot, "fdot", REAL8, 0, OPTIONAL, "First spindown" ) == XLAL_SUCCESS, XLAL_EFUNC );
181
182 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_AlphaWeight, "AlphaWeight", REAL8, 0, OPTIONAL, "sky Alpha for weight calculation" ) == XLAL_SUCCESS, XLAL_EFUNC );
183 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_DeltaWeight, "DeltaWeight", REAL8, 0, OPTIONAL, "sky Delta for weight calculation" ) == XLAL_SUCCESS, XLAL_EFUNC );
184
185 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_nfSizeCylinder, "nfSizeCylinder", INT4, 0, OPTIONAL, "Size of cylinder of PHMDs" ) == XLAL_SUCCESS, XLAL_EFUNC );
186
187 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_blocksRngMed, "blocksRngMed", INT4, 0, OPTIONAL, "Running Median block size" ) == XLAL_SUCCESS, XLAL_EFUNC );
188 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_maxBinsClean, "maxBinsClean", INT4, 0, OPTIONAL, "Maximum number of bins in cleaning" ) == XLAL_SUCCESS, XLAL_EFUNC );
189 XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_outfile, "outfile", STRING, 0, OPTIONAL, "output file name" ) == XLAL_SUCCESS, XLAL_EFUNC );
190
191
192 /* read all command line variables */
193 BOOLEAN should_exit = 0;
195 if ( should_exit ) {
196 exit( 1 );
197 }
198
199 /* very basic consistency checks on user input */
200 if ( uvar_fStart < 0 ) {
201 fprintf( stderr, "start frequency must be positive\n" );
202 exit( 1 );
203 }
204
205 if ( uvar_fSearchBand < 0 ) {
206 fprintf( stderr, "search frequency band must be positive\n" );
207 exit( 1 );
208 }
209
210 if ( uvar_peakThreshold < 0 ) {
211 fprintf( stderr, "peak selection threshold must be positive\n" );
212 exit( 1 );
213 }
214
215
216 /***** start main calculations *****/
217
218
219
220 /* read sft Files and set up weights */
221 {
222 /* new SFT I/O data types */
223 SFTCatalog *catalog = NULL;
224 static SFTConstraints constraints;
225
226 REAL8 doppWings, f_min, f_max;
227
228 /* set detector constraint */
229 constraints.detector = NULL;
230
232 XLALGPSSetREAL8( &startTimeGPS, uvar_startTime );
233 constraints.minStartTime = &startTimeGPS;
234 }
235
237 XLALGPSSetREAL8( &endTimeGPS, uvar_endTime );
238 constraints.maxStartTime = &endTimeGPS;
239 }
240
241 if ( XLALUserVarWasSet( &uvar_timeStampsFile ) ) {
242 XLAL_CHECK_MAIN( ( inputTimeStampsVector = XLALReadTimestampsFile( uvar_timeStampsFile ) ) != NULL, XLAL_EFUNC );
243 constraints.timestamps = inputTimeStampsVector;
244 }
245
246 /* get sft catalog */
247 XLAL_CHECK_MAIN( ( catalog = XLALSFTdataFind( uvar_sftDir, &constraints ) ) != NULL, XLAL_EFUNC );
248 if ( ( catalog == NULL ) || ( catalog->length == 0 ) ) {
249 fprintf( stderr, "Unable to match any SFTs with pattern '%s'\n", uvar_sftDir );
250 exit( 1 );
251 }
252
253 /* now we can free the inputTimeStampsVector */
254 if ( XLALUserVarWasSet( &uvar_timeStampsFile ) ) {
255 XLALDestroyTimestampVector( inputTimeStampsVector );
256 }
257
258 /* get some sft parameters */
259 mObsCoh = catalog->length; /* number of sfts */
260 deltaF = catalog->data->header.deltaF; /* frequency resolution */
261 timeBase = 1.0 / deltaF; /* coherent integration time */
262 // unused: INT8 f0Bin = floor( uvar_fStart * timeBase + 0.5); /* initial search frequency */
263 // unused: INT4 length = uvar_fSearchBand * timeBase; /* total number of search bins - 1 */
264
265 /* catalog is ordered in time so we can get start, end time and tObs*/
266 firstTimeStamp = catalog->data[0].header.epoch;
267 // unused: LIGOTimeGPS lastTimeStamp = catalog->data[mObsCoh - 1].header.epoch;
268
269 /* allocate memory for velocity vector */
270 velV.length = mObsCoh;
271 velV.data = NULL;
272 velV.data = ( REAL8Cart3Coor * )LALCalloc( mObsCoh, sizeof( REAL8Cart3Coor ) );
273
274 /* allocate memory for timestamps vector */
275 timeV.length = mObsCoh;
276 timeV.data = NULL;
277 timeV.data = ( LIGOTimeGPS * )LALCalloc( mObsCoh, sizeof( LIGOTimeGPS ) );
278
279 /* allocate memory for vector of time differences from start */
280 timeDiffV.length = mObsCoh;
281 timeDiffV.data = NULL;
282 timeDiffV.data = ( REAL8 * )LALCalloc( mObsCoh, sizeof( REAL8 ) );
283
284 /* add wings for Doppler modulation and running median block size*/
285 doppWings = ( uvar_fStart + uvar_fSearchBand ) * VTOT;
286 f_min = uvar_fStart - doppWings - ( uvar_blocksRngMed + uvar_nfSizeCylinder ) * deltaF;
287 f_max = uvar_fStart + uvar_fSearchBand + doppWings + ( uvar_blocksRngMed + uvar_nfSizeCylinder ) * deltaF;
288
289 /* read sft files making sure to add extra bins for running median */
290 /* read the sfts */
291 XLAL_CHECK_MAIN( ( inputSFTs = XLALLoadMultiSFTs( catalog, f_min, f_max ) ) != NULL, XLAL_EFUNC );
292
293
294 /* clean sfts if required */
295 if ( XLALUserVarWasSet( &uvar_linefiles ) ) {
296 RandomParams *randPar = NULL;
297 FILE *fpRand = NULL;
298 INT4 seed, ranCount;
299
300 if ( ( fpRand = fopen( "/dev/urandom", "r" ) ) == NULL ) {
301 fprintf( stderr, "Error in opening /dev/urandom" );
302 exit( 1 );
303 }
304
305 if ( ( ranCount = fread( &seed, sizeof( seed ), 1, fpRand ) ) != 1 ) {
306 fprintf( stderr, "Error in getting random seed" );
307 exit( 1 );
308 }
309
310 LAL_CALL( LALCreateRandomParams( &status, &randPar, seed ), &status );
311
312 LAL_CALL( LALRemoveKnownLinesInMultiSFTVector( &status, inputSFTs, uvar_maxBinsClean, uvar_blocksRngMed, uvar_linefiles, randPar ), &status );
313
314 LAL_CALL( LALDestroyRandomParams( &status, &randPar ), &status );
315 fclose( fpRand );
316 } /* end cleaning */
317
318
319 /* SFT info -- assume all SFTs have same length */
320 numifo = inputSFTs->length;
321 binsSFT = inputSFTs->data[0]->data->data->length;
322 sftFminBin = ( INT4 ) floor( inputSFTs->data[0]->data[0].f0 * timeBase + 0.5 );
323
324
325 XLALDestroySFTCatalog( catalog );
326
327 } /* end of sft reading block */
328
329
330
331 /* get detector velocities weights vector, and timestamps */
332 {
333 MultiNoiseWeights *multweight = NULL;
334 MultiPSDVector *multPSD = NULL;
335 UINT4 iIFO, iSFT, j;
336
337 /* get ephemeris */
338 XLAL_CHECK_MAIN( ( edat = XLALInitBarycenter( uvar_earthEphemeris, uvar_sunEphemeris ) ) != NULL, XLAL_EFUNC );
339
340
341 /* normalize sfts */
342 XLAL_CHECK_MAIN( ( multPSD = XLALNormalizeMultiSFTVect( inputSFTs, uvar_blocksRngMed, NULL ) ) != NULL, XLAL_EFUNC );
343
344 /* set up weights */
345 weightsV.length = mObsCoh;
346 weightsV.data = ( REAL8 * )LALCalloc( 1, mObsCoh * sizeof( REAL8 ) );
347
348 /* initialize all weights to unity */
350
351 /* compute multi noise weights if required */
352 if ( uvar_weighNoise ) {
353 XLAL_CHECK_MAIN( ( multweight = XLALComputeMultiNoiseWeights( multPSD, uvar_blocksRngMed, 0 ) ) != NULL, XLAL_EFUNC );
354 }
355
356 /* we are now done with the psd */
357 XLALDestroyMultiPSDVector( multPSD );
358
359 /* get information about all detectors including velocity and timestamps */
360 /* note that this function returns the velocity at the
361 mid-time of the SFTs -- should not make any difference */
362 const REAL8 tOffset = 0.5 / inputSFTs->data[0]->data[0].deltaF;
363 XLAL_CHECK_MAIN( ( mdetStates = XLALGetMultiDetectorStatesFromMultiSFTs( inputSFTs, edat, tOffset ) ) != NULL, XLAL_EFUNC );
364
365 /* copy the timestamps, weights, and velocity vector */
366 for ( j = 0, iIFO = 0; iIFO < numifo; iIFO++ ) {
367
368 numsft = mdetStates->data[iIFO]->length;
369
370 for ( iSFT = 0; iSFT < numsft; iSFT++, j++ ) {
371
372 velV.data[j].x = mdetStates->data[iIFO]->data[iSFT].vDetector[0];
373 velV.data[j].y = mdetStates->data[iIFO]->data[iSFT].vDetector[1];
374 velV.data[j].z = mdetStates->data[iIFO]->data[iSFT].vDetector[2];
375
376 if ( uvar_weighNoise ) {
377 weightsV.data[j] = multweight->data[iIFO]->data[iSFT];
378 }
379
380 /* mid time of sfts */
381 timeV.data[j] = mdetStates->data[iIFO]->data[iSFT].tGPS;
382
383 } /* loop over SFTs */
384
385 } /* loop over IFOs */
386
387 if ( uvar_weighNoise ) {
389 }
390
391 /* compute the time difference relative to startTime for all SFTs */
392 for ( j = 0; j < mObsCoh; j++ ) {
393 timeDiffV.data[j] = XLALGPSDiff( timeV.data + j, &firstTimeStamp );
394 }
395
396 if ( uvar_weighNoise ) {
397 XLALDestroyMultiNoiseWeights( multweight );
398 }
399
400 } /* end block for noise weights, velocity and time */
401
402
403
404 /* calculate amplitude modulation weights if required */
405 if ( uvar_weighAM ) {
406 MultiAMCoeffs *multiAMcoef = NULL;
407 UINT4 iIFO, iSFT;
408 SkyPosition skypos;
409
410 /* get the amplitude modulation coefficients */
411 skypos.longitude = uvar_AlphaWeight;
412 skypos.latitude = uvar_DeltaWeight;
414 XLAL_CHECK_MAIN( ( multiAMcoef = XLALComputeMultiAMCoeffs( mdetStates, NULL, skypos ) ) != NULL, XLAL_EFUNC );
415
416 /* loop over the weights and multiply them by the appropriate
417 AM coefficients */
418 for ( k = 0, iIFO = 0; iIFO < numifo; iIFO++ ) {
419
420 numsft = mdetStates->data[iIFO]->length;
421
422 for ( iSFT = 0; iSFT < numsft; iSFT++, k++ ) {
423
424 REAL8 a, b;
425
426 a = multiAMcoef->data[iIFO]->a->data[iSFT];
427 b = multiAMcoef->data[iIFO]->b->data[iSFT];
428 weightsV.data[k] *= ( a * a + b * b );
429 } /* loop over SFTs */
430 } /* loop over IFOs */
431
433
434 XLALDestroyMultiAMCoeffs( multiAMcoef );
435 } /* end AM weights calculation */
436
437
438 /* misc. memory allocations */
439
440 /* memory for one spindown */
441 pulsarTemplate.spindown.length = 1;
442 pulsarTemplate.spindown.data = NULL;
443 pulsarTemplate.spindown.data = ( REAL8 * )LALMalloc( sizeof( REAL8 ) );
444
445 /* copy template parameters */
446 pulsarTemplate.spindown.data[0] = uvar_fdot;
447 pulsarTemplate.f0 = uvar_Freq;
448 pulsarTemplate.latitude = uvar_Delta;
449 pulsarTemplate.longitude = uvar_Alpha;
450
451 /* memory for f(t) vector */
452 foft.length = mObsCoh;
453 foft.data = NULL;
454 foft.data = ( REAL8 * )LALMalloc( mObsCoh * sizeof( REAL8 ) );
455
456 /* memory for peakgram */
457 pg1.length = binsSFT;
458 pg1.data = NULL;
459 pg1.data = ( UCHAR * )LALCalloc( binsSFT, sizeof( UCHAR ) );
460
461
462 /* block for calculating peakgram and number count */
463 {
464 UINT4 iIFO, iSFT;
465 INT4 ind;
466 REAL8 sumWeightSquare;
467 SFTtype *sft;
468
469 /* compute mean and sigma for noise only */
470 /* first calculate the sum of the weights squared */
471 sumWeightSquare = 0.0;
472 for ( k = 0; k < ( INT4 )mObsCoh; k++ ) {
473 sumWeightSquare += weightsV.data[k] * weightsV.data[k];
474 }
475
476 /* probability of selecting a peak expected mean and standard deviation for noise only */
477 alphaPeak = exp( - uvar_peakThreshold );
478 meanN = mObsCoh * alphaPeak;
479 sigmaN = sqrt( sumWeightSquare * alphaPeak * ( 1.0 - alphaPeak ) );
480
481
482 /* the received frequency as a function of time */
483 LAL_CALL( ComputeFoft( &status, &foft, &pulsarTemplate, &timeDiffV, &velV, timeBase ), &status );
484
485 numberCount = 0;
486
487 /* loop over SFT, generate peakgram and get number count */
488 UINT4 j;
489 for ( j = 0, iIFO = 0; iIFO < numifo; iIFO++ ) {
490 numsft = mdetStates->data[iIFO]->length;
491
492 for ( iSFT = 0; iSFT < numsft; iSFT++, j++ ) {
493
494 sft = inputSFTs->data[iIFO]->data + iSFT;
495
496 LAL_CALL( SFTtoUCHARPeakGram( &status, &pg1, sft, uvar_peakThreshold ), &status );
497
498 ind = floor( foft.data[j] * timeBase - sftFminBin + 0.5 );
499
500 numberCount += pg1.data[ind] * weightsV.data[j];
501
502 } /* loop over SFTs */
503 } /* loop over IFOs */
504
505 }
506
507
508 {
509 FILE *fp = NULL;
510 fp = fopen( "./tempout", "w" );
511 fprintf( fp, "%g %g %g\n", ( numberCount - meanN ) / sigmaN, meanN, sigmaN );
512 fprintf( stdout, "%g %g %g\n", ( numberCount - meanN ) / sigmaN, meanN, sigmaN );
513 fclose( fp );
514 }
515
516
517 /* free memory */
518 LALFree( pulsarTemplate.spindown.data );
519 LALFree( timeV.data );
520 LALFree( timeDiffV.data );
521 LALFree( foft.data );
522 LALFree( velV.data );
523
524 LALFree( weightsV.data );
525
527
529
530 XLALDestroyMultiSFTVector( inputSFTs );
531 LALFree( pg1.data );
532
534
536
537 if ( lalDebugLevel ) {
539 }
540
541 return status.statusCode;
542}
543
544
545
546
547
550 REAL8Vector *foft,
551 HoughTemplate *pulsarTemplate,
552 REAL8Vector *timeDiffV,
554 REAL8 timeBase )
555{
556
557 INT4 mObsCoh;
558 REAL8 f0new, vcProdn, timeDiffN;
559 INT4 f0newBin;
560 REAL8 sourceDelta, sourceAlpha, cosDelta;
561 INT4 j, i, nspin, factorialN;
562 REAL8Cart3Coor sourceLocation;
563
564 /* --------------------------------------------- */
567
568 /* Make sure the arguments are not NULL: */
573
577
578 sourceDelta = pulsarTemplate->latitude;
579 sourceAlpha = pulsarTemplate->longitude;
580 cosDelta = cos( sourceDelta );
581
582 sourceLocation.x = cosDelta * cos( sourceAlpha );
583 sourceLocation.y = cosDelta * sin( sourceAlpha );
584 sourceLocation.z = sin( sourceDelta );
585
586 mObsCoh = foft->length;
587 nspin = pulsarTemplate->spindown.length;
588
589 for ( j = 0; j < mObsCoh; ++j ) { /* loop for all different time stamps */
590 vcProdn = velV->data[j].x * sourceLocation.x
591 + velV->data[j].y * sourceLocation.y
592 + velV->data[j].z * sourceLocation.z;
593 f0new = pulsarTemplate->f0;
594 factorialN = 1;
595 timeDiffN = timeDiffV->data[j];
596
597 for ( i = 0; i < nspin; ++i ) { /* loop for spin-down values */
598 factorialN *= ( i + 1 );
599 f0new += pulsarTemplate->spindown.data[i] * timeDiffN / factorialN;
600 timeDiffN *= timeDiffN;
601 }
602 f0newBin = floor( f0new * timeBase + 0.5 );
603 foft->data[j] = f0newBin * ( 1.0 + vcProdn ) / timeBase;
604 }
605
607 /* normal exit */
608 RETURN( status );
609}
Header file for non-demodulated Hough search.
#define DRIVEHOUGHCOLOR_MSGENULL
#define DRIVEHOUGHCOLOR_ENULL
#define LAL_INT4_MAX
void REPORTSTATUS(LALStatus *status)
lal_errhandler_t lal_errhandler
int LAL_ERR_EXIT(LALStatus *stat, const char *func, const char *file, const int line, volatile const char *id)
#define LAL_CALL(function, statusptr)
int j
int k
void LALCheckMemoryLeaks(void)
#define LALCalloc(m, n)
#define LALMalloc(n)
#define LALFree(p)
const LALVCSInfoList lalPulsarVCSInfoList
NULL-terminated list of VCS and build information for LALPulsar and its dependencies
#define ATTATCHSTATUSPTR(statusptr)
#define ASSERT(assertion, statusptr, code, mesg)
#define DETATCHSTATUSPTR(statusptr)
#define INITSTATUS(statusptr)
#define RETURN(statusptr)
void SFTtoUCHARPeakGram(LALStatus *status, UCHARPeakGram *pg, const SFTtype *sft, REAL8 thr)
Convert a normalized sft into a peakgram by selecting bins in which power exceeds a threshold.
Definition: PeakSelect.c:371
#define STRING(a)
#define fprintf
int main(int argc, char *argv[])
BOOLEAN uvar_printSigma
#define F0
#define NFSIZE
#define THRESHOLD
void ComputeFoft(LALStatus *status, REAL8Vector *foft, HoughTemplate *pulsarTemplate, REAL8Vector *timeDiffV, REAL8Cart3CoorVector *velV, REAL8 timeBase)
#define EARTHEPHEMERIS
BOOLEAN uvar_printTemplates
#define TRUE
BOOLEAN uvar_printStats
BOOLEAN uvar_printMaps
BOOLEAN uvar_printEvents
#define BLOCKSRNGMED
#define FBAND
#define SUNEPHEMERIS
#define MAXFILENAMELENGTH
void XLALDestroyMultiDetectorStateSeries(MultiDetectorStateSeries *mdetStates)
Helper function to get rid of a multi-IFO DetectorStateSeries Note, this is "NULL-robust" in the sens...
MultiDetectorStateSeries * XLALGetMultiDetectorStatesFromMultiSFTs(const MultiSFTVector *multiSFTs, const EphemerisData *edat, REAL8 tOffset)
Get the 'detector state' (ie detector-tensor, position, velocity, etc) for the given multi-vector of ...
EphemerisData * XLALInitBarycenter(const CHAR *earthEphemerisFile, const CHAR *sunEphemerisFile)
XLAL interface to reading ephemeris files 'earth' and 'sun', and return ephemeris-data in old backwar...
void XLALDestroyEphemerisData(EphemerisData *edat)
Destructor for EphemerisData struct, NULL robust.
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
unsigned char UCHAR
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
#define lalDebugLevel
void LALHOUGHNormalizeWeights(LALStatus *status, REAL8Vector *weightV)
Normalizes weight factors so that their sum is N.
Definition: DriveHough.c:668
void LALHOUGHInitializeWeights(LALStatus *status, REAL8Vector *weightV)
Initializes weight factors to unity.
Definition: DriveHough.c:633
#define VTOT
Total detector velocity/c TO BE CHANGED DEPENDING ON DETECTOR.
Definition: LUT.h:207
MultiPSDVector * XLALNormalizeMultiSFTVect(MultiSFTVector *multsft, UINT4 blockSize, const MultiNoiseFloor *assumeSqrtSX)
Function for normalizing a multi vector of SFTs in a multi IFO search and returns the running-median ...
void XLALDestroyMultiPSDVector(MultiPSDVector *multvect)
Destroy a multi PSD-vector.
Definition: PSDutils.c:102
MultiNoiseWeights * XLALComputeMultiNoiseWeights(const MultiPSDVector *rngmed, UINT4 blocksRngMed, UINT4 excludePercentile)
Computes weight factors arising from MultiSFTs with different noise floors.
Definition: PSDutils.c:285
void XLALDestroyMultiNoiseWeights(MultiNoiseWeights *weights)
Destroy a MultiNoiseWeights object.
Definition: PSDutils.c:172
void LALCreateRandomParams(LALStatus *status, RandomParams **params, INT4 seed)
static const INT4 a
void LALDestroyRandomParams(LALStatus *status, RandomParams **params)
void LALRemoveKnownLinesInMultiSFTVector(LALStatus *status, MultiSFTVector *MultiSFTVect, INT4 width, INT4 window, LALStringVector *linefiles, RandomParams *randPar)
top level function to remove lines from a multi sft vector given a list of files containing list of k...
Definition: lib/SFTClean.c:953
void XLALDestroySFTCatalog(SFTCatalog *catalog)
Free an 'SFT-catalogue'.
Definition: SFTcatalog.c:329
MultiSFTVector * XLALLoadMultiSFTs(const SFTCatalog *inputCatalog, REAL8 fMin, REAL8 fMax)
Function to load a catalog of SFTs from possibly different detectors.
Definition: SFTfileIO.c:416
void XLALDestroyMultiSFTVector(MultiSFTVector *multvect)
Destroy a multi SFT-vector.
Definition: SFTtypes.c:424
SFTCatalog * XLALSFTdataFind(const CHAR *file_pattern, const SFTConstraints *constraints)
Find the list of SFTs matching the file_pattern and satisfying the given constraints,...
Definition: SFTcatalog.c:71
LIGOTimeGPSVector * XLALReadTimestampsFile(const CHAR *fname)
backwards compatible wrapper to XLALReadTimestampsFileConstrained() without GPS-time constraints
void XLALDestroyTimestampVector(LIGOTimeGPSVector *vect)
De-allocate a LIGOTimeGPSVector.
Definition: SFTtimestamps.c:69
COORDINATESYSTEM_EQUATORIAL
int XLALUserVarReadAllInput(BOOLEAN *should_exit, int argc, char *argv[], const LALVCSInfoList vcs_list)
void XLALDestroyUserVars(void)
#define XLALRegisterNamedUvar(cvar, name, type, option, category,...)
int XLALUserVarWasSet(const void *cvar)
#define XLAL_CHECK_MAIN(assertion,...)
XLAL_SUCCESS
XLAL_EFUNC
LIGOTimeGPS * XLALGPSSetREAL8(LIGOTimeGPS *epoch, REAL8 t)
REAL8 XLALGPSDiff(const LIGOTimeGPS *t1, const LIGOTimeGPS *t0)
int deltaF
REAL8 uvar_fdot
CHAR * uvar_sftDir
REAL8 uvar_startTime
INT4 uvar_blocksRngMed
REAL8 uvar_endTime
REAL4Vector * b
(weighted) per-SFT antenna-pattern function
Definition: LALComputeAM.h:65
REAL4Vector * a
(weighted) per-SFT antenna-pattern function
Definition: LALComputeAM.h:64
COMPLEX8Sequence * data
COMPLEX8FrequencySeries * data
Pointer to the data array.
REAL8 vDetector[3]
Cart.
LIGOTimeGPS tGPS
GPS timestamps corresponding to this entry.
DetectorState * data
array of DetectorState entries
UINT4 length
total number of entries
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
REAL8Vector spindown
A vector of 'timestamps' of type LIGOTimeGPS.
Definition: SFTfileIO.h:188
LIGOTimeGPS * data
array of timestamps
Definition: SFTfileIO.h:193
UINT4 length
number of timestamps
Definition: SFTfileIO.h:192
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
Multi-IFO time-series of DetectorStates.
DetectorStateSeries ** data
vector of pointers to DetectorStateSeries
One noise-weight (number) per SFT (therefore indexed over IFOs and SFTs.
Definition: PSDutils.h:71
REAL8Vector ** data
weights-vector for each detector
Definition: PSDutils.h:76
A collection of PSD vectors – one for each IFO in a multi-IFO search.
Definition: PSDutils.h:62
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
REAL4 * data
Three dimensional Cartessian coordinates.
Definition: LUT.h:308
REAL8 y
Definition: LUT.h:310
REAL8 x
Definition: LUT.h:309
REAL8 z
Definition: LUT.h:311
REAL8Cart3Coor * data
x.y.z
UINT4 length
number of elements
REAL8 * data
An "SFT-catalogue": a vector of SFTdescriptors, as returned by XLALSFTdataFind()
Definition: SFTfileIO.h:238
SFTDescriptor * data
array of data-entries describing matched SFTs
Definition: SFTfileIO.h:243
UINT4 length
number of SFTs in catalog
Definition: SFTfileIO.h:242
'Constraints' for SFT-matching: which detector, within which time-stretch and which timestamps exactl...
Definition: SFTfileIO.h:212
CHAR * detector
2-char channel-prefix describing the detector (eg 'H1', 'H2', 'L1', 'G1' etc)
Definition: SFTfileIO.h:213
LIGOTimeGPSVector * timestamps
list of timestamps
Definition: SFTfileIO.h:216
LIGOTimeGPS * maxStartTime
only include SFTs whose epoch is < maxStartTime
Definition: SFTfileIO.h:215
LIGOTimeGPS * minStartTime
only include SFTs whose epoch is >= minStartTime
Definition: SFTfileIO.h:214
SFTtype header
SFT-header info.
Definition: SFTfileIO.h:228
REAL8 longitude
REAL8 latitude
CoordinateSystem system
Explicit peakgram structure – 1 if power in bin is above threshold and 0 if below.
Definition: PeakSelect.h:114
UCHAR * data
pointer to the data {0,1}
Definition: PeakSelect.h:120
INT4 length
number of elements in data
Definition: PeakSelect.h:118
double f_min
double f_max