Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-b246709
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Weave.c
Go to the documentation of this file.
1//
2// Copyright (C) 2016, 2017 Karl Wette
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_Weave
23///
24
25#include "config.h"
26
27#include "Weave.h"
28#include "SetupData.h"
29#include "SearchIteration.h"
30#include "ComputeResults.h"
31#include "CacheResults.h"
32#include "OutputResults.h"
33#include "SearchTiming.h"
34
35#ifdef LALPULSAR_CUDA_ENABLED
36#include <cuda.h>
37#include <cuda_runtime_api.h>
38#endif
39
40#include <lal/LogPrintf.h>
41#include <lal/UserInput.h>
42#include <lal/Random.h>
43
44int main( int argc, char *argv[] )
45{
46
47 // Set help information
48 lalUserVarHelpBrief = "search for gravitational-wave pulsars";
49
50 ////////// Parse user input //////////
51
52 // Optional arguments to XLALCreateFstatInput(), used to initialise some user input variables
54
55 // Initialise user input variables
56 struct uvar_type {
57 BOOLEAN validate_sft_files, interpolation, lattice_rand_offset, mean2F_hgrm, segment_info, simulate_search, time_search, cache_all_gc, strict_spindown_bounds;
58 CHAR *setup_file, *sft_files, *output_file, *ckpt_output_file;
59 LALStringVector *sft_timestamps_files, *sft_noise_sqrtSX, *injections, *Fstat_assume_sqrtSX, *lrs_oLGX;
60 REAL8 sft_timebase, semi_max_mismatch, coh_max_mismatch, ckpt_output_period, ckpt_output_exit, lrs_Fstar0sc, nc_2Fth;
61 REAL8Range alpha, delta, freq, f1dot, f2dot, f3dot, f4dot;
62 REAL8Vector *random_injection;
63 UINT4 sky_patch_count, sky_patch_index, freq_partitions, f1dot_partitions, Fstat_run_med_window, Fstat_Dterms, toplist_limit, rand_seed, cache_max_size;
64 int lattice, Fstat_method, Fstat_SSB_precision, toplists, extra_statistics, recalc_statistics;
65 } uvar_struct = {
66 .Fstat_Dterms = Fstat_opt_args.Dterms,
67 .Fstat_SSB_precision = Fstat_opt_args.SSBprec,
68 .Fstat_method = FMETHOD_RESAMP_BEST,
69 .Fstat_run_med_window = Fstat_opt_args.runningMedianWindow,
70 .alpha = {0, LAL_TWOPI},
71 .delta = {-LAL_PI_2, LAL_PI_2},
72 .freq_partitions = 1,
73 .f1dot_partitions = 1,
74 .interpolation = 1,
75 .lattice = TILING_LATTICE_ANSTAR,
76 .toplist_limit = 1000,
77 .toplists = WEAVE_STATISTIC_MEAN2F,
78 .extra_statistics = WEAVE_STATISTIC_NONE,
79 .recalc_statistics = WEAVE_STATISTIC_NONE,
80 .nc_2Fth = 5.2,
81 };
82 struct uvar_type *const uvar = &uvar_struct;
83
84 // Register user input variables:
85 //
86 // - General
87 //
89 setup_file, STRING, 'S', REQUIRED,
90 "Setup file generated by lalpulsar_WeaveSetup; the segment list, parameter-space metrics, and other required data. "
91 );
93 output_file, STRING, 'o', REQUIRED,
94 "Output file which stores all quantities computed by lalpulsar_Weave. "
95 );
96 //
97 // - SFT input/generation and signal generation
98 //
99 lalUserVarHelpOptionSubsection = "SFT input/generation and signal generation";
101 sft_files, STRING, 'I', NODEFAULT,
102 "Pattern matching the SFT files to be analysed. Possibilities are:\n"
103 " - '<SFT file>;<SFT file>;...', where <SFT file> may contain wildcards\n - 'list:<file containing list of SFT files>'"
104 );
106 validate_sft_files, BOOLEAN, 'V', DEVELOPER,
107 "Validate the checksums of the SFTs matched by " UVAR_STR( sft_files ) " before loading them into memory. "
108 );
110 sft_timebase, REAL8, 't', NODEFAULT,
111 "Generate SFTs with this timebase instead of loading from files. "
112 );
114 sft_timestamps_files, STRINGVector, 'T', DEVELOPER,
115 "Files containing timestamps for the generated SFTs; if not given, SFTs with contiguous timestamps are generated. "
116 "Arguments correspond to the detectors in the setup file given by " UVAR_STR( setup_file )
117 "; for example, if the setup file was created with " UVAR_STR( detectors ) " set to 'H1,L1', then an argument of "
118 "'t1.txt,t2.txt' to this option will read H1 timestamps from the file 't1.txt', and L1 timestamps from the file 't2.txt'. "
119 "The timebase of the generated SFTs is specified by " UVAR_STR( sft_timebase ) ". "
120 );
122 sft_noise_sqrtSX, STRINGVector, 'p', NODEFAULT,
123 "Inject fake Gaussian noise with these amplitude spectral densities [sqrt(Sh)] into the generated SFTs. "
124 "Arguments correspond to the detectors in the setup file given by " UVAR_STR( setup_file )
125 "; for example, if the setup file was created with " UVAR_STR( detectors ) " set to 'H1,L1', then an argument of "
126 "'1.2,3.4' to this option will generate H1 SFTs with a noise sqrt(Sh) of 1.2, and L1 SFTs with a noise sqrt(Sh) of 3.4. "
127 );
129 injections, STRINGVector, 'J', NODEFAULT,
131 );
133 random_injection, REAL8Vector, 'H', NODEFAULT,
134 "Inject a simulated signal, with phase parameters drawn randomly from the search parameter space, "
135 "and with a randomly-generated polarisation angle (psi) and initial phase (phi0). "
136 "The amplitude of the injection is specified by either 1 or 2 arguments to " UVAR_STR( random_injection ) ":\n"
137 " - h0: generate a random cosine of inclination angle (cosi), then set aPlus=0.5*h0*(1+cosi^2), aCross=h0*cosi;\n"
138 " - aPlus,aCross: use the given plus- and cross-polarisation amplitudes."
139 );
140
141 //
142 // - Search parameter space
143 //
144 lalUserVarHelpOptionSubsection = "Search parameter space";
146 alpha, RAJRange, 'a', NODEFAULT,
147 "Search parameter space in right ascension. "
148 "If not specified, an all-sky search is performed; otherwise " UVAR_STR( delta ) " must also be specified. "
149 "Range for a partial-sky search is limited to PI radians. "
150 );
152 delta, DECJRange, 'd', NODEFAULT,
153 "Search parameter space in declination. "
154 "If not specified, an all-sky search is performed; otherwise " UVAR_STR( alpha ) " must also be specified. "
155 );
157 sky_patch_count, UINT4, 'K', DEVELOPER,
158 "Divide the entire sky into this number of ~equal-template-count patches. "
159 "Requires " UVAR_STR( sky_patch_index ) ", mutually exclusive with " UVAR_STR2AND( alpha, delta ) ". "
160 );
162 sky_patch_index, UINT4, 'k', DEVELOPER,
163 "Search the sky patch given by this index, from zero to one less than " UVAR_STR( sky_patch_count ) ". "
164 "Requires " UVAR_STR( sky_patch_count ) ", mutually exclusive with " UVAR_STR2AND( alpha, delta ) ". "
165 );
167 freq, REAL8Range, 'f', REQUIRED,
168 "Search parameter space in frequency, in Hertz. "
169 );
171 freq_partitions, UINT4, 'F', DEVELOPER,
172 "Internally divide the frequency parameter space into this number of ~equal-width partitions. "
173 );
175 f1dot, REAL8Range, '1', OPTIONAL,
176 "Search parameter space in first spindown, in Hertz/second. "
177 );
179 f1dot_partitions, UINT4, '!', DEVELOPER,
180 "Internally divide the first spindown parameter space into this number of ~equal-width partitions. "
181 );
183 f2dot, REAL8Range, '2', OPTIONAL,
184 "Search parameter space in second spindown, in Hertz/second^2. "
185 );
187 f3dot, REAL8Range, '3', DEVELOPER,
188 "Search parameter space in third spindown, in Hertz/second^3. "
189 );
191 f4dot, REAL8Range, '4', DEVELOPER,
192 "Search parameter space in fourth spindown, in Hertz/second^4. "
193 "(Just in case a nearby supernova goes off!) "
194 );
196 strict_spindown_bounds, BOOLEAN, '0', DEVELOPER,
197 "Do not add padding to spindown parameter space bounds."
198 );
199 //
200 // - Lattice tiling setup
201 //
202 lalUserVarHelpOptionSubsection = "Lattice tiling setup";
204 semi_max_mismatch, REAL8, 's', REQUIRED,
205 "Maximum metric mismatch of the lattice tiling on which semicoherent quantities are computed, e.g. F-statistics averaged over segments. "
206 );
208 coh_max_mismatch, REAL8, 'c', NODEFAULT,
209 "Maximum metric mismatch of the per-segment lattice tilings on which coherent quantities are computed, e.g. coherent F-statistics. "
210 "If the search setup contains only 1 segment, then this option must not be specified. "
211 );
213 interpolation, BOOLEAN, 'i', OPTIONAL,
214 "If TRUE, perform interpolation from the semicoherent lattice tiling to the per-segment coherent lattice tilings; "
215 UVAR_STR( coh_max_mismatch ) " must also be specified in this case. "
216 "If FALSE, turn off interpolation and use the same lattice tiling for both semicoherent and coherent computations; "
217 UVAR_STR( coh_max_mismatch ) " must not be specified in this case. "
218 );
220 lattice, UserEnum, &TilingLatticeChoices, 'l', DEVELOPER,
221 "Type of lattice used to generate the lattice tilings. "
222 );
224 lattice_rand_offset, BOOLEAN, 'j', DEVELOPER,
225 "Offset the physical parameter-space origin of the lattice tilings by a random fraction of the lattice step size. "
226 "This is important when performing mismatch studies to ensure that the mismatch distribution is fully sampled. "
227 );
228 //
229 // - F-statistic computation
230 //
231 lalUserVarHelpOptionSubsection = "F-statistic computation";
233 Fstat_method, UserEnum, XLALFstatMethodChoices(), 'm', DEVELOPER,
234 "Method used to calculate the F-statistic. "
235 );
237 Fstat_run_med_window, UINT4, 'w', DEVELOPER,
238 "Size of the running median window used to normalise SFTs and compute noise weight. "
239 );
241 Fstat_assume_sqrtSX, STRINGVector, 'q', DEVELOPER,
242 "Assume that the noise in the SFTs have known amplitude spectral densities [sqrt(Sh)], which are given by the arguments to "
243 "this option, and normalise the SFTs by these given sqrt(Sh). "
244 "Arguments correspond to the detectors in the setup file given by " UVAR_STR( setup_file )
245 "; for example, if the setup file was created with " UVAR_STR( detectors ) " set to 'H1,L1', then an argument of "
246 "'3.2,4.3' to this option will assume that H1 SFTs contain noise with a sqrt(Sh) of 3.2, and L1 SFTs contain noise with a sqrt(Sh) of 4.3. "
247 "If this option is not given, the SFTs are normalised using noise sqrt(Sh) estimated from the SFTs themselves. "
248 );
250 Fstat_Dterms, UINT4, 0, DEVELOPER,
251 "Number of Dirichlet kernel terms to use in computing the F-statistic. May not be available for all F-statistic methods. "
252 );
254 Fstat_SSB_precision, UserEnum, &SSBprecisionChoices, 0, DEVELOPER,
255 "Precision in calculating the barycentric transformation. "
256 );
257 //
258 // Various statistics input arguments
259 //
260 lalUserVarHelpOptionSubsection = "Various statistics input arguments";
262 lrs_Fstar0sc, REAL8, 0, OPTIONAL,
263 "(Semi-coherent) transition-scale parameter 'Fstar0sc' (=Nseg*Fstar0coh) for B_S/GL.. family of statistics."
264 );
266 lrs_oLGX, STRINGVector, 0, OPTIONAL,
267 "Per-detector line-vs-Gauss prior odds 'oLGX' (Defaults to oLGX=1/Ndet) for B_S/GL.. family of statistics."
268 );
270 nc_2Fth, REAL8, 0, OPTIONAL,
271 "Number count: per-segment 2F threshold value."
272 );
273 //
274 // - Output control
275 //
276 lalUserVarHelpOptionSubsection = "Output control";
278 toplist_limit, UINT4, 'n', OPTIONAL,
279 "Maximum number of candidates to return in an output toplist; if 0, all candidates are returned. "
280 );
282 toplists, UserFlag, &WeaveToplistChoices, 'L', OPTIONAL,
283 "Sets which combination of toplists to return in the output file given by " UVAR_STR( output_file ) ":\n"
285 );
287 mean2F_hgrm, BOOLEAN, 0, DEVELOPER,
288 "Output a histogram of all mean multi-Fstatistics computed by the search. "
289 );
291 extra_statistics, UserFlag, &WeaveStatisticChoices, 'E', OPTIONAL,
292 "Sets which ('stage 0') extra statistics to compute and return in the output file given by " UVAR_STR( output_file ) ":\n"
294 );
295
297 recalc_statistics, UserFlag, &WeaveStatisticChoices, 'R', OPTIONAL,
298 "Sets which extra *recalc* statistics to compute on final toplist without interpolation. See " UVAR_STR( extra_statistics ) " for statistics descriptions.\n"
299 );
300
302 segment_info, BOOLEAN, 'M', DEVELOPER,
303 "Output various information regarding the segment list, e.g. number of SFTs within each segment. "
304 );
305 //
306 // - Checkpointing
307 //
308 lalUserVarHelpOptionSubsection = "Checkpointing";
310 ckpt_output_file, STRING, 'C', DEVELOPER,
311 "File to which to periodically write checkpoints of output results. "
312 );
314 ckpt_output_period, REAL8, 'z', DEVELOPER,
315 "Write checkpoints of output results after this time period (in seconds) has elapsed. "
316 );
318 ckpt_output_exit, REAL8, 0, DEVELOPER,
319 "Write a checkpoint of output results after this fraction of the search has been completed, then exit. "
320 "Arguments to this option must be in the range [0,1]. "
321 "(This option is only really useful for testing the checkpointing feature.) "
322 );
323 //
324 // - Esoterica
325 //
326 lalUserVarHelpOptionSubsection = "Esoterica";
328 rand_seed, UINT4, 'e', DEVELOPER,
329 "Random seed used to initialise random number generators. "
330 );
332 simulate_search, BOOLEAN, 0, DEVELOPER,
333 "Simulate search; perform all search actions apart from computing any results. "
334 "If SFT parameters (i.e. " UVAR_STR( sft_files ) " or " UVAR_STR( sft_timebase ) ") are supplied, simulate search with full memory allocation, i.e. with F-statistic input data, cached coherent results, etc. "
335 "Otherwise, perform search with minimal memory allocation, i.e. do not allocate memory for any data or results. "
336 );
338 time_search, BOOLEAN, 0, DEVELOPER,
339 "Collect and output detailed timing information from various stages of the search pipeline. "
340 );
342 cache_max_size, UINT4, 0, DEVELOPER,
343 "Limit the size of the internal caches, used to store intermediate results, to this number of items per segment. "
344 "If zero, the caches will grow in size to store all items that are still required. "
345 "Has no effect when performing a fully-coherent single-segment search, or a non-interpolating search. "
346 );
348 cache_all_gc, BOOLEAN, 0, DEVELOPER,
349 "If TRUE, try to instead remove as many items as possible, provided that they are no longer required. "
350 "If FALSE, whenever an item is added to the internal caches, at most one item that may no longer be required is removed. "
351 "Has no effect when performing a fully-coherent single-segment search, or a non-interpolating search. "
352 );
353
354 // Parse user input
355 XLAL_CHECK_MAIN( xlalErrno == 0, XLAL_EFUNC, "A call to XLALRegisterUvarMember() failed" );
356 BOOLEAN should_exit = 0;
358 if ( should_exit ) {
359 return EXIT_FAILURE;
360 }
361
362 // Check user input:
363 //
364 // - General
365 //
366
367 //
368 // - SFT input/generation and signal generation
369 //
370 XLALUserVarCheck( &should_exit,
371 uvar->simulate_search || UVAR_SET2( sft_files, sft_timebase ) == 1,
372 "Exactly one of " UVAR_STR2OR( sft_files, sft_timebase ) " must be specified" );
373 XLALUserVarCheck( &should_exit,
374 !UVAR_SET( sft_files ) || !UVAR_ALLSET3( sft_timebase, sft_timestamps_files, sft_noise_sqrtSX ),
375 UVAR_STR( sft_files ) " are mutually exclusive with " UVAR_STR3AND( sft_timebase, sft_timestamps_files, sft_noise_sqrtSX ) );
376 XLALUserVarCheck( &should_exit,
377 !UVAR_SET( validate_sft_files ) || UVAR_SET( sft_files ),
378 UVAR_STR( validate_sft_files ) " requires " UVAR_STR( sft_files ) );
379 XLALUserVarCheck( &should_exit,
380 !UVAR_SET( sft_timebase ) || uvar->sft_timebase > 0,
381 UVAR_STR( sft_timebase ) " must be strictly positive" );
382 XLALUserVarCheck( &should_exit,
383 !UVAR_ALLSET2( injections, random_injection ),
384 UVAR_STR( injections ) " and " UVAR_STR( random_injection ) " are mutually exclusive" );
385 XLALUserVarCheck( &should_exit,
386 !UVAR_SET( random_injection ) || uvar->random_injection->length <= 2,
387 UVAR_STR( random_injection ) " must be passed either 1 or 2 arguments" );
388 //
389 // - Search parameter space
390 //
391 XLALUserVarCheck( &should_exit,
392 -LAL_PI_2 <= uvar->delta[0] && uvar->delta[1] <= LAL_PI_2,
393 UVAR_STR( delta ) " must be within range [-PI/2,PI/2]" );
394 XLALUserVarCheck( &should_exit,
395 !UVAR_SET2( sky_patch_count, sky_patch_index ) || !UVAR_ALLSET2( alpha, delta ),
396 UVAR_STR2AND( sky_patch_count, sky_patch_index ) " are mutually exclusive with " UVAR_STR2AND( alpha, delta ) );
397 XLALUserVarCheck( &should_exit,
398 UVAR_SET2( sky_patch_count, sky_patch_index ) != 1,
399 UVAR_STR( sky_patch_count ) " requires " UVAR_STR( sky_patch_index ) " and vice versa" );
400 XLALUserVarCheck( &should_exit,
401 !UVAR_SET( sky_patch_index ) || uvar->sky_patch_index < uvar->sky_patch_count,
402 UVAR_STR( sky_patch_index ) " must be positive and strictly less than " UVAR_STR( sky_patch_count ) );
403 XLALUserVarCheck( &should_exit,
404 uvar->freq_partitions > 0,
405 UVAR_STR( freq_partitions ) " must be strictly positive" );
406 XLALUserVarCheck( &should_exit,
407 uvar->f1dot_partitions > 0,
408 UVAR_STR( freq_partitions ) " must be strictly positive" );
409 XLALUserVarCheck( &should_exit,
410 !UVAR_SET( f1dot ) || UVAR_SET( freq ),
411 UVAR_STR( freq ) " must be specified if " UVAR_STR( f1dot ) " is specified" );
412 XLALUserVarCheck( &should_exit,
413 !UVAR_SET( f2dot ) || UVAR_SET( f1dot ),
414 UVAR_STR( f1dot ) " must be specified if " UVAR_STR( f2dot ) " is specified" );
415 XLALUserVarCheck( &should_exit,
416 !UVAR_SET( f3dot ) || UVAR_SET( f2dot ),
417 UVAR_STR( f2dot ) " must be specified if " UVAR_STR( f3dot ) " is specified" );
418 XLALUserVarCheck( &should_exit,
419 !UVAR_SET( f4dot ) || UVAR_SET( f3dot ),
420 UVAR_STR( f3dot ) " must be specified if " UVAR_STR( f4dot ) " is specified" );
421 //
422 // - Lattice tiling setup
423 //
424 XLALUserVarCheck( &should_exit,
425 uvar->semi_max_mismatch > 0,
426 UVAR_STR( semi_max_mismatch ) " must be strictly positive" );
427 XLALUserVarCheck( &should_exit,
428 !UVAR_SET( coh_max_mismatch ) || uvar->coh_max_mismatch > 0,
429 UVAR_STR( coh_max_mismatch ) " must be strictly positive" );
430 //
431 // - F-statistic computation
432 //
433 XLALUserVarCheck( &should_exit,
434 uvar->Fstat_SSB_precision < SSBPREC_LAST,
435 UVAR_STR( Fstat_SSB_precision ) " must be in range [0,%u)", SSBPREC_LAST );
436 //
437 // - Output control
438 //
439
440 //
441 // - Checkpointing
442 //
443 XLALUserVarCheck( &should_exit,
444 !UVAR_ALLSET2( ckpt_output_period, ckpt_output_exit ),
445 UVAR_STR2AND( ckpt_output_period, ckpt_output_exit ) " are mutually exclusive" );
446 XLALUserVarCheck( &should_exit,
447 !UVAR_SET( ckpt_output_period ) || uvar->ckpt_output_period > 0,
448 UVAR_STR( ckpt_output_period ) " must be strictly positive" );
449 XLALUserVarCheck( &should_exit,
450 !UVAR_SET( ckpt_output_exit ) || ( 0 <= uvar->ckpt_output_exit && uvar->ckpt_output_exit <= 1 ),
451 UVAR_STR( ckpt_output_exit ) " must be in range [0,1]" );
452 //
453 // - Esoterica
454 //
455 XLALUserVarCheck( &should_exit,
456 !UVAR_ALLSET2( time_search, simulate_search ),
457 UVAR_STR2AND( time_search, simulate_search ) " are mutually exclusive" );
458 XLALUserVarCheck( &should_exit,
459 !UVAR_ALLSET2( time_search, ckpt_output_file ),
460 UVAR_STR2AND( time_search, ckpt_output_file ) " are mutually exclusive" );
461
462 // Exit if required
463 if ( should_exit ) {
464 return EXIT_FAILURE;
465 }
466 LogPrintf( LOG_NORMAL, "Parsed user input successfully\n" );
467
468 // Allocate random number generator
469 RandomParams *rand_par = XLALCreateRandomParams( uvar->rand_seed );
470 XLAL_CHECK_MAIN( rand_par != NULL, XLAL_EFUNC );
471
472 ////////// Load setup data //////////
473
474 // Initialise setup data
475 WeaveSetupData XLAL_INIT_DECL( setup );
476
477 {
478 // Open setup file
479 LogPrintf( LOG_NORMAL, "Opening setup file '%s' for reading ...\n", uvar->setup_file );
480 FITSFile *file = XLALFITSFileOpenRead( uvar->setup_file );
481 XLAL_CHECK_MAIN( file != NULL, XLAL_EFUNC );
482
483 // Read setup data
485
486 // Close output file
488 LogPrintf( LOG_NORMAL, "Closed setup file '%s'\n", uvar->setup_file );
489 }
490
491 // Print reference time
492 LogPrintf( LOG_NORMAL, "Setup file reference time = %" LAL_GPS_FORMAT "\n", LAL_GPS_PRINT( setup.ref_time ) );
493
494 // Number of detectors and segments
495 const UINT4 ndetectors = setup.detectors->length;
496 const UINT4 nsegments = setup.segments->length;
497
498 // Concatenate list of detector into a string
499 char *setup_detectors_string = XLALConcatStringVector( setup.detectors, "," );
500 XLAL_CHECK_MAIN( setup_detectors_string != NULL, XLAL_EFUNC );
501
502 // Compute segment list range
503 LIGOTimeGPS segments_start, segments_end;
504 XLAL_CHECK_MAIN( XLALSegListRange( setup.segments, &segments_start, &segments_end ) == XLAL_SUCCESS, XLAL_EFUNC );
505 LogPrintf( LOG_NORMAL, "Setup file segment list range = [%" LAL_GPS_FORMAT ", %" LAL_GPS_FORMAT "] GPS, segment count = %u\n", LAL_GPS_PRINT( segments_start ), LAL_GPS_PRINT( segments_end ), nsegments );
506
507 ////////// Set up calculation of various requested output statistics //////////
508
509 WeaveStatisticsParams *statistics_params = XLALCalloc( 1, sizeof( *statistics_params ) );
510 XLAL_CHECK_MAIN( statistics_params != NULL, XLAL_ENOMEM );
511 statistics_params->detectors = XLALCopyStringVector( setup.detectors );
512 XLAL_CHECK_MAIN( statistics_params->detectors != NULL, XLAL_EFUNC );
513 statistics_params->nsegments = nsegments;
514 statistics_params->ref_time = setup.ref_time;
515
516 //
517 // Figure out which statistics need to be computed, and when, in order to
518 // produce all the requested toplist-statistics in the "main loop" and all remaining
519 // extra-statistics in the "completion loop" after the toplist has been computed
520 // work out dependency-map for different statistics sets: toplist-ranking, output, total set of dependencies in main/completion loop ...
521 XLAL_CHECK_MAIN( XLALWeaveStatisticsParamsSetDependencyMap( statistics_params, uvar->toplists, uvar->extra_statistics, uvar->recalc_statistics ) == XLAL_SUCCESS, XLAL_EFUNC );
522
523 // Prepare memory for coherent F-stat argument setups
524 statistics_params->coh_input = XLALCalloc( nsegments, sizeof( statistics_params->coh_input[0] ) );
525 XLAL_CHECK_MAIN( statistics_params->coh_input != NULL, XLAL_ENOMEM );
526
527 // ---------- prepare setup for line-robust statistics if requested ----------
528 if ( statistics_params->all_statistics_to_compute & ( WEAVE_STATISTIC_BSGL | WEAVE_STATISTIC_BSGLtL | WEAVE_STATISTIC_BtSGLtL ) ) {
529 REAL4 *oLGX_p = NULL;
531 if ( uvar->lrs_oLGX != NULL ) {
532 XLAL_CHECK_MAIN( uvar->lrs_oLGX->length == ndetectors, XLAL_EINVAL, "length(lrs-oLGX) = %d must equal number of detectors (%d)'\n", uvar->lrs_oLGX->length, ndetectors );
533 XLAL_CHECK_MAIN( XLALParseLinePriors( &oLGX[0], uvar->lrs_oLGX ) == XLAL_SUCCESS, XLAL_EFUNC );
534 oLGX_p = &oLGX[0];
535 }
536 const BOOLEAN useLogCorrection = 0;
537 statistics_params->BSGL_setup = XLALCreateBSGLSetup( ndetectors, uvar->lrs_Fstar0sc, oLGX_p, useLogCorrection, nsegments );
538 XLAL_CHECK_MAIN( statistics_params->BSGL_setup != NULL, XLAL_EFUNC );
539 }
540 // set number-count threshold
541 statistics_params->nc_2Fth = uvar->nc_2Fth;
542
543 // If asked to return a histogram of mean multi-Fstatistics, ensure they're computed
544 if ( uvar->mean2F_hgrm ) {
545 statistics_params->mainloop_statistics |= WEAVE_STATISTIC_MEAN2F;
546 }
547
548 ////////// Set up lattice tilings //////////
549
550 // Check interpolation/maximum mismatch options are consistent with the type of search being performed
551 if ( nsegments == 1 ) {
552 XLALUserVarCheck( &should_exit,
553 !UVAR_SET( coh_max_mismatch ),
554 UVAR_STR( coh_max_mismatch ) " must not be specified if setup file '%s' contains only 1 segment", uvar->setup_file );
555 XLALUserVarCheck( &should_exit,
556 !UVAR_SET( interpolation ) || !uvar->interpolation,
557 UVAR_STR( interpolation ) " must either be FALSE or not specified if setup file '%s' contains only 1 segment", uvar->setup_file );
558 } else if ( uvar->interpolation ) {
559 XLALUserVarCheck( &should_exit,
560 UVAR_SET( coh_max_mismatch ),
561 UVAR_STR( coh_max_mismatch ) " must be specified if " UVAR_STR( interpolation ) " is true" );
562 } else {
563 XLALUserVarCheck( &should_exit,
564 !UVAR_SET( coh_max_mismatch ),
565 UVAR_STR( coh_max_mismatch ) " must not be set if " UVAR_STR( interpolation ) " is false" );
566 }
567 if ( should_exit ) {
568 return EXIT_FAILURE;
569 }
570
571 // Decide which mismatches to use, depending of whether this is an interpolating or non-interpolating search
572 const BOOLEAN interpolation = ( nsegments > 1 ) ? uvar->interpolation : 0;
573 const double semi_max_mismatch = uvar->semi_max_mismatch;
574 const double coh_max_mismatch = interpolation ? uvar->coh_max_mismatch : semi_max_mismatch;
575 if ( nsegments == 1 ) {
576 LogPrintf( LOG_NORMAL, "Performing a fully-coherent single-segment search with maximum (semicoherent) mismatch = %.15g\n", semi_max_mismatch );
577 } else if ( !interpolation ) {
578 LogPrintf( LOG_NORMAL, "Performing a non-interpolating search with maximum (semicoherent) mismatch = %.15g\n", semi_max_mismatch );
579 } else {
580 LogPrintf( LOG_NORMAL, "Performing an interpolating search with maximum semicoherent mismatch = %.15g, maximum coherent mismatch = %.15g\n", semi_max_mismatch, coh_max_mismatch );
581 }
582
583 // Scale metrics to fiducial frequency, given by maximum search frequency
585 LogPrintf( LOG_NORMAL, "Metric fiducial frequency set to maximum search frequency = %.15g Hz\n", uvar->freq[1] );
586
587 // Equalise metric frequency spacing, given the specified maximum mismatches
588 XLAL_CHECK_MAIN( XLALEqualizeReducedSuperskyMetricsFreqSpacing( setup.metrics, coh_max_mismatch, semi_max_mismatch ) == XLAL_SUCCESS, XLAL_EFUNC );
589
590 // Store user input spindown ranges in an array for ease of use
591 const double uvarspins[][2] = {
592 { uvar->f1dot[0], uvar->f1dot[1] },
593 { uvar->f2dot[0], uvar->f2dot[1] },
594 { uvar->f3dot[0], uvar->f3dot[1] },
595 { uvar->f4dot[0], uvar->f4dot[1] }
596 };
597
598 // Check that metrics computed in setup file have sufficient spindown dimensions to cover user input
599 size_t nmetricspins = 0;
600 XLAL_CHECK_MAIN( XLALSuperskyMetricsDimensions( setup.metrics, &nmetricspins ) == XLAL_SUCCESS, XLAL_EFUNC );
601 const size_t nuvarspins = XLAL_NUM_ELEM( uvarspins );
602 XLAL_CHECK_MAIN( nmetricspins <= nuvarspins, XLAL_EINVAL, "Number of spindowns from metrics (%zu) computed in setup file '%s' must be <= %zu", nmetricspins, uvar->setup_file, nuvarspins );
603 const size_t ninputspins = UVAR_SET4( f1dot, f2dot, f3dot, f4dot );
604 XLAL_CHECK_MAIN( ninputspins <= nmetricspins, XLAL_EINVAL, "Number of spindowns from user input (%zu) must be <= number of spindowns from metrics (%zu) computed in setup file '%s'", ninputspins, nmetricspins, uvar->setup_file );
605
606 // Number of parameter-space dimensions: 2 for sky + 1 for frequency + 'nmetricspins' for spindowns
607 const size_t ndim = 2 + 1 + nmetricspins;
608
609 // Number of coherent parameter-space tilings
610 // - If performing a fully-coherent search (i.e. of a single segment), we only need the semicoherent
611 // tiling; otherwise we need coherent tilings for each segment, plus the semicoherent tiling
612 const size_t ncohtiles = ( nsegments > 1 ) ? nsegments : 0;
613
614 // Total number of parameter-space tilings; always number of coherent tilings plus the semicoherent tiling
615 const size_t ntiles = ncohtiles + 1;
616
617 // Index of semicoherent tiling in arrays; always the last element
618 const size_t isemi = ntiles - 1;
619
620 // Create parameter-space tilings
621 LatticeTiling *XLAL_INIT_DECL( tiling, [ntiles] );
622 for ( size_t i = 0; i < ntiles; ++i ) {
624 XLAL_CHECK_MAIN( tiling[i] != NULL, XLAL_EFUNC );
625 }
626
627 // Create arrays to store the appropriate parameter-space metrics for each tiling
628 gsl_matrix *XLAL_INIT_DECL( rssky_metric, [ntiles] );
629 SuperskyTransformData *XLAL_INIT_DECL( rssky_transf, [ntiles] );
630 for ( size_t i = 0; i < nsegments; ++i ) {
631 rssky_metric[i] = setup.metrics->coh_rssky_metric[i];
632 rssky_transf[i] = setup.metrics->coh_rssky_transf[i];
633 }
634 rssky_metric[isemi] = setup.metrics->semi_rssky_metric;
635 rssky_transf[isemi] = setup.metrics->semi_rssky_transf;
636
637 // Create arrays to store the range of physical coordinates covered by each tiling
638 const PulsarDopplerParams *XLAL_INIT_DECL( min_phys, [ntiles] );
639 const PulsarDopplerParams *XLAL_INIT_DECL( max_phys, [ntiles] );
640
641 //
642 // Set up semicoherent lattice tiling
643 //
644
645 // Set sky semicoherent parameter-space bounds
646 // - Compute area of sky semicoherent parameter space for later output
647 double semi_sky_area = 0.0;
648 if ( UVAR_SET( sky_patch_count ) ) {
649 XLAL_CHECK_MAIN( XLALSetSuperskyEqualAreaSkyBounds( tiling[isemi], rssky_metric[isemi], semi_max_mismatch, uvar->sky_patch_count, uvar->sky_patch_index ) == XLAL_SUCCESS, XLAL_EFUNC );
650 LogPrintf( LOG_NORMAL, "Search sky parameter space sky patch = %u of %u\n", uvar->sky_patch_index, uvar->sky_patch_count );
651 semi_sky_area = 4.0 * LAL_PI / uvar->sky_patch_count;
652 } else {
653 XLAL_CHECK_MAIN( XLALSetSuperskyPhysicalSkyBounds( tiling[isemi], rssky_metric[isemi], rssky_transf[isemi], uvar->alpha[0], uvar->alpha[1], uvar->delta[0], uvar->delta[1] ) == XLAL_SUCCESS, XLAL_EFUNC );
654 LogPrintf( LOG_NORMAL, "Search sky parameter space right ascension = [%.15g, %.15g] rad\n", uvar->alpha[0], uvar->alpha[1] );
655 LogPrintf( LOG_NORMAL, "Search sky parameter space declination = [%.15g, %.15g] rad\n", uvar->delta[0], uvar->delta[1] );
656 semi_sky_area = ( uvar->alpha[1] - uvar->alpha[0] ) * ( sin( uvar->delta[1] ) - sin( uvar->delta[0] ) );
657 }
658
659 // Set frequency/spindown semicoherent parameter-space bounds
660 XLAL_CHECK_MAIN( XLALSetSuperskyPhysicalSpinBound( tiling[isemi], rssky_transf[isemi], 0, uvar->freq[0], uvar->freq[1] ) == XLAL_SUCCESS, XLAL_EFUNC );
661 LogPrintf( LOG_NORMAL, "Search frequency parameter space = [%.15g, %.15g] Hz\n", uvar->freq[0], uvar->freq[1] );
662 for ( size_t s = 1; s <= nmetricspins; ++s ) {
663 XLAL_CHECK_MAIN( XLALSetSuperskyPhysicalSpinBound( tiling[isemi], rssky_transf[isemi], s, uvarspins[s - 1][0], uvarspins[s - 1][1] ) == XLAL_SUCCESS, XLAL_EFUNC );
664 XLAL_CHECK_MAIN( XLALSetSuperskyPhysicalSpinBoundPadding( tiling[isemi], rssky_transf[isemi], s, !uvar->strict_spindown_bounds ) == XLAL_SUCCESS, XLAL_EFUNC );
665 LogPrintf( LOG_NORMAL, "Search %zu-order spindown parameter space = [%.15g, %.15g] Hz/s^%zu %s padding\n", s, uvarspins[s - 1][0], uvarspins[s - 1][1], s, uvar->strict_spindown_bounds ? "without" : "with" );
666 }
667
668 // Add random offsets to physical origin of semicoherent lattice tiling, if requested
669 if ( UVAR_SET( lattice_rand_offset ) ) {
671 }
672
673 // Set semicoherent parameter-space lattice and metric
674 XLAL_CHECK_MAIN( XLALSetTilingLatticeAndMetric( tiling[isemi], uvar->lattice, rssky_metric[isemi], semi_max_mismatch ) == XLAL_SUCCESS, XLAL_EFUNC );
675
676 // Print number of (tiled) parameter-space dimensions
677 LogPrintf( LOG_NORMAL, "Number of (tiled) parameter-space dimensions = %zu (%zu)\n", ndim, XLALTiledLatticeTilingDimensions( tiling[isemi] ) );
678
679 // Register callback to compute range of physical coordinates covered by semicoherent parameter space
680 XLAL_CHECK_MAIN( XLALRegisterSuperskyLatticePhysicalRangeCallback( tiling[isemi], rssky_transf[isemi], &min_phys[isemi], &max_phys[isemi] ) == XLAL_SUCCESS, XLAL_EFUNC );
681
682 // Register callbacks to compute, for each coherent tiling, range of coherent reduced supersky coordinates which enclose semicoherent parameter space
683 // - Arrays are of length 'ntiles' since 'ncohtiles' will be zero for a fully-coherent search
684 const gsl_vector *XLAL_INIT_DECL( coh_min_rssky, [ntiles] );
685 const gsl_vector *XLAL_INIT_DECL( coh_max_rssky, [ntiles] );
686 for ( size_t i = 0; i < ncohtiles; ++i ) {
687 XLAL_CHECK_MAIN( XLALRegisterSuperskyLatticeSuperskyRangeCallback( tiling[isemi], rssky_transf[isemi], rssky_transf[i], &coh_min_rssky[i], &coh_max_rssky[i] ) == XLAL_SUCCESS, XLAL_EFUNC );
688 }
689
690 // Iterate over semicoherent tiling and perform callback actions
691 LogPrintf( LOG_NORMAL, "Setting up semicoherent lattice tiling ...\n" );
693 LogPrintf( LOG_NORMAL, "Finished setting up semicoherent lattice tiling\n" );
694
695 // Output number of points in semicoherent tiling
696 {
697 const LatticeTilingStats *stats = XLALLatticeTilingStatistics( tiling[isemi], ndim - 1 );
698 XLAL_CHECK_MAIN( stats != NULL, XLAL_EFUNC );
699 LogPrintf( LOG_NORMAL, "Number of semicoherent templates = %" LAL_UINT8_FORMAT "\n", stats->total_points );
700 }
701
702 // Get frequency spacing used by parameter-space tiling
703 // - XLALEqualizeReducedSuperskyMetricsFreqSpacing() ensures this is the same for all segments
704 const double dfreq = XLALLatticeTilingStepSize( tiling[isemi], ndim - 1 );
705
706 //
707 // Set up coherent lattice tilings
708 //
709 LogPrintf( LOG_NORMAL, "Setting up coherent lattice tilings ...\n" );
710 for ( size_t i = 0; i < ncohtiles; ++i ) {
711
712 // Set coherent parameter-space bounds which enclose semicoherent parameter space
713 XLAL_CHECK_MAIN( XLALSetSuperskyRangeBounds( tiling[i], coh_min_rssky[i], coh_max_rssky[i] ) == XLAL_SUCCESS, XLAL_EFUNC );
714
715 // Add random offsets to physical origin of coherent lattice tiling, if requested
716 if ( UVAR_SET( lattice_rand_offset ) ) {
718 }
719
720 // Ensure that coherent and semicoherent lattice tilings have the same tiled/non-tiled dimensions
722
723 // Set coherent parameter-space lattice and metric
724 XLAL_CHECK_MAIN( XLALSetTilingLatticeAndMetric( tiling[i], uvar->lattice, rssky_metric[i], coh_max_mismatch ) == XLAL_SUCCESS, XLAL_EFUNC );
725
726 // Register callback to compute range of physical coordinates covered by coherent parameter space
727 XLAL_CHECK_MAIN( XLALRegisterSuperskyLatticePhysicalRangeCallback( tiling[i], rssky_transf[i], &min_phys[i], &max_phys[i] ) == XLAL_SUCCESS, XLAL_EFUNC );
728
729 // Iterate over coherent tiling and perform callback actions
731
732 }
733 LogPrintf( LOG_NORMAL, "Finished setting up coherent lattice tilings\n" );
734
735 ////////// Load input data //////////
736
737 // Load or generate SFTs, unless search is being simulated
738 SFTCatalog *sft_catalog = NULL;
739 if ( UVAR_SET( sft_files ) ) {
740
741 // Load SFT catalog from files given by 'sft_files'
742 LogPrintf( LOG_NORMAL, "Loading SFTs matching '%s' into catalog ...\n", uvar->sft_files );
743 sft_catalog = XLALSFTdataFind( uvar->sft_files, NULL );
744 XLAL_CHECK_MAIN( sft_catalog != NULL, XLAL_EFUNC );
745 XLAL_CHECK_MAIN( sft_catalog->length > 0, XLAL_EFUNC );
746 LogPrintf( LOG_NORMAL, "Loaded SFT catalog from SFTs matching '%s'\n", uvar->sft_files );
747
748 // Validate checksums of SFTs, if requested
749 if ( uvar->validate_sft_files ) {
750 LogPrintf( LOG_NORMAL, "Validating SFTs ...\n" );
751 BOOLEAN crc_check = 0;
752 XLAL_CHECK_MAIN( XLALCheckCRCSFTCatalog( &crc_check, sft_catalog ) == XLAL_SUCCESS, XLAL_EFUNC );
753 XLAL_CHECK_MAIN( crc_check, XLAL_EFUNC, "Failed to validate checksums of SFTs" );
754 LogPrintf( LOG_NORMAL, "Finished validating SFTs\n" );
755 }
756
757 } else if ( UVAR_SET( sft_timebase ) ) {
758
759 // Create timestamps for generated SFTs
760 MultiLIGOTimeGPSVector *sft_timestamps = NULL;
761 if ( UVAR_SET( sft_timestamps_files ) ) {
762
763 // Check that the number of SFT timestamp files is consistent with the number of detectors
764 XLAL_CHECK_MAIN( uvar->sft_timestamps_files->length == ndetectors, XLAL_EINVAL, "Number SFT timestamp files (%i) is inconsistent with number of detectors (%i) in setup file '%s'", uvar->sft_timestamps_files->length, ndetectors, uvar->setup_file );
765
766 // Load SFT timestamps from files given by 'sft_timestamps_files'
767 sft_timestamps = XLALReadMultiTimestampsFiles( uvar->sft_timestamps_files );
768 XLAL_CHECK_MAIN( sft_timestamps != NULL, XLAL_EFUNC );
769 for ( size_t i = 0; i < ndetectors; ++i ) {
770 sft_timestamps->data[i]->deltaT = uvar->sft_timebase;
771 LogPrintf( LOG_NORMAL, "Loaded SFT timestamps for detector '%s' from file '%s'\n", setup.detectors->data[i], uvar->sft_timestamps_files->data[i] );
772 }
773
774 } else {
775
776 // Generate identical SFT timestamps for each detector, starting from beginning of segment list, with timebase given by 'sft_timebase'
777 sft_timestamps = XLALMakeMultiTimestamps( segments_start, XLALGPSDiff( &segments_end, &segments_start ), uvar->sft_timebase, 0, ndetectors );
778 XLAL_CHECK_MAIN( sft_timestamps != NULL, XLAL_EFUNC );
779 LogPrintf( LOG_NORMAL, "Generated SFT timestamps for %i detectors, timebase = %.15g sec\n", ndetectors, uvar->sft_timebase );
780
781 }
782
783 // Generate SFT catalog for detectors 'sft_detectors' and timestamps 'sft_timestamps'
784 sft_catalog = XLALMultiAddToFakeSFTCatalog( sft_catalog, setup.detectors, sft_timestamps );
785 XLAL_CHECK_MAIN( sft_catalog != NULL, XLAL_EFUNC );
786 XLAL_CHECK_MAIN( sft_catalog->length > 0, XLAL_EFUNC );
787
788 // Cleanup
789 XLALDestroyMultiTimestamps( sft_timestamps );
790
791 }
792 if ( sft_catalog != NULL ) {
793
794 // Check that all SFT catalog detectors were included in setup file
795 LALStringVector *sft_catalog_detectors = XLALListIFOsInCatalog( sft_catalog );
796 XLAL_CHECK_MAIN( sft_catalog_detectors != NULL, XLAL_EFUNC );
797 char *sft_catalog_detectors_string = XLALConcatStringVector( sft_catalog_detectors, "," );
798 XLAL_CHECK_MAIN( sft_catalog_detectors_string != NULL, XLAL_EFUNC );
799 XLAL_CHECK_MAIN( strcmp( sft_catalog_detectors_string, setup_detectors_string ) == 0, XLAL_EINVAL, "List of detectors '%s' in SFT catalog differs from list of detectors '%s' in setup file '%s'", sft_catalog_detectors_string, setup_detectors_string, uvar->setup_file );
800
801 // Log number of SFTs, both in total and for each detector
802 MultiSFTCatalogView *sft_catalog_view = XLALGetMultiSFTCatalogView( sft_catalog );
803 XLAL_CHECK_MAIN( sft_catalog_view != NULL, XLAL_EINVAL );
804 XLAL_CHECK_MAIN( sft_catalog_view->length > 0, XLAL_EFUNC );
805 LogPrintf( LOG_NORMAL, "Using %u SFTs in total", sft_catalog->length );
806 for ( size_t j = 0; j < sft_catalog_view->length; ++j ) {
807 XLAL_CHECK_MAIN( sft_catalog_view->data[j].length > 0, XLAL_EFUNC );
808 char *detector_name = XLALGetChannelPrefix( sft_catalog_view->data[j].data[0].header.name );
809 XLAL_CHECK_MAIN( detector_name != NULL, XLAL_EFUNC );
810 XLAL_CHECK_MAIN( XLALFindStringInVector( detector_name, setup.detectors ) >= 0, XLAL_EFAILED );
811 LogPrintfVerbatim( LOG_NORMAL, ", %u SFTs from detector '%s'", sft_catalog_view->data[j].length, detector_name );
812 XLALFree( detector_name );
813 }
815
816 // Cleanup
817 XLALDestroyStringVector( sft_catalog_detectors );
818 XLALFree( sft_catalog_detectors_string );
819 XLALDestroyMultiSFTCatalogView( sft_catalog_view );
820
821 }
822
823 // Decide on search simulation level
824 const WeaveSimulationLevel simulation_level = uvar->simulate_search ? ( WEAVE_SIMULATE | ( sft_catalog == NULL ? WEAVE_SIMULATE_MIN_MEM : 0 ) ) : 0;
825 if ( simulation_level & WEAVE_SIMULATE ) {
826 LogPrintf( LOG_NORMAL, "Simulating search with %s memory usage; no results will be computed\n", simulation_level & WEAVE_SIMULATE_MIN_MEM ? "minimal" : "full" );
827 }
828
829 // Set up injections
831 if ( UVAR_SET( injections ) ) {
832
833 // Parse signal injection string
834 injections = XLALPulsarParamsFromUserInput( uvar->injections, &setup.ref_time );
836
837 } else if ( UVAR_SET( random_injection ) ) {
838
839 // Create injection parameters vector
842 PulsarParams *const inj = &injections->data[0];
843 strcpy( inj->name, "random_injection" );
844
845 // Draw random phase parameters from search parameter space
846 {
847 double XLAL_INIT_DECL( random_point, [ndim] );
848 gsl_matrix_view random_point_matrix = gsl_matrix_view_array( random_point, ndim, 1 );
849 XLAL_CHECK_MAIN( XLALRandomLatticeTilingPoints( tiling[isemi], 0, rand_par, &random_point_matrix.matrix ) == XLAL_SUCCESS, XLAL_EFUNC );
850 gsl_vector_view random_point_vector = gsl_vector_view_array( random_point, ndim );
851 XLAL_CHECK_MAIN( XLALConvertSuperskyToPhysicalPoint( &inj->Doppler, &random_point_vector.vector, NULL, rssky_transf[isemi] ) == XLAL_SUCCESS, XLAL_EFUNC );
852 }
853
854 // Randomly generate and set amplitude parameters
855 inj->Amp.psi = LAL_TWOPI * XLALUniformDeviate( rand_par );
856 inj->Amp.phi0 = LAL_TWOPI * XLALUniformDeviate( rand_par );
857 if ( uvar->random_injection->length == 1 ) {
858 const REAL8 h0 = uvar->random_injection->data[0];
859 const REAL8 cosi = -1 + 2 * XLALUniformDeviate( rand_par );
860 inj->Amp.aPlus = 0.5 * h0 * ( 1.0 + cosi * cosi );
861 inj->Amp.aCross = h0 * cosi;
862 } else {
863 inj->Amp.aPlus = uvar->random_injection->data[0];
864 inj->Amp.aCross = uvar->random_injection->data[1];
865 }
866
867 }
868
869 // Set F-statistic optional arguments
870 Fstat_opt_args.randSeed = uvar->rand_seed;
871 Fstat_opt_args.SSBprec = uvar->Fstat_SSB_precision;
872 Fstat_opt_args.Dterms = uvar->Fstat_Dterms;
873 Fstat_opt_args.runningMedianWindow = uvar->Fstat_run_med_window;
874 Fstat_opt_args.FstatMethod = uvar->Fstat_method;
875 Fstat_opt_args.injectSources = injections;
876 Fstat_opt_args.prevInput = NULL;
877 Fstat_opt_args.collectTiming = uvar->time_search;
878
879 // Load input data required for computing coherent results
880 const LALStringVector *sft_noise_sqrtSX = UVAR_SET( sft_noise_sqrtSX ) ? uvar->sft_noise_sqrtSX : NULL;
881 const LALStringVector *Fstat_assume_sqrtSX = UVAR_SET( Fstat_assume_sqrtSX ) ? uvar->Fstat_assume_sqrtSX : NULL;
882 LogPrintf( LOG_NORMAL, "Loading input data for coherent results ...\n" );
883 XLAL_INIT_MEM( statistics_params->n2F_det );
884 for ( size_t i = 0; i < nsegments; ++i ) {
885 statistics_params->coh_input[i] = XLALWeaveCohInputCreate( setup.detectors, simulation_level, sft_catalog, i, &setup.segments->segs[i], min_phys[i], max_phys[i], dfreq, setup.ephemerides, sft_noise_sqrtSX, Fstat_assume_sqrtSX, &Fstat_opt_args, statistics_params, 0 );
886 XLAL_CHECK_MAIN( statistics_params->coh_input[i] != NULL, XLAL_EFUNC );
887 }
888 if ( !( simulation_level & WEAVE_SIMULATE_MIN_MEM ) && ( statistics_params->mainloop_statistics & WEAVE_STATISTIC_COH2F_DET ) ) {
889 for ( size_t i = 0; i < ndetectors; ++i ) {
890 XLAL_CHECK_MAIN( 0 < statistics_params->n2F_det[i] && statistics_params->n2F_det[i] <= nsegments, XLAL_EFAILED, "Invalid number of per-detector F-statistics (%u) for detector '%s'", statistics_params->n2F_det[i], setup.detectors->data[i] );
891 }
892 }
893
894 LogPrintf( LOG_NORMAL, "Finished loading input data for coherent results\n" );
895
896 // Create caches to store intermediate results from coherent parameter-space tilings
897 // - If no interpolation, caching is not required so reduce maximum cache size to 1
898 WeaveCache *XLAL_INIT_DECL( coh_cache, [nsegments] );
899 for ( size_t i = 0; i < nsegments; ++i ) {
900 const size_t cache_max_size = interpolation ? uvar->cache_max_size : 1;
901 const BOOLEAN cache_all_gc = interpolation ? uvar->cache_all_gc : 0;
902 coh_cache[i] = XLALWeaveCacheCreate( tiling[i], interpolation, rssky_transf[i], rssky_transf[isemi], statistics_params->coh_input[i], cache_max_size, cache_all_gc );
903 XLAL_CHECK_MAIN( coh_cache[i] != NULL, XLAL_EFUNC );
904 }
905
906 ////////// Perform search //////////
907
908 // Create iterator over the main loop search parameter space
909 WeaveSearchIterator *main_loop_itr = XLALWeaveMainLoopSearchIteratorCreate( tiling[isemi], uvar->freq_partitions, uvar->f1dot_partitions );
910 XLAL_CHECK_MAIN( main_loop_itr != NULL, XLAL_EFUNC );
911
912 // Create storage for cache queries for coherent results in each segment
913 WeaveCacheQueries *queries = XLALWeaveCacheQueriesCreate( tiling[isemi], rssky_transf[isemi], dfreq, nsegments, uvar->freq_partitions );
914 XLAL_CHECK_MAIN( queries != NULL, XLAL_EFUNC );
915
916 // Pointer to final semicoherent results
917 WeaveSemiResults *semi_res = NULL;
918
919 // Create output results structure
920 WeaveOutputResults *out = XLALWeaveOutputResultsCreate( &setup.ref_time, ninputspins, statistics_params, uvar->toplist_limit, uvar->mean2F_hgrm );
921 XLAL_CHECK_MAIN( out != NULL, XLAL_EFUNC );
922
923 // Create search timing structure
924 WeaveSearchTiming *tim = XLALWeaveSearchTimingCreate( uvar->time_search, statistics_params );
925 XLAL_CHECK_MAIN( tim != NULL, XLAL_EFUNC );
926
927 // Number of times output results have been restored from a checkpoint
928 UINT4 ckpt_output_count = 0;
929
930 // Try to restore output results from a checkpoint file, if given
931 if ( UVAR_SET( ckpt_output_file ) ) {
932
933 // Try to open output checkpoint file
934 LogPrintf( LOG_NORMAL, "Trying to open output checkpoint file '%s' for reading ...\n", uvar->ckpt_output_file );
935 int errnum = 0;
936 FITSFile *file = NULL;
937 XLAL_TRY( file = XLALFITSFileOpenRead( uvar->ckpt_output_file ), errnum );
938 if ( errnum == XLAL_ENOENT ) {
939 LogPrintf( LOG_NORMAL, "Output checkpoint file '%s' does not exist; no checkpoint will be loaded\n", uvar->ckpt_output_file );
940 } else {
941 XLAL_CHECK_MAIN( errnum == 0 && file != NULL, XLAL_EFUNC );
942 LogPrintf( LOG_NORMAL, "Output checkpoint file '%s' exists; checkpoint will be loaded\n", uvar->ckpt_output_file );
943
944 // Read number of times output results have been restored from a checkpoint
945 XLAL_CHECK_MAIN( XLALFITSHeaderReadUINT4( file, "ckptcnt", &ckpt_output_count ) == XLAL_SUCCESS, XLAL_EFUNC );
946 XLAL_CHECK_MAIN( ckpt_output_count > 0, XLAL_EIO, "Invalid output checkpoint file '%s'", uvar->ckpt_output_file );
947
948 // Read maximum size of toplists
949 UINT4 toplist_limit = 0;
950 XLAL_CHECK( XLALFITSHeaderReadUINT4( file, "toplimit", &toplist_limit ) == XLAL_SUCCESS, XLAL_EFUNC );
951
952 // Read output results
953 XLAL_CHECK_MAIN( XLALWeaveOutputResultsReadAppend( file, &out, toplist_limit ) == XLAL_SUCCESS, XLAL_EFUNC, "Invalid output checkpoint file '%s'", uvar->ckpt_output_file );
954
955 // Restore state of main loop iterator
956 XLAL_CHECK_MAIN( XLALWeaveSearchIteratorRestore( main_loop_itr, file ) == XLAL_SUCCESS, XLAL_EFUNC, "Invalid output checkpoint file '%s'", uvar->ckpt_output_file );
957
958 // Close output checkpoint file
960 LogPrintf( LOG_NORMAL, "Closed output checkpoint file '%s'\n", uvar->ckpt_output_file );
961
962 }
963
964 }
965
966 // Start timing main search loop
968
969 // Elapsed wall time at which search was last checkpointed
970 double wall_ckpt_elapsed = 0;
971
972 // Elapsed wall time at which progress was last printed, and interval at which to print progress
973 double wall_prog_elapsed = 0;
974 double wall_prog_period = 5.0;
975
976 // Whether to print predicted remaining time, and previous prediction for total elapsed time
977 BOOLEAN wall_prog_remain_print = 0;
978 double wall_prog_total_prev = 0;
979
980 // Print initial progress
981 LogPrintf( LOG_NORMAL, "Starting main loop at %.3g%% complete, peak memory %.1fMB\n", XLALWeaveSearchIteratorProgress( main_loop_itr ), XLALGetPeakHeapUsageMB() );
982
983 // Begin main loop
984 BOOLEAN search_complete = 0;
985 while ( !search_complete ) {
986
987 // Switch timing section
989
990 // Get next semicoherent frequency block
991 // - Exit main loop if iteration is complete
992 // - Expire cache items if requested by iterator
993 BOOLEAN expire_cache = 0;
994 UINT8 semi_index = 0;
995 const gsl_vector *semi_rssky = NULL;
996 INT4 semi_left = 0;
997 INT4 semi_right = 0;
998 UINT4 freq_partition_index = 0;
999 XLAL_CHECK( XLALWeaveSearchIteratorNext( main_loop_itr, &search_complete, &expire_cache, &semi_index, &semi_rssky, &semi_left, &semi_right, &freq_partition_index ) == XLAL_SUCCESS, XLAL_EFUNC );
1000 if ( search_complete ) {
1002 break;
1003 } else if ( expire_cache ) {
1004 for ( size_t i = 0; i < nsegments; ++i ) {
1006 }
1007 }
1008
1009 // Switch timing section
1011
1012 // Initialise cache queries
1013 XLAL_CHECK_MAIN( XLALWeaveCacheQueriesInit( queries, semi_index, semi_rssky, semi_left, semi_right, freq_partition_index ) == XLAL_SUCCESS, XLAL_EFUNC );
1014
1015 // Query for coherent results for each segment
1016 for ( size_t i = 0; i < nsegments; ++i ) {
1017 XLAL_CHECK_MAIN( XLALWeaveCacheQuery( coh_cache[i], queries, i ) == XLAL_SUCCESS, XLAL_EFUNC );
1018 }
1019
1020 // Finalise cache queries
1022 UINT4 semi_nfreqs = 0;
1023 XLAL_CHECK_MAIN( XLALWeaveCacheQueriesFinal( queries, &semi_phys, &semi_nfreqs ) == XLAL_SUCCESS, XLAL_EFUNC );
1024 if ( semi_nfreqs == 0 ) {
1026 continue;
1027 }
1028
1029 // Switch timing section
1031
1032 // Retrieve coherent results from each segment
1033 const WeaveCohResults *XLAL_INIT_DECL( coh_res, [nsegments] );
1034 UINT8 XLAL_INIT_DECL( coh_index, [nsegments] );
1035 UINT4 XLAL_INIT_DECL( coh_offset, [nsegments] );
1036 for ( size_t i = 0; i < nsegments; ++i ) {
1037 XLAL_CHECK_MAIN( XLALWeaveCacheRetrieve( coh_cache[i], queries, i, &coh_res[i], &coh_index[i], &coh_offset[i], tim ) == XLAL_SUCCESS, XLAL_EFUNC );
1038 XLAL_CHECK_MAIN( coh_res[i] != NULL, XLAL_EFUNC );
1039 }
1040
1041 // Switch timing section
1043
1044 // Initialise semicoherent results
1045 XLAL_CHECK_MAIN( XLALWeaveSemiResultsInit( &semi_res, simulation_level, ndetectors, nsegments, semi_index, &semi_phys, dfreq, semi_nfreqs, statistics_params ) == XLAL_SUCCESS, XLAL_EFUNC );
1046
1047 // Add coherent results to semicoherent results
1048 XLAL_CHECK_MAIN( XLALWeaveSemiResultsComputeSegs( semi_res, nsegments, coh_res, coh_index, coh_offset, tim ) == XLAL_SUCCESS, XLAL_EFUNC );
1049
1050 // Switch timing section
1052
1053 // Compute all toplist-ranking semicoherent results
1055
1056 // Switch timing section
1058
1059 // Add semicoherent results to output
1060 XLAL_CHECK_MAIN( XLALWeaveOutputResultsAdd( out, semi_res, semi_nfreqs ) == XLAL_SUCCESS, XLAL_EFUNC );
1061
1062 // Switch timing section
1064
1065 // Main iterator percentage complete
1066 const REAL4 prog_per_cent = XLALWeaveSearchIteratorProgress( main_loop_itr );
1067
1068 // Current elapsed wall and CPU times
1069 double wall_elapsed = 0, cpu_elapsed = 0;
1070 XLAL_CHECK_MAIN( XLALWeaveSearchTimingElapsed( tim, &wall_elapsed, &cpu_elapsed ) == XLAL_SUCCESS, XLAL_EFUNC );
1071
1072 // Print iteration progress, if required
1073 if ( wall_elapsed - wall_prog_elapsed >= wall_prog_period ) {
1074
1075 // Print progress
1076 LogPrintf( LOG_NORMAL, "%s at %.3g%% complete", simulation_level & WEAVE_SIMULATE ? "Simulation" : "Search", prog_per_cent );
1077
1078 // Print elapsed time
1079 LogPrintfVerbatim( LOG_NORMAL, ", elapsed %.1f sec", wall_elapsed );
1080
1081 // Print remaining time, if it can be reliably predicted
1082 const double wall_prog_remain = XLALWeaveSearchIteratorRemainingTime( main_loop_itr, wall_elapsed );
1083 const double wall_prog_total = wall_elapsed + wall_prog_remain;
1084 if ( wall_prog_remain_print || fabs( wall_prog_total - wall_prog_total_prev ) <= 0.1 * wall_prog_total_prev ) {
1085 LogPrintfVerbatim( LOG_NORMAL, ", remaining ~%.1f sec", wall_prog_remain );
1086 wall_prog_remain_print = 1; // Always print remaining time once it can be reliably predicted
1087 } else {
1088 wall_prog_total_prev = wall_prog_total;
1089 }
1090
1091 // Print CPU usage
1092 LogPrintfVerbatim( LOG_NORMAL, ", CPU %.1f%%", 100.0 * cpu_elapsed / wall_elapsed );
1093
1094 // Print memory usage
1095 LogPrintfVerbatim( LOG_NORMAL, ", peak memory %.1fMB", XLALGetPeakHeapUsageMB() );
1096#ifdef LALPULSAR_CUDA_ENABLED
1097 if ( Fstat_opt_args.FstatMethod == FMETHOD_RESAMP_CUDA ) {
1098 size_t CUDA_free_mem_B = 0;
1099 size_t CUDA_tot_mem_B = 0;
1100 XLAL_CHECK_MAIN( cudaMemGetInfo( &CUDA_free_mem_B, &CUDA_tot_mem_B ) == cudaSuccess, XLAL_EERR );
1101 XLAL_CHECK_MAIN( CUDA_free_mem_B <= CUDA_tot_mem_B, XLAL_EERR );
1102 const size_t CUDA_used_mem_B = CUDA_tot_mem_B - CUDA_free_mem_B;
1103 const REAL4 CUDA_used_mem_MB = CUDA_used_mem_B / ( 1024.0 * 1024.0 );
1104 const REAL4 CUDA_tot_mem_MB = CUDA_tot_mem_B / ( 1024.0 * 1024.0 );
1105 LogPrintfVerbatim( LOG_NORMAL, ", CUDA memory %.1f/%.1fMB", CUDA_used_mem_MB, CUDA_tot_mem_MB );
1106 }
1107#endif
1108
1109 // Print mean maximum size obtained by caches
1110 {
1111 REAL4 cache_mean_max_size = 0;
1112 XLAL_CHECK_MAIN( XLALWeaveGetCacheMeanMaxSize( &cache_mean_max_size, nsegments, coh_cache ) == XLAL_SUCCESS, XLAL_EFUNC );
1113 LogPrintfVerbatim( LOG_NORMAL, ", cache max ~%0.1f", cache_mean_max_size );
1114 }
1115
1116 // Finish progress printing
1118
1119 // Update elapsed wall time at which progress was last printed, and increase interval at which to print progress
1120 wall_prog_elapsed = wall_elapsed;
1121 wall_prog_period = GSL_MIN( 1200, wall_prog_period * 1.5 );
1122
1123 }
1124
1125 // Checkpoint output results, if required
1126 if ( UVAR_SET( ckpt_output_file ) ) {
1127
1128 // Switch timing section
1130
1131 // Decide whether to checkpoint output results
1132 const BOOLEAN do_ckpt_output_period = UVAR_SET( ckpt_output_period ) && wall_elapsed - wall_ckpt_elapsed >= uvar->ckpt_output_period;
1133 const BOOLEAN do_ckpt_output_exit = UVAR_SET( ckpt_output_exit ) && prog_per_cent >= 100.0 * uvar->ckpt_output_exit;
1134 if ( do_ckpt_output_period || do_ckpt_output_exit ) {
1135
1136 // Open output checkpoint file
1137 FITSFile *file = XLALFITSFileOpenWrite( uvar->ckpt_output_file );
1138 XLAL_CHECK_MAIN( file != NULL, XLAL_EFUNC );
1141
1142 // Write number of times output results have been restored from a checkpoint
1143 ++ckpt_output_count;
1144 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT4( file, "ckptcnt", ckpt_output_count, "number of checkpoints" ) == XLAL_SUCCESS, XLAL_EFUNC );
1145
1146 // Write output results
1148
1149 // Save state of main loop iterator
1151
1152 // Close output checkpoint file
1154
1155 // Print progress
1156 LogPrintf( LOG_NORMAL, "Wrote output checkpoint to file '%s' at %.3g%% complete, elapsed %.1f sec\n", uvar->ckpt_output_file, prog_per_cent, wall_elapsed );
1157
1158 // Exit main loop, if checkpointing was triggered by 'do_ckpt_output_exit'
1159 if ( do_ckpt_output_exit ) {
1160 LogPrintf( LOG_NORMAL, "Exiting main seach loop after writing output checkpoint\n" );
1162 break;
1163 }
1164
1165 // Update elapsed wall time at which search was last checkpointed
1166 wall_ckpt_elapsed = wall_elapsed;
1167
1168 }
1169
1170 // Switch timing section
1172
1173 }
1174
1175 } // End of main loop
1176
1177 // Clear all cache items from memory
1178 for ( size_t i = 0; i < nsegments; ++i ) {
1180 }
1181
1182 // Print progress
1183 double wall_main = 0, cpu_main = 0;
1184 XLAL_CHECK_MAIN( XLALWeaveSearchTimingElapsed( tim, &wall_main, &cpu_main ) == XLAL_SUCCESS, XLAL_EFUNC );
1185 LogPrintf( LOG_NORMAL, "Finished main loop at %.3g%% complete, main-loop time %.1f sec, CPU %.1f%%, peak memory %.1fMB\n", XLALWeaveSearchIteratorProgress( main_loop_itr ), wall_main, 100.0 * cpu_main / wall_main, XLALGetPeakHeapUsageMB() );
1186
1187 // Prepare completion-loop calculations:
1188 // if any 'recalc' (= 'stage 1') statistics have been requested: we'll need 'Demod' Fstatistic setups
1189 // so if the 'stage 0' calculation used Demod => nothing to do, if it used 'Resamp' then re-compute the setups
1190 if ( ( uvar->recalc_statistics != WEAVE_STATISTIC_NONE ) ) {
1191 statistics_params->coh_input_recalc = XLALCalloc( nsegments, sizeof( statistics_params->coh_input_recalc[0] ) );
1192 XLAL_CHECK_MAIN( statistics_params->coh_input_recalc != NULL, XLAL_ENOMEM );
1193 FstatOptionalArgs Fstat_opt_args_recalc = Fstat_opt_args;
1194 Fstat_opt_args_recalc.FstatMethod = FMETHOD_DEMOD_BEST;
1195 Fstat_opt_args_recalc.prevInput = NULL;
1196 XLAL_INIT_MEM( statistics_params->n2F_det );
1197 for ( size_t i = 0; i < nsegments; ++i ) {
1198 statistics_params->coh_input_recalc[i] = XLALWeaveCohInputCreate( setup.detectors, simulation_level, sft_catalog, i, &setup.segments->segs[i], min_phys[i], max_phys[i], 0, setup.ephemerides, sft_noise_sqrtSX, Fstat_assume_sqrtSX, &Fstat_opt_args_recalc, statistics_params, 1 );
1199 XLAL_CHECK_MAIN( statistics_params->coh_input_recalc[i] != NULL, XLAL_EFUNC );
1200 }
1201 if ( !( simulation_level & WEAVE_SIMULATE_MIN_MEM ) && ( statistics_params->completionloop_statistics[1] & WEAVE_STATISTIC_COH2F_DET ) ) {
1202 for ( size_t i = 0; i < ndetectors; ++i ) {
1203 XLAL_CHECK_MAIN( 0 < statistics_params->n2F_det[i] && statistics_params->n2F_det[i] <= nsegments, XLAL_EFAILED, "Invalid number of per-detector F-statistics (%u) for detector '%s'", statistics_params->n2F_det[i], setup.detectors->data[i] );
1204 }
1205 }
1206 }
1207
1208 // Switch timing section
1210
1211 // Completion loop: compute all extra statistics that weren't required in the main loop
1213
1214 // Switch timing section
1216
1217 // Stop timing main search loop, and get total wall and CPU times
1218 double wall_total = 0, cpu_total = 0;
1219 XLAL_CHECK_MAIN( XLALWeaveSearchTimingStop( tim, &wall_total, &cpu_total ) == XLAL_SUCCESS, XLAL_EFUNC );
1220 LogPrintf( LOG_NORMAL, "Finished completion-loop, total time %.1f sec, CPU %.1f%%, peak memory %.1fMB\n", wall_total, 100.0 * cpu_total / wall_total, XLALGetPeakHeapUsageMB() );
1221
1222 ////////// Output search results //////////
1223
1224 if ( search_complete ) {
1225
1226 // Open output file
1227 LogPrintf( LOG_NORMAL, "Opening output file '%s' for writing ...\n", uvar->output_file );
1228 FITSFile *file = XLALFITSFileOpenWrite( uvar->output_file );
1229 XLAL_CHECK_MAIN( file != NULL, XLAL_EFUNC );
1232
1233 // Write number of times output results were restored from a checkpoint
1234 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT4( file, "numckpt", ckpt_output_count, "number of checkpoints" ) == XLAL_SUCCESS, XLAL_EFUNC );
1235
1236 // Write list of detectors
1237 XLAL_CHECK_MAIN( XLALFITSHeaderWriteStringVector( file, "detect", setup.detectors, "setup detectors" ) == XLAL_SUCCESS, XLAL_EFUNC );
1238
1239 // Write number of segments
1240 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT4( file, "nsegment", nsegments, "number of segments" ) == XLAL_SUCCESS, XLAL_EFUNC );
1241
1242 // Write frequency spacing
1243 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, "dfreq", dfreq, "frequency spacing" ) == XLAL_SUCCESS, XLAL_EFUNC );
1244
1245 // Write semicoherent parameter-space bounds
1246 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, "semiparam skyarea [sr]", semi_sky_area, "area of sky parameter space" ) == XLAL_SUCCESS, XLAL_EFUNC );
1247 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, "semiparam minfreq [Hz]", uvar->freq[0], "minimum frequency range" ) == XLAL_SUCCESS, XLAL_EFUNC );
1248 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, "semiparam maxfreq [Hz]", uvar->freq[1], "maximum frequency range" ) == XLAL_SUCCESS, XLAL_EFUNC );
1249 for ( size_t s = 1; s <= ninputspins; ++s ) {
1250 char keyword[64];
1251 char comment[64];
1252 snprintf( keyword, sizeof( keyword ), "semiparam minf%zudot [Hz/s^%zu]", s, s );
1253 snprintf( comment, sizeof( comment ), "minimum %zu-order spindown range", s );
1254 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, keyword, uvarspins[s - 1][0], comment ) == XLAL_SUCCESS, XLAL_EFUNC );
1255 snprintf( keyword, sizeof( keyword ), "semiparam maxf%zudot [Hz/s^%zu]", s, s );
1256 snprintf( comment, sizeof( comment ), "maximum %zu-order spindown range", s );
1257 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, keyword, uvarspins[s - 1][1], comment ) == XLAL_SUCCESS, XLAL_EFUNC );
1258 }
1259
1260 // Write cumulative number of semicoherent templates in each dimension
1261 for ( size_t i = 0; i < ndim; ++i ) {
1262 char keyword[64];
1263 const LatticeTilingStats *semi_stats = XLALLatticeTilingStatistics( tiling[isemi], i );
1264 XLAL_CHECK_MAIN( semi_stats != NULL, XLAL_EFUNC );
1265 XLAL_CHECK_MAIN( semi_stats->name != NULL, XLAL_EFUNC );
1266 XLAL_CHECK_MAIN( semi_stats->total_points > 0, XLAL_EFUNC );
1267 snprintf( keyword, sizeof( keyword ), "nsemitmpl %s", semi_stats->name );
1268 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT8( file, keyword, semi_stats->total_points, "cumulative number of semicoherent templates" ) == XLAL_SUCCESS, XLAL_EFUNC );
1269 }
1270
1271 // Write number of computed coherent results, and number of coherent and semicoherent templates
1272 UINT8 coh_nres = 0, coh_ntmpl = 0, semi_ntmpl = 0;
1273 XLAL_CHECK_MAIN( XLALWeaveCacheQueriesGetCounts( queries, &coh_nres, &coh_ntmpl, &semi_ntmpl ) == XLAL_SUCCESS, XLAL_EFUNC );
1274 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT8( file, "ncohres", coh_nres, "number of computed coherent results" ) == XLAL_SUCCESS, XLAL_EFUNC );
1275 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT8( file, "ncohtpl", coh_ntmpl, "number of coherent templates" ) == XLAL_SUCCESS, XLAL_EFUNC );
1276 XLAL_CHECK_MAIN( XLALFITSHeaderWriteUINT8( file, "nsemitpl", semi_ntmpl, "number of semicoherent templates" ) == XLAL_SUCCESS, XLAL_EFUNC );
1277
1278 // Write peak memory usage
1279 XLAL_CHECK_MAIN( XLALFITSHeaderWriteREAL8( file, "peakmem [MB]", XLALGetPeakHeapUsageMB(), "peak memory usage" ) == XLAL_SUCCESS, XLAL_EFUNC );
1280
1281 // Write timing information
1283
1284 // Write various information from coherent input data
1285 XLAL_CHECK_MAIN( XLALWeaveCohInputWriteInfo( file, nsegments, statistics_params->coh_input ) == XLAL_SUCCESS, XLAL_EFUNC );
1286
1287 // Write various information from caches
1288 XLAL_CHECK_MAIN( XLALWeaveCacheWriteInfo( file, nsegments, coh_cache ) == XLAL_SUCCESS, XLAL_EFUNC );
1289
1290 // Write search results, unless search is being simulated
1291 if ( simulation_level == 0 ) {
1293 }
1294
1295 // Write list of injections
1296 if ( injections != NULL ) {
1298 }
1299
1300 // Write various segment information from coherent input data
1301 if ( uvar->segment_info ) {
1302 XLAL_CHECK_MAIN( XLALWeaveCohInputWriteSegInfo( file, nsegments, statistics_params->coh_input ) == XLAL_SUCCESS, XLAL_EFUNC );
1303 }
1304
1305 // Close output file
1307 LogPrintf( LOG_NORMAL, "Closed output file '%s'\n", uvar->output_file );
1308
1309 }
1310
1311 ////////// Cleanup memory and exit //////////
1312
1313 // Cleanup memory from search timing
1315
1316 // Cleanup memory from output results
1318
1319 // Cleanup memory from semicoherent results
1320 XLALWeaveSemiResultsDestroy( semi_res );
1321
1322 // Cleanup memory from parameter-space iteration
1323 XLALWeaveSearchIteratorDestroy( main_loop_itr );
1324
1325 // Cleanup memory from computing 'stage 0' coherent results
1327 for ( size_t i = 0; i < nsegments; ++i ) {
1328 XLALWeaveCacheDestroy( coh_cache[i] );
1329 }
1330
1331 // Cleanup memory from loading input data
1332 XLALDestroySFTCatalog( sft_catalog );
1333
1334 // Cleanup memory from lattice tilings
1335 for ( size_t i = 0; i < ntiles; ++i ) {
1337 }
1338
1339 // Cleanup memory from setup data
1340 XLALWeaveSetupDataClear( &setup );
1341 XLALFree( setup_detectors_string );
1342
1343 // Cleanup random number generator
1344 XLALDestroyRandomParams( rand_par );
1345
1346 // Cleanup injection parameters vector
1348
1349 // Cleanup memory from user input
1351
1352 // Check for memory leaks
1354
1355 if ( search_complete ) {
1356 LogPrintf( LOG_NORMAL, "Finished successfully!\n" );
1357 } else {
1358 LogPrintf( LOG_NORMAL, "Finished but search not completed!\n" );
1359 }
1360
1361 return EXIT_SUCCESS;
1362
1363}
1364
1365// Local Variables:
1366// c-file-style: "linux"
1367// c-basic-offset: 2
1368// End:
WeaveCache * XLALWeaveCacheCreate(const LatticeTiling *coh_tiling, const BOOLEAN interpolation, const SuperskyTransformData *coh_rssky_transf, const SuperskyTransformData *semi_rssky_transf, WeaveCohInput *coh_input, const UINT4 max_size, const BOOLEAN all_gc)
Create a cache.
Definition: CacheResults.c:616
void XLALWeaveCacheQueriesDestroy(WeaveCacheQueries *queries)
Destroy storage for a series of cache queries.
Definition: CacheResults.c:383
int XLALWeaveCacheQuery(const WeaveCache *cache, WeaveCacheQueries *queries, const UINT4 query_index)
Query a cache for the results nearest to a given coherent point.
Definition: CacheResults.c:449
int XLALWeaveCacheQueriesFinal(WeaveCacheQueries *queries, PulsarDopplerParams *semi_phys, UINT4 *semi_nfreqs)
Finalise a series of cache queries.
Definition: CacheResults.c:515
int XLALWeaveCacheQueriesInit(WeaveCacheQueries *queries, const UINT8 semi_index, const gsl_vector *semi_rssky, const INT4 semi_left, const INT4 semi_right, const UINT4 freq_partition_index)
Initialise a series of cache queries.
Definition: CacheResults.c:404
int XLALWeaveCacheExpire(WeaveCache *cache)
Expire all items in the cache.
Definition: CacheResults.c:833
int XLALWeaveGetCacheMeanMaxSize(REAL4 *cache_mean_max_size, const size_t ncache, WeaveCache *const *cache)
Determine the mean maximum size obtained by caches.
Definition: CacheResults.c:783
void XLALWeaveCacheDestroy(WeaveCache *cache)
Destroy a cache.
Definition: CacheResults.c:766
int XLALWeaveCacheRetrieve(WeaveCache *cache, const WeaveCacheQueries *queries, const UINT4 query_index, const WeaveCohResults **coh_res, UINT8 *coh_index, UINT4 *coh_offset, WeaveSearchTiming *tim)
Retrieve coherent results for a given query, or compute new coherent results if not found.
Definition: CacheResults.c:874
WeaveCacheQueries * XLALWeaveCacheQueriesCreate(const LatticeTiling *semi_tiling, const SuperskyTransformData *semi_rssky_transf, const double dfreq, const UINT4 nqueries, const UINT4 nfreq_partitions)
Create storage for a series of cache queries.
Definition: CacheResults.c:305
int XLALWeaveCacheQueriesGetCounts(const WeaveCacheQueries *queries, UINT8 *coh_nres, UINT8 *coh_ntmpl, UINT8 *semi_ntmpl)
Get number of computed coherent results, and number of coherent and semicoherent templates.
Definition: CacheResults.c:580
int XLALWeaveCacheClear(WeaveCache *cache)
Clear all items in the cache from memory.
Definition: CacheResults.c:852
int XLALWeaveCacheWriteInfo(FITSFile *file, const size_t ncache, WeaveCache *const *cache)
Write various information from caches to a FITS file.
Definition: CacheResults.c:809
Module which caches computed coherent results.
int XLALWeaveSemiResultsComputeMain(WeaveSemiResults *semi_res, WeaveSearchTiming *tim)
Compute all remaining toplist-ranking semicoherent statistics (ie 'mainloop-statistics').
int XLALWeaveCohInputWriteInfo(FITSFile *file, const size_t ncoh_input, WeaveCohInput *const *coh_input)
Write various information from coherent input data to a FITS file.
int XLALWeaveCohInputWriteSegInfo(FITSFile *file, const size_t ncoh_input, WeaveCohInput *const *coh_input)
Write various segment information from coherent input data to a FITS file.
int XLALWeaveSemiResultsInit(WeaveSemiResults **semi_res, const WeaveSimulationLevel simulation_level, const UINT4 ndetectors, const UINT4 nsegments, const UINT8 semi_index, const PulsarDopplerParams *semi_phys, const double dfreq, const UINT4 semi_nfreqs, const WeaveStatisticsParams *statistics_params)
Create and initialise semicoherent results.
void XLALWeaveSemiResultsDestroy(WeaveSemiResults *semi_res)
Destroy final semicoherent results.
WeaveCohInput * XLALWeaveCohInputCreate(const LALStringVector *setup_detectors, const WeaveSimulationLevel simulation_level, const SFTCatalog *sft_catalog, const UINT4 segment_index, const LALSeg *segment, const PulsarDopplerParams *min_phys, const PulsarDopplerParams *max_phys, const double dfreq, const EphemerisData *ephemerides, const LALStringVector *sft_noise_sqrtSX, const LALStringVector *Fstat_assume_sqrtSX, FstatOptionalArgs *Fstat_opt_args, WeaveStatisticsParams *statistics_params, BOOLEAN recalc_stage)
Create coherent input data.
int XLALWeaveSemiResultsComputeSegs(WeaveSemiResults *semi_res, const UINT4 nsegments, const WeaveCohResults **coh_res, const UINT8 *coh_index, const UINT4 *coh_offset, WeaveSearchTiming *tim)
Add a new set of coherent results to the semicoherent results.
Module which computes coherent and semicoherent results.
int XLALFITSHeaderWriteUINT4(FITSFile UNUSED *file, const CHAR UNUSED *key, const UINT4 UNUSED value, const CHAR UNUSED *comment)
Definition: FITSFileIO.c:762
int XLALFITSHeaderWriteStringVector(FITSFile UNUSED *file, const CHAR UNUSED *key, const LALStringVector UNUSED *values, const CHAR UNUSED *comment)
Definition: FITSFileIO.c:1412
FITSFile * XLALFITSFileOpenWrite(const CHAR UNUSED *file_name)
Definition: FITSFileIO.c:263
int XLALFITSFileWriteVCSInfo(FITSFile UNUSED *file, const LALVCSInfoList UNUSED vcs_list)
Definition: FITSFileIO.c:466
int XLALFITSFileWriteUVarCmdLine(FITSFile UNUSED *file)
Definition: FITSFileIO.c:508
FITSFile * XLALFITSFileOpenRead(const CHAR UNUSED *file_name)
Definition: FITSFileIO.c:320
int XLALFITSHeaderWriteUINT8(FITSFile UNUSED *file, const CHAR UNUSED *key, const UINT8 UNUSED value, const CHAR UNUSED *comment)
Definition: FITSFileIO.c:827
void XLALFITSFileClose(FITSFile UNUSED *file)
Definition: FITSFileIO.c:245
int XLALFITSHeaderWriteREAL8(FITSFile UNUSED *file, const CHAR UNUSED *key, const REAL8 UNUSED value, const CHAR UNUSED *comment)
Definition: FITSFileIO.c:1150
int XLALFITSHeaderReadUINT4(FITSFile UNUSED *file, const CHAR UNUSED *key, UINT4 UNUSED *value)
Definition: FITSFileIO.c:797
int j
void LALCheckMemoryLeaks(void)
const LALVCSInfoList lalPulsarVCSInfoList
NULL-terminated list of VCS and build information for LALPulsar and its dependencies
static double double delta
int XLALWeaveOutputResultsAdd(WeaveOutputResults *out, const WeaveSemiResults *semi_res, const UINT4 semi_nfreqs)
Add semicoherent results to output.
int XLALWeaveOutputResultsWrite(FITSFile *file, const WeaveOutputResults *out)
Write output results to a FITS file.
int XLALWeaveOutputResultsReadAppend(FITSFile *file, WeaveOutputResults **out, UINT4 toplist_limit)
Read results from a FITS file and append to new/existing output results.
void XLALWeaveOutputResultsDestroy(WeaveOutputResults *out)
Free output results.
int XLALWeaveOutputResultsCompletionLoop(WeaveOutputResults *out)
Compute all the missing 'completion-loop' statistics for all toplist entries.
WeaveOutputResults * XLALWeaveOutputResultsCreate(const LIGOTimeGPS *ref_time, const size_t nspins, WeaveStatisticsParams *statistics_params, const UINT4 toplist_limit, const BOOLEAN mean2F_hgrm)
Create output results.
Module which handles the output results.
#define STRING(a)
int XLALWeaveSearchIteratorNext(WeaveSearchIterator *itr, BOOLEAN *iteration_complete, BOOLEAN *expire_cache, UINT8 *semi_index, const gsl_vector **semi_rssky, INT4 *semi_left, INT4 *semi_right, UINT4 *repetition_index)
Advance to next state of iterator.
REAL8 XLALWeaveSearchIteratorProgress(const WeaveSearchIterator *itr)
Return progress of iterator as a percentage.
void XLALWeaveSearchIteratorDestroy(WeaveSearchIterator *itr)
Destroy iterator.
int XLALWeaveSearchIteratorSave(const WeaveSearchIterator *itr, FITSFile *file)
Save state of iterator to a FITS file.
WeaveSearchIterator * XLALWeaveMainLoopSearchIteratorCreate(const LatticeTiling *semi_tiling, const UINT4 freq_partitions, const UINT4 f1dot_partitions)
Create iterator over the main loop search parameter space.
REAL8 XLALWeaveSearchIteratorRemainingTime(const WeaveSearchIterator *itr, const REAL8 elapsed_time)
Return estimate of time remaining for iteration to complete, assuming a equal dstribution in computat...
int XLALWeaveSearchIteratorRestore(WeaveSearchIterator *itr, FITSFile *file)
Restore state of iterator from a FITS file.
Module which implements iterators over search parameter spaces.
int XLALWeaveSearchTimingElapsed(WeaveSearchTiming *tim, double *wall_elapsed, double *cpu_elapsed)
Return elapsed wall and CPU times since start of search timing.
Definition: SearchTiming.c:212
int XLALWeaveSearchTimingWriteInfo(FITSFile *file, const WeaveSearchTiming *tim, const WeaveCacheQueries *queries)
Write information from search timing to a FITS file.
Definition: SearchTiming.c:360
int XLALWeaveSearchTimingStop(WeaveSearchTiming *tim, double *wall_total, double *cpu_total)
Stop timing of search.
Definition: SearchTiming.c:242
int XLALWeaveSearchTimingSection(WeaveSearchTiming *tim, const WeaveSearchTimingSection prev_section, const WeaveSearchTimingSection next_section)
Change the search section currently being timed.
Definition: SearchTiming.c:286
void XLALWeaveSearchTimingDestroy(WeaveSearchTiming *tim)
Destroy a search timing structure.
Definition: SearchTiming.c:163
int XLALWeaveSearchTimingStart(WeaveSearchTiming *tim)
Start timing of search.
Definition: SearchTiming.c:175
const char * comment
Definition: SearchTiming.c:94
WeaveSearchTiming * XLALWeaveSearchTimingCreate(const BOOLEAN detailed_timing, const WeaveStatisticsParams *statistics_params)
Create a search timing structure.
Definition: SearchTiming.c:140
Module which collects search timings and builds a timing model.
@ WEAVE_SEARCH_TIMING_COH
Computation of coherent results section.
Definition: SearchTiming.h:46
@ WEAVE_SEARCH_TIMING_OUTPUT
Result output section.
Definition: SearchTiming.h:52
@ WEAVE_SEARCH_TIMING_SEMISEG
Computation of per-segment semicoherent results section.
Definition: SearchTiming.h:48
@ WEAVE_SEARCH_TIMING_ITER
Parameter space iteration section.
Definition: SearchTiming.h:42
@ WEAVE_SEARCH_TIMING_QUERY
Cache queries section.
Definition: SearchTiming.h:44
@ WEAVE_SEARCH_TIMING_CKPT
Checkpointing section.
Definition: SearchTiming.h:54
@ WEAVE_SEARCH_TIMING_OTHER
Unaccounted section.
Definition: SearchTiming.h:58
@ WEAVE_SEARCH_TIMING_CMPL
Completion-loop section.
Definition: SearchTiming.h:56
@ WEAVE_SEARCH_TIMING_SEMI
Computation of semicoherent results section.
Definition: SearchTiming.h:50
void XLALWeaveSetupDataClear(WeaveSetupData *setup)
Free contents of setup data.
Definition: SetupData.c:40
int XLALWeaveSetupDataRead(FITSFile *file, WeaveSetupData *setup)
Read setup data from a FITS file.
Definition: SetupData.c:91
Module which handles the setup data.
int main(int argc, char *argv[])
Definition: Weave.c:44
@ WEAVE_SIMULATE
Simulate search (implicitly with full memory allocation)
Definition: Weave.h:82
@ WEAVE_SIMULATE_MIN_MEM
Simulate search with minimal memory allocation.
Definition: Weave.h:84
enum tagWeaveSimulationLevel WeaveSimulationLevel
Definition: Weave.h:60
int XLALWeaveStatisticsParamsSetDependencyMap(WeaveStatisticsParams *statistics_params, const WeaveStatisticType toplist_stats, const WeaveStatisticType extra_output_stats, const WeaveStatisticType recalc_stats)
Fill StatisticsParams logic for given toplist and extra-output stats.
const char *const WeaveToplistHelpString
User input help string for toplist ranking statistics.
const UserChoices WeaveToplistChoices
User input choices for toplist ranking statistics.
const char *const WeaveStatisticHelpString
User input help string for all supported statistics.
const UserChoices WeaveStatisticChoices
User input choices for all supported statistics.
@ WEAVE_STATISTIC_COH2F_DET
Per segment per-detector F-statistic.
@ WEAVE_STATISTIC_BtSGLtL
(transient-)line robust log10(B_tS/GLtL) statistic
@ WEAVE_STATISTIC_MEAN2F
Multi-detector average (over segments) F-statistic.
@ WEAVE_STATISTIC_BSGLtL
(transient-)line robust log10(B_S/GLtL) statistic
@ WEAVE_STATISTIC_BSGL
Line-robust log10(B_S/GL) statistic.
@ WEAVE_STATISTIC_NONE
PulsarParamsVector * XLALCreatePulsarParamsVector(UINT4 numPulsars)
Create zero-initialized PulsarParamsVector for numPulsars.
PulsarParamsVector * XLALPulsarParamsFromUserInput(const LALStringVector *UserInput, const LIGOTimeGPS *refTimeDef)
Function to determine the PulsarParamsVector input from a user-input defining CW sources.
const char *const InjectionSourcesHelpString
void XLALDestroyPulsarParamsVector(PulsarParamsVector *ppvect)
Destructor for PulsarParamsVector type.
int XLALFITSWritePulsarParamsVector(FITSFile *file, const CHAR *tableName, const PulsarParamsVector *list)
Write a PulsarParamsVector to a FITS file.
const FstatOptionalArgs FstatOptionalArgsDefaults
Global initializer for setting FstatOptionalArgs to default values.
Definition: ComputeFstat.c:90
const UserChoices * XLALFstatMethodChoices(void)
Return pointer to a static array of all (available) FstatMethodType choices.
@ FMETHOD_DEMOD_BEST
Demod: best guess of the fastest available hotloop
Definition: ComputeFstat.h:121
@ FMETHOD_RESAMP_BEST
Resamp: best guess of the fastest available implementation
Definition: ComputeFstat.h:125
@ FMETHOD_RESAMP_CUDA
Resamp: CUDA resampling
Definition: ComputeFstat.h:124
#define LAL_GPS_FORMAT
#define LAL_GPS_PRINT(gps)
struct tagFITSFile FITSFile
Representation of a FITS file.
Definition: FITSFileIO.h:54
#define LAL_PI_2
#define LAL_PI
#define LAL_TWOPI
unsigned char BOOLEAN
#define XLAL_NUM_ELEM(x)
#define XLAL_INIT_MEM(x)
uint64_t UINT8
double REAL8
#define XLAL_INIT_DECL(var,...)
char CHAR
uint32_t UINT4
int32_t INT4
float REAL4
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
#define LAL_UINT8_FORMAT
int XLALSetTiledLatticeDimensionsFromTiling(LatticeTiling *tiling, const LatticeTiling *ref_tiling)
Set the tiled (i.e.
REAL8 XLALLatticeTilingStepSize(const LatticeTiling *tiling, const size_t dim)
Return the step size of the lattice tiling in a given dimension, or 0 for non-tiled dimensions.
int XLALSetTilingLatticeAndMetric(LatticeTiling *tiling, const TilingLattice lattice, const gsl_matrix *metric, const double max_mismatch)
Set the tiling lattice, parameter-space metric, and maximum prescribed mismatch.
size_t XLALTiledLatticeTilingDimensions(const LatticeTiling *tiling)
Return the number of tiled dimensions of the lattice tiling.
int XLALPerformLatticeTilingCallbacks(const LatticeTiling *tiling)
Perform all registered lattice tiling callbacks.
void XLALDestroyLatticeTiling(LatticeTiling *tiling)
Destroy a lattice tiling.
LatticeTiling * XLALCreateLatticeTiling(const size_t ndim)
Create a new lattice tiling.
const LatticeTilingStats * XLALLatticeTilingStatistics(const LatticeTiling *tiling, const size_t dim)
Return statistics related to the number/value of lattice tiling points in a dimension.
const UserChoices TilingLatticeChoices
Static array of all :tagTilingLattice choices, for use by the UserInput module parsing routines.
int XLALSetLatticeTilingRandomOriginOffsets(LatticeTiling *tiling, RandomParams *rng)
Offset the physical parameter-space origin of the lattice tiling by a random fraction of the lattice ...
int XLALRandomLatticeTilingPoints(const LatticeTiling *tiling, const double scale, RandomParams *rng, gsl_matrix *random_points)
Generate random points within the parameter space of the lattice tiling.
@ TILING_LATTICE_ANSTAR
An-star ( ) lattice.
Definition: LatticeTiling.h:65
BSGLSetup * XLALCreateBSGLSetup(const UINT4 numDetectors, const REAL4 Fstar0sc, const REAL4 oLGX[PULSAR_MAX_DETECTORS], const BOOLEAN useLogCorrection, const UINT4 numSegments)
int XLALParseLinePriors(REAL4 oLGX[PULSAR_MAX_DETECTORS], const LALStringVector *oLGX_string)
Parse string-vectors (typically input by user) of N per-detector line-to-Gaussian prior ratios to a ...
REAL8 XLALGetPeakHeapUsageMB(void)
void void LogPrintfVerbatim(LogLevel_t, const char *format,...) _LAL_GCC_PRINTF_FORMAT_(2
void LogPrintf(LogLevel_t, const char *format,...) _LAL_GCC_PRINTF_FORMAT_(2
LOG_NORMAL
#define PULSAR_MAX_DETECTORS
maximal number of detectors we can handle (for static arrays of detector quantities)
RandomParams * XLALCreateRandomParams(INT4 seed)
REAL4 XLALUniformDeviate(RandomParams *params)
void XLALDestroyRandomParams(RandomParams *params)
void XLALDestroySFTCatalog(SFTCatalog *catalog)
Free an 'SFT-catalogue'.
Definition: SFTcatalog.c:329
void XLALDestroyMultiSFTCatalogView(MultiSFTCatalogView *multiView)
Destroys a MultiSFTCatalogView, without freeing the original catalog that the 'view' was referring to...
Definition: SFTcatalog.c:496
MultiLIGOTimeGPSVector * XLALMakeMultiTimestamps(LIGOTimeGPS tStart, REAL8 Tspan, REAL8 Tsft, REAL8 Toverlap, UINT4 numDet)
Same as XLALMakeTimestamps() just for several detectors, additionally specify the number of detectors...
CHAR * XLALGetChannelPrefix(const CHAR *name)
Find the valid CW detector prefix.
Definition: SFTnaming.c:203
SFTCatalog * XLALMultiAddToFakeSFTCatalog(SFTCatalog *catalog, const LALStringVector *detectors, const MultiLIGOTimeGPSVector *timestamps)
Multi-detector and multi-timestamp wrapper of XLALAddToFakeSFTCatalog().
Definition: SFTcatalog.c:749
int XLALCheckCRCSFTCatalog(BOOLEAN *crc_check, SFTCatalog *catalog)
This function reads in the SFTs in the catalog and validates their CRC64 checksums.
Definition: SFTcatalog.c:524
MultiSFTCatalogView * XLALGetMultiSFTCatalogView(const SFTCatalog *catalog)
Return a MultiSFTCatalogView generated from an input SFTCatalog.
Definition: SFTcatalog.c:380
LALStringVector * XLALListIFOsInCatalog(const SFTCatalog *catalog)
Return a sorted string vector listing the unique IFOs in the given catalog.
Definition: SFTcatalog.c:562
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
void XLALDestroyMultiTimestamps(MultiLIGOTimeGPSVector *multiTS)
Destroy a MultiLIGOTimeGPSVector timestamps vector.
MultiLIGOTimeGPSVector * XLALReadMultiTimestampsFiles(const LALStringVector *fnames)
backwards compatible wrapper to XLALReadMultiTimestampsFilesConstrained() without GPS-time constraint...
const UserChoices SSBprecisionChoices
Static array of all SSBprecision choices, for use by the UserInput module parsing routines.
Definition: SSBtimes.c:41
@ SSBPREC_LAST
end marker
Definition: SSBtimes.h:50
int XLALSegListRange(const LALSegList *seglist, LIGOTimeGPS *start, LIGOTimeGPS *end)
void XLALDestroyStringVector(LALStringVector *vect)
char * XLALConcatStringVector(const LALStringVector *strings, const char *sep)
INT4 XLALFindStringInVector(const char *needle, const LALStringVector *haystack)
LALStringVector * XLALCopyStringVector(const LALStringVector *vect)
int XLALSetSuperskyPhysicalSkyBounds(LatticeTiling *tiling, gsl_matrix *rssky_metric, SuperskyTransformData *rssky_transf, const double alpha1, const double alpha2, const double delta1, const double delta2)
Set parameter-space bounds on the physical sky position for a lattice tiling using the reduced super...
int XLALSetSuperskyEqualAreaSkyBounds(LatticeTiling *tiling, const gsl_matrix *rssky_metric, const double max_mismatch, const UINT4 patch_count, const UINT4 patch_index)
Divide the reduced supersky parameter space into patch_count equal-area patches, and set parameter-sp...
int XLALSuperskyMetricsDimensions(const SuperskyMetrics *metrics, size_t *spindowns)
Return dimensions of the supersky metrics.
int XLALEqualizeReducedSuperskyMetricsFreqSpacing(SuperskyMetrics *metrics, const double coh_max_mismatch, const double semi_max_mismatch)
Project and rescale the reduced supersky metrics in the frequency dimension, such that all reduced su...
int XLALScaleSuperskyMetricsFiducialFreq(SuperskyMetrics *metrics, const double new_fiducial_freq)
Scale all supersky metrics and their coordinate transform data to a new fiducial frequency.
int XLALConvertSuperskyToPhysicalPoint(PulsarDopplerParams *out_phys, const gsl_vector *in_rssky, const gsl_vector *ref_rssky, const SuperskyTransformData *rssky_transf)
Convert a point from supersky to physical coordinates.
int XLALSetSuperskyRangeBounds(LatticeTiling *tiling, const gsl_vector *min_rssky, const gsl_vector *max_rssky)
Set parameter-space bounds on an entire lattice tiling given minimum and maximum ranges in reduced su...
int XLALRegisterSuperskyLatticePhysicalRangeCallback(LatticeTiling *tiling, const SuperskyTransformData *rssky_transf, const PulsarDopplerParams **min_phys, const PulsarDopplerParams **max_phys)
Register a lattice tiling callback function which computes the physical range covered by a reduced su...
int XLALSetSuperskyPhysicalSpinBoundPadding(LatticeTiling *tiling, const SuperskyTransformData *rssky_transf, const size_t s, const bool padding)
Set parameter-space bound padding on the physical frequency/spindowns for a lattice tiling using the...
int XLALSetSuperskyPhysicalSpinBound(LatticeTiling *tiling, const SuperskyTransformData *rssky_transf, const size_t s, const double bound1, const double bound2)
Set parameter-space bounds on the physical frequency/spindowns for a lattice tiling using the reduce...
int XLALRegisterSuperskyLatticeSuperskyRangeCallback(LatticeTiling *tiling, const SuperskyTransformData *rssky_transf, const SuperskyTransformData *rssky2_transf, const gsl_vector **min_rssky2, const gsl_vector **max_rssky2)
Register a lattice tiling callback function which computes the range covered by a reduced supersky la...
const char * lalUserVarHelpOptionSubsection
const char * lalUserVarHelpBrief
#define UVAR_STR3AND(n1, n2, n3)
#define UVAR_STR2AND(n1, n2)
int XLALUserVarReadAllInput(BOOLEAN *should_exit, int argc, char *argv[], const LALVCSInfoList vcs_list)
void XLALDestroyUserVars(void)
#define UVAR_ALLSET2(n1, n2)
#define UVAR_SET4(n1, n2, n3, n4)
#define XLALRegisterUvarMember(name, type, option, category,...)
#define UVAR_STR(n)
void XLALUserVarCheck(BOOLEAN *should_exit, const int assertion, const CHAR *fmt,...) _LAL_GCC_PRINTF_FORMAT_(3
#define UVAR_STR2OR(n1, n2)
#define UVAR_ALLSET3(n1, n2, n3)
#define UVAR_SET(n)
#define XLALRegisterUvarAuxDataMember(name, type, cdata, option, category,...)
#define UVAR_SET2(n1, n2)
REAL8 REAL8Range[2]
#define xlalErrno
#define XLAL_CHECK(assertion,...)
#define XLAL_CHECK_MAIN(assertion,...)
#define XLAL_TRY(statement, errnum)
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_ENOENT
XLAL_EFUNC
XLAL_EERR
XLAL_EIO
XLAL_EINVAL
XLAL_EFAILED
REAL8 XLALGPSDiff(const LIGOTimeGPS *t1, const LIGOTimeGPS *t0)
injections
output_file
out
LALDetector detectors[LAL_NUM_DETECTORS]
double alpha
CHAR name[LALNameLength]
Struct of optional 'advanced level' and (potentially method-specific) arguments to be passed to the ...
Definition: ComputeFstat.h:137
UINT4 randSeed
Random-number seed value used in case of fake Gaussian noise generation (injectSqrtSX)
Definition: ComputeFstat.h:138
PulsarParamsVector * injectSources
Vector of parameters of CW signals to simulate and inject.
Definition: ComputeFstat.h:143
UINT4 runningMedianWindow
If SFT noise weights are calculated from the SFTs, the running median window length to use.
Definition: ComputeFstat.h:141
FstatInput * prevInput
An FstatInput structure from a previous call to XLALCreateFstatInput(); may contain common workspace ...
Definition: ComputeFstat.h:146
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
SSBprecision SSBprec
Barycentric transformation precision.
Definition: ComputeFstat.h:139
REAL8 deltaT
'length' of each timestamp (e.g.
Definition: SFTfileIO.h:194
Statistics related to the number/value of lattice tiling points in a dimension.
UINT8 total_points
Total number of points up to this dimension.
const char * name
Name of parameter-space dimension.
A collection of (multi-IFO) LIGOTimeGPSVector time-stamps vectors.
Definition: SFTfileIO.h:198
LIGOTimeGPSVector ** data
timestamps vector for each ifo
Definition: SFTfileIO.h:203
A multi-SFT-catalogue "view": a multi-IFO vector of SFT-catalogs.
Definition: SFTfileIO.h:255
SFTCatalog * data
array of SFT-catalog pointers
Definition: SFTfileIO.h:260
UINT4 length
number of detectors
Definition: SFTfileIO.h:259
REAL8 aCross
Signal amplitude (cross)
REAL8 psi
polarization angle psi
REAL8 aPlus
Signal amplitude (plus)
REAL8 phi0
initial signal-phase (at some reference time)
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
Type defining the parameters of a pulsar-source of CW Gravitational waves.
PulsarAmplitudeParams Amp
'Amplitude-parameters': h0, cosi, phi0, psi
CHAR name[LALNameLength]
'name' for this sources, can be an empty string
PulsarDopplerParams Doppler
'Phase-evolution parameters': {skypos, fkdot, orbital params }
Straightforward vector type of N PulsarParams structs.
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
SFTtype header
SFT-header info.
Definition: SFTfileIO.h:228
Input data segment info.
UserInput_t uvar_struct
Definition: sw_inj_frames.c:93