static int MismatchSquareTest( const char *lattice_name, const double freqband, const double f1dotband, const double f2dotband, const UINT8 total_ref, const double mism_hist_ref[MISM_HIST_BINS] ) { // Create lattice tiling LatticeTiling *tiling = XLALCreateLatticeTiling(3); XLAL_CHECK(tiling != NULL, XLAL_EFUNC); // Add bounds const double fndot[3] = {100, 0, 0}; const double fndotband[3] = {freqband, f1dotband, f2dotband}; for (size_t i = 0; i < 3; ++i) { printf("Bounds: f%zudot=%0.3g, f%zudotband=%0.3g\n", i, fndot[i], i, fndotband[i]); XLAL_CHECK(XLALSetLatticeTilingConstantBound(tiling, i, fndot[i], fndot[i] + fndotband[i]) == XLAL_SUCCESS, XLAL_EFUNC); } // Set metric to the spindown metric const double max_mismatch = 0.3; gsl_matrix *GAMAT(metric, 3, 3); for (size_t i = 0; i < metric->size1; ++i) { for (size_t j = i; j < metric->size2; ++j) { const double Tspan = 432000; const double metric_i_j_num = 4.0 * LAL_PI * LAL_PI * pow(Tspan, i + j + 2) * (i + 1) * (j + 1); const double metric_i_j_denom = LAL_FACT[i + 1] * LAL_FACT[j + 1] * (i + 2) * (j + 2) * (i + j + 3); gsl_matrix_set(metric, i, j, metric_i_j_num / metric_i_j_denom); gsl_matrix_set(metric, j, i, gsl_matrix_get(metric, i, j)); } } printf("Lattice type: %s\n", lattice_name); XLAL_CHECK(XLALSetTilingLatticeAndMetric(tiling, lattice_name, metric, max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch test XLAL_CHECK(MismatchTest(tiling, metric, max_mismatch, total_ref, mism_hist_ref) == XLAL_SUCCESS, XLAL_EFUNC); return XLAL_SUCCESS; }
static int MismatchAgeBrakeTest( const char *lattice_name, const double freq, const double freqband, const UINT8 total_ref, const double mism_hist_ref[MISM_HIST_BINS] ) { // Create lattice tiling LatticeTiling *tiling = XLALCreateLatticeTiling(3); XLAL_CHECK(tiling != NULL, XLAL_EFUNC); // Add bounds printf("Bounds: freq=%0.3g, freqband=%0.3g\n", freq, freqband); XLAL_CHECK(XLALSetLatticeTilingConstantBound(tiling, 0, freq, freq + freqband) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetLatticeTilingF1DotAgeBrakingBound(tiling, 0, 1, 1e11, 2, 5) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetLatticeTilingF2DotBrakingBound(tiling, 0, 1, 2, 2, 5) == XLAL_SUCCESS, XLAL_EFUNC); // Set metric to the spindown metric const double max_mismatch = 0.3; gsl_matrix *GAMAT(metric, 3, 3); for (size_t i = 0; i < metric->size1; ++i) { for (size_t j = i; j < metric->size2; ++j) { const double Tspan = 1036800; const double metric_i_j_num = 4.0 * LAL_PI * LAL_PI * pow(Tspan, i + j + 2) * (i + 1) * (j + 1); const double metric_i_j_denom = LAL_FACT[i + 1] * LAL_FACT[j + 1] * (i + 2) * (j + 2) * (i + j + 3); gsl_matrix_set(metric, i, j, metric_i_j_num / metric_i_j_denom); gsl_matrix_set(metric, j, i, gsl_matrix_get(metric, i, j)); } } printf("Lattice type: %s\n", lattice_name); XLAL_CHECK(XLALSetTilingLatticeAndMetric(tiling, lattice_name, metric, max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch test XLAL_CHECK(MismatchTest(tiling, metric, max_mismatch, total_ref, mism_hist_ref) == XLAL_SUCCESS, XLAL_EFUNC); return XLAL_SUCCESS; }
static int BasicTest( size_t n, const int bound_on_0, const int bound_on_1, const int bound_on_2, const int bound_on_3, const char *lattice_name, const UINT8 total_ref_0, const UINT8 total_ref_1, const UINT8 total_ref_2, const UINT8 total_ref_3 ) { const int bound_on[4] = {bound_on_0, bound_on_1, bound_on_2, bound_on_3}; const UINT8 total_ref[4] = {total_ref_0, total_ref_1, total_ref_2, total_ref_3}; // Create lattice tiling LatticeTiling *tiling = XLALCreateLatticeTiling(n); XLAL_CHECK(tiling != NULL, XLAL_EFUNC); // Add bounds for (size_t i = 0; i < n; ++i) { XLAL_CHECK(bound_on[i] == 0 || bound_on[i] == 1, XLAL_EFAILED); XLAL_CHECK(XLALSetLatticeTilingConstantBound(tiling, i, 0.0, bound_on[i] * pow(100.0, 1.0/n)) == XLAL_SUCCESS, XLAL_EFUNC); } // Set metric to the Lehmer matrix const double max_mismatch = 0.3; { gsl_matrix *GAMAT(metric, n, n); for (size_t i = 0; i < n; ++i) { for (size_t j = 0; j < n; ++j) { const double ii = i+1, jj = j+1; gsl_matrix_set(metric, i, j, jj >= ii ? ii/jj : jj/ii); } } XLAL_CHECK(XLALSetTilingLatticeAndMetric(tiling, lattice_name, metric, max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); GFMAT(metric); printf("Number of (tiled) dimensions: %zu (%zu)\n", XLALTotalLatticeTilingDimensions(tiling), XLALTiledLatticeTilingDimensions(tiling)); printf(" Bounds: %i %i %i %i\n", bound_on_0, bound_on_1, bound_on_2, bound_on_3); printf(" Lattice type: %s\n", lattice_name); } // Create lattice tiling locator LatticeTilingLocator *loc = XLALCreateLatticeTilingLocator(tiling); XLAL_CHECK(loc != NULL, XLAL_EFUNC); if (lalDebugLevel & LALINFOBIT) { printf(" Index trie:\n"); XLAL_CHECK(XLALPrintLatticeTilingIndexTrie(loc, stdout) == XLAL_SUCCESS, XLAL_EFUNC); } for (size_t i = 0; i < n; ++i) { // Create lattice tiling iterator and locator over 'i+1' dimensions LatticeTilingIterator *itr = XLALCreateLatticeTilingIterator(tiling, i+1); XLAL_CHECK(itr != NULL, XLAL_EFUNC); // Count number of points const UINT8 total = XLALTotalLatticeTilingPoints(itr); printf("Number of lattice points in %zu dimensions: %" LAL_UINT8_FORMAT "\n", i+1, total); XLAL_CHECK(imaxabs(total - total_ref[i]) <= 1, XLAL_EFUNC, "ERROR: |total - total_ref[%zu]| = |%" LAL_UINT8_FORMAT " - %" LAL_UINT8_FORMAT "| > 1", i, total, total_ref[i]); for (UINT8 k = 0; XLALNextLatticeTilingPoint(itr, NULL) > 0; ++k) { const UINT8 itr_index = XLALCurrentLatticeTilingIndex(itr); XLAL_CHECK(k == itr_index, XLAL_EFUNC, "ERROR: k = %" LAL_UINT8_FORMAT " != %" LAL_UINT8_FORMAT " = itr_index", k, itr_index); } XLAL_CHECK(XLALResetLatticeTilingIterator(itr) == XLAL_SUCCESS, XLAL_EFUNC); // Check tiling statistics printf(" Check tiling statistics ..."); for (size_t j = 0; j < n; ++j) { const LatticeTilingStats *stats = XLALLatticeTilingStatistics(tiling, j); XLAL_CHECK(stats != NULL, XLAL_EFUNC); XLAL_CHECK(imaxabs(stats->total_points - total_ref[j]) <= 1, XLAL_EFAILED, "\n " "ERROR: |total - total_ref[%zu]| = |%" LAL_UINT8_FORMAT " - %" LAL_UINT8_FORMAT "| > 1", j, stats->total_points, total_ref[j]); XLAL_CHECK(stats->min_points <= stats->avg_points, XLAL_EFAILED, "\n " "ERROR: min_points = %" LAL_INT4_FORMAT " > %g = avg_points", stats->min_points, stats->avg_points); XLAL_CHECK(stats->max_points >= stats->avg_points, XLAL_EFAILED, "\n " "ERROR: max_points = %" LAL_INT4_FORMAT " < %g = avg_points", stats->max_points, stats->avg_points); } printf(" done\n"); // Get all points gsl_matrix *GAMAT(points, n, total); XLAL_CHECK(XLALNextLatticeTilingPoints(itr, &points) == (int)total, XLAL_EFUNC); XLAL_CHECK(XLALNextLatticeTilingPoint(itr, NULL) == 0, XLAL_EFUNC); // Get nearest points to each template, check for consistency printf(" Testing XLALNearestLatticeTiling{Point|Block}() ..."); gsl_vector *GAVEC(nearest, n); UINT8Vector *nearest_indexes = XLALCreateUINT8Vector(n); XLAL_CHECK(nearest_indexes != NULL, XLAL_ENOMEM); for (UINT8 k = 0; k < total; ++k) { gsl_vector_const_view point_view = gsl_matrix_const_column(points, k); const gsl_vector *point = &point_view.vector; XLAL_CHECK(XLALNearestLatticeTilingPoint(loc, point, nearest, nearest_indexes) == XLAL_SUCCESS, XLAL_EFUNC); gsl_vector_sub(nearest, point); double err = gsl_blas_dasum(nearest) / n; XLAL_CHECK(err < 1e-6, XLAL_EFAILED, "\n " "ERROR: err = %e < 1e-6", err); XLAL_CHECK(nearest_indexes->data[i] == k, XLAL_EFAILED, "\n " "ERROR: nearest_indexes[%zu] = %" LAL_UINT8_FORMAT " != %" LAL_UINT8_FORMAT "\n", i, nearest_indexes->data[i], k); if (0 < i) { const LatticeTilingStats *stats = XLALLatticeTilingStatistics(tiling, i); UINT8 nearest_index = 0; UINT4 nearest_left = 0, nearest_right = 0; XLAL_CHECK(XLALNearestLatticeTilingBlock(loc, point, i, nearest, &nearest_index, &nearest_left, &nearest_right) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(nearest_index == nearest_indexes->data[i-1], XLAL_EFAILED, "\n " "ERROR: nearest_index = %" LAL_UINT8_FORMAT " != %" LAL_UINT8_FORMAT "\n", nearest_index, nearest_indexes->data[i-1]); UINT4 nearest_len = 1 + nearest_left + nearest_right; XLAL_CHECK(nearest_len <= stats->max_points, XLAL_EFAILED, "\n " "ERROR: nearest_len = %i > %i = stats[%zu]->max_points\n", nearest_len, stats->max_points, i); } if (i+1 < n) { const LatticeTilingStats *stats = XLALLatticeTilingStatistics(tiling, i+1); UINT8 nearest_index = 0; UINT4 nearest_left = 0, nearest_right = 0; XLAL_CHECK(XLALNearestLatticeTilingBlock(loc, point, i+1, nearest, &nearest_index, &nearest_left, &nearest_right) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(nearest_index == nearest_indexes->data[i], XLAL_EFAILED, "\n " "ERROR: nearest_index = %" LAL_UINT8_FORMAT " != %" LAL_UINT8_FORMAT "\n", nearest_index, nearest_indexes->data[i]); UINT4 nearest_len = 1 + nearest_left + nearest_right; XLAL_CHECK(nearest_len <= stats->max_points, XLAL_EFAILED, "\n " "ERROR: nearest_len = %i > %i = stats[%zu]->max_points\n", nearest_len, stats->max_points, i+1); } } printf(" done\n"); // Cleanup XLALDestroyLatticeTilingIterator(itr); GFMAT(points); GFVEC(nearest); XLALDestroyUINT8Vector(nearest_indexes); } for (size_t i = 0; i < n; ++i) { // Create alternating lattice tiling iterator over 'i+1' dimensions LatticeTilingIterator *itr_alt = XLALCreateLatticeTilingIterator(tiling, i+1); XLAL_CHECK(itr_alt != NULL, XLAL_EFUNC); XLAL_CHECK(XLALSetLatticeTilingAlternatingIterator(itr_alt, true) == XLAL_SUCCESS, XLAL_EFUNC); // Count number of points, check for consistency with non-alternating count UINT8 total = 0; while (XLALNextLatticeTilingPoint(itr_alt, NULL) > 0) ++total; XLAL_CHECK(imaxabs(total - total_ref[i]) <= 1, XLAL_EFUNC, "ERROR: alternating |total - total_ref[%zu]| = |%" LAL_UINT8_FORMAT " - %" LAL_UINT8_FORMAT "| > 1", i, total, total_ref[i]); // Cleanup XLALDestroyLatticeTilingIterator(itr_alt); } // Cleanup XLALDestroyLatticeTiling(tiling); XLALDestroyLatticeTilingLocator(loc); LALCheckMemoryLeaks(); printf("\n"); fflush(stdout); return XLAL_SUCCESS; }
static int SuperskyTest( const double T, const double max_mismatch, const char *lattice_name, const UINT8 patch_count, const double freq, const double freqband, const UINT8 total_ref, const double mism_hist_ref[MISM_HIST_BINS] ) { // Create lattice tiling LatticeTiling *tiling = XLALCreateLatticeTiling(3); XLAL_CHECK(tiling != NULL, XLAL_EFUNC); // Compute reduced supersky metric const double Tspan = T * 86400; LIGOTimeGPS ref_time; XLALGPSSetREAL8(&ref_time, 900100100); LALSegList segments; { XLAL_CHECK(XLALSegListInit(&segments) == XLAL_SUCCESS, XLAL_EFUNC); LALSeg segment; LIGOTimeGPS start_time = ref_time, end_time = ref_time; XLALGPSAdd(&start_time, -0.5 * Tspan); XLALGPSAdd(&end_time, 0.5 * Tspan); XLAL_CHECK(XLALSegSet(&segment, &start_time, &end_time, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSegListAppend(&segments, &segment) == XLAL_SUCCESS, XLAL_EFUNC); } MultiLALDetector detectors = { .length = 1, .sites = { lalCachedDetectors[LAL_LLO_4K_DETECTOR] } }; EphemerisData *edat = XLALInitBarycenter(TEST_DATA_DIR "earth00-19-DE405.dat.gz", TEST_DATA_DIR "sun00-19-DE405.dat.gz"); XLAL_CHECK(edat != NULL, XLAL_EFUNC); SuperskyMetrics *metrics = XLALComputeSuperskyMetrics(0, &ref_time, &segments, freq, &detectors, NULL, DETMOTION_SPIN | DETMOTION_PTOLEORBIT, edat); XLAL_CHECK(metrics != NULL, XLAL_EFUNC); gsl_matrix *rssky_metric = metrics->semi_rssky_metric, *rssky_transf = metrics->semi_rssky_transf; metrics->semi_rssky_metric = metrics->semi_rssky_transf = NULL; XLALDestroySuperskyMetrics(metrics); XLALSegListClear(&segments); XLALDestroyEphemerisData(edat); // Add bounds printf("Bounds: supersky, sky patch 0/%" LAL_UINT8_FORMAT ", freq=%0.3g, freqband=%0.3g\n", patch_count, freq, freqband); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSkyPatch(tiling, rssky_metric, rssky_transf, patch_count, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSpinBound(tiling, rssky_transf, 0, freq, freq + freqband) == XLAL_SUCCESS, XLAL_EFUNC); GFMAT(rssky_transf); // Set metric printf("Lattice type: %s\n", lattice_name); XLAL_CHECK(XLALSetTilingLatticeAndMetric(tiling, lattice_name, rssky_metric, max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch test XLAL_CHECK(MismatchTest(tiling, rssky_metric, max_mismatch, total_ref, mism_hist_ref) == XLAL_SUCCESS, XLAL_EFUNC); return XLAL_SUCCESS; } static int MultiSegSuperskyTest(void) { printf("Performing multiple-segment tests ...\n"); // Compute reduced supersky metrics const double Tspan = 86400; LIGOTimeGPS ref_time; XLALGPSSetREAL8(&ref_time, 900100100); LALSegList segments; { XLAL_CHECK(XLALSegListInit(&segments) == XLAL_SUCCESS, XLAL_EFUNC); LALSeg segment; { LIGOTimeGPS start_time = ref_time, end_time = ref_time; XLALGPSAdd(&start_time, -4 * Tspan); XLALGPSAdd(&end_time, -3 * Tspan); XLAL_CHECK(XLALSegSet(&segment, &start_time, &end_time, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSegListAppend(&segments, &segment) == XLAL_SUCCESS, XLAL_EFUNC); } { LIGOTimeGPS start_time = ref_time, end_time = ref_time; XLALGPSAdd(&start_time, -0.5 * Tspan); XLALGPSAdd(&end_time, 0.5 * Tspan); XLAL_CHECK(XLALSegSet(&segment, &start_time, &end_time, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSegListAppend(&segments, &segment) == XLAL_SUCCESS, XLAL_EFUNC); } { LIGOTimeGPS start_time = ref_time, end_time = ref_time; XLALGPSAdd(&start_time, 3.5 * Tspan); XLALGPSAdd(&end_time, 4.5 * Tspan); XLAL_CHECK(XLALSegSet(&segment, &start_time, &end_time, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSegListAppend(&segments, &segment) == XLAL_SUCCESS, XLAL_EFUNC); } } MultiLALDetector detectors = { .length = 1, .sites = { lalCachedDetectors[LAL_LLO_4K_DETECTOR] } }; EphemerisData *edat = XLALInitBarycenter(TEST_DATA_DIR "earth00-19-DE405.dat.gz", TEST_DATA_DIR "sun00-19-DE405.dat.gz"); XLAL_CHECK(edat != NULL, XLAL_EFUNC); SuperskyMetrics *metrics = XLALComputeSuperskyMetrics(1, &ref_time, &segments, 50, &detectors, NULL, DETMOTION_SPIN | DETMOTION_PTOLEORBIT, edat); XLAL_CHECK(metrics != NULL, XLAL_EFUNC); // Project and rescale semicoherent metric to give equal frequency spacings const double coh_max_mismatch = 0.2, semi_max_mismatch = 0.4; XLAL_CHECK(XLALEqualizeReducedSuperskyMetricsFreqSpacing(metrics, coh_max_mismatch, semi_max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); // Create lattice tilings LatticeTiling *coh_tiling[metrics->num_segments]; for (size_t n = 0; n < metrics->num_segments; ++n) { coh_tiling[n] = XLALCreateLatticeTiling(4); XLAL_CHECK(coh_tiling[n] != NULL, XLAL_EFUNC); } LatticeTiling *semi_tiling = XLALCreateLatticeTiling(4); XLAL_CHECK(semi_tiling != NULL, XLAL_EFUNC); // Add bounds for (size_t n = 0; n < metrics->num_segments; ++n) { XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSkyPatch(coh_tiling[n], metrics->coh_rssky_metric[n], metrics->coh_rssky_transf[n], 1, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSpinBound(coh_tiling[n], metrics->coh_rssky_transf[n], 0, 50, 50 + 1e-4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSpinBound(coh_tiling[n], metrics->coh_rssky_transf[n], 1, 0, -5e-10) == XLAL_SUCCESS, XLAL_EFUNC); } XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSkyPatch(semi_tiling, metrics->semi_rssky_metric, metrics->semi_rssky_transf, 1, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSpinBound(semi_tiling, metrics->semi_rssky_transf, 0, 50, 50 + 1e-4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK(XLALSetSuperskyLatticeTilingPhysicalSpinBound(semi_tiling, metrics->semi_rssky_transf, 1, 0, -5e-10) == XLAL_SUCCESS, XLAL_EFUNC); // Set metric for (size_t n = 0; n < metrics->num_segments; ++n) { XLAL_CHECK(XLALSetTilingLatticeAndMetric(coh_tiling[n], "Ans", metrics->coh_rssky_metric[n], coh_max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); } XLAL_CHECK(XLALSetTilingLatticeAndMetric(semi_tiling, "Ans", metrics->semi_rssky_metric, semi_max_mismatch) == XLAL_SUCCESS, XLAL_EFUNC); // Check lattice step sizes in frequency const size_t ifreq = 3; const double semi_dfreq = XLALLatticeTilingStepSizes(semi_tiling, ifreq); for (size_t n = 0; n < metrics->num_segments; ++n) { const double coh_dfreq = XLALLatticeTilingStepSizes(coh_tiling[n], ifreq); const double tol = 1e-8; XLAL_CHECK(fabs(coh_dfreq - semi_dfreq) < tol * semi_dfreq, XLAL_EFAILED, " ERROR: semi_dfreq=%0.15e, coh_dfreq[%zu]=%0.15e, |coh_dfreq - semi_dfreq| >= %g * semi_dfreq", semi_dfreq, n, coh_dfreq, tol); } // Check computation of spindown range for coherent tilings for (size_t n = 0; n < metrics->num_segments; ++n) { PulsarSpinRange spin_range; XLAL_CHECK(XLALSuperskyLatticePulsarSpinRange(&spin_range, coh_tiling[n], metrics->coh_rssky_transf[n]) == XLAL_SUCCESS, XLAL_EFUNC); } // Cleanup for (size_t n = 0; n < metrics->num_segments; ++n) { XLALDestroyLatticeTiling(coh_tiling[n]); } XLALDestroyLatticeTiling(semi_tiling); XLALDestroySuperskyMetrics(metrics); XLALSegListClear(&segments); XLALDestroyEphemerisData(edat); LALCheckMemoryLeaks(); printf("\n"); fflush(stdout); return XLAL_SUCCESS; } int main(void) { // Perform basic tests XLAL_CHECK_MAIN(BasicTest(1, 0, 0, 0, 0, "Zn" , 1, 1, 1, 1) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(1, 1, 1, 1, 1, "Ans", 93, 0, 0, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(1, 1, 1, 1, 1, "Zn" , 93, 0, 0, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(2, 0, 0, 0, 0, "Ans", 1, 1, 1, 1) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(2, 1, 1, 1, 1, "Ans", 12, 144, 0, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(2, 1, 1, 1, 1, "Zn" , 13, 190, 0, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(3, 0, 0, 0, 0, "Zn" , 1, 1, 1, 1) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(3, 1, 1, 1, 1, "Ans", 8, 46, 332, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(3, 1, 1, 1, 1, "Zn" , 8, 60, 583, 0) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 0, 0, 0, "Ans", 1, 1, 1, 1) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 0, 0, 1, "Ans", 1, 1, 1, 4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 0, 1, 0, "Ans", 1, 1, 4, 4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 0, 1, 1, "Ans", 1, 1, 4, 20) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 1, 0, 0, "Ans", 1, 4, 4, 4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 1, 0, 1, "Ans", 1, 5, 5, 25) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 1, 1, 0, "Ans", 1, 5, 24, 24) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 0, 1, 1, 1, "Ans", 1, 5, 20, 115) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 0, 0, 0, "Ans", 4, 4, 4, 4) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 0, 0, 1, "Ans", 5, 5, 5, 23) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 0, 1, 0, "Ans", 5, 5, 23, 23) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 0, 1, 1, "Ans", 6, 6, 24, 139) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 1, 0, 0, "Ans", 5, 25, 25, 25) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 1, 0, 1, "Ans", 6, 30, 30, 162) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 1, 1, 0, "Ans", 6, 27, 151, 151) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 1, 1, 1, "Ans", 6, 30, 145, 897) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(BasicTest(4, 1, 1, 1, 1, "Zn" , 7, 46, 287, 2543) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch tests with a square parameter space XLAL_CHECK_MAIN(MismatchSquareTest("Zn", 0.03, 0, 0, 21460, Z1_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchSquareTest("Zn", 2e-4, -2e-9, 0, 23763, Z2_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchSquareTest("Zn", 1e-4, -1e-9, 1e-17, 19550, Z3_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchSquareTest("Ans", 0.03, 0, 0, 21460, A1s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchSquareTest("Ans", 2e-4, -2e-9, 0, 18283, A2s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchSquareTest("Ans", 1e-4, -2e-9, 2e-17, 20268, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch tests with an age--braking index parameter space XLAL_CHECK_MAIN(MismatchAgeBrakeTest("Ans", 100, 4.0e-5, 37872, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchAgeBrakeTest("Ans", 200, 1.5e-5, 37232, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(MismatchAgeBrakeTest("Ans", 300, 1.0e-5, 37022, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); // Perform mismatch tests with the reduced supersky parameter space and metric XLAL_CHECK_MAIN(SuperskyTest(1.1, 0.8, "Ans", 1, 50, 2.0e-5, 20548, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(SuperskyTest(1.5, 0.8, "Ans", 3, 50, 2.0e-5, 20202, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN(SuperskyTest(2.5, 0.8, "Ans", 17, 50, 2.0e-5, 29147, A3s_mism_hist) == XLAL_SUCCESS, XLAL_EFUNC); // Perform tests with the reduced supersky parameter space metric and multiple segments XLAL_CHECK_MAIN(MultiSegSuperskyTest() == XLAL_SUCCESS, XLAL_EFUNC); return EXIT_SUCCESS; }
/** * Set up a full multi-dimensional grid-scan. * Currently this only emulates a 'factored' grid-scan with 'sky x Freq x f1dot ...' , but * keeps all details within the DopplerScan module for future extension to real multidimensional * grids. * * NOTE: Use 'XLALNextDopplerPos()' to step through this template grid. * */ DopplerFullScanState * XLALInitDopplerFullScan ( const DopplerFullScanInit *init /**< [in] initialization parameters */ ) { XLAL_CHECK_NULL ( init != NULL, XLAL_EINVAL ); DopplerFullScanState *thisScan; XLAL_CHECK_NULL ( (thisScan = LALCalloc (1, sizeof(*thisScan) )) != NULL, XLAL_ENOMEM ); thisScan->gridType = init->gridType; /* store the user-input spinRange (includes refTime) in DopplerFullScanState */ thisScan->spinRange.refTime = init->searchRegion.refTime; memcpy ( thisScan->spinRange.fkdot, init->searchRegion.fkdot, sizeof(PulsarSpins) ); memcpy ( thisScan->spinRange.fkdotBand, init->searchRegion.fkdotBand, sizeof(PulsarSpins) ); // check that some old metric-codes aren't used with refTime!=startTime, which they don't handle correctly switch ( thisScan->gridType ) { case GRID_METRIC: case GRID_METRIC_SKYFILE: case GRID_SPINDOWN_SQUARE: /* square parameter space */ case GRID_SPINDOWN_AGEBRK: /* age-braking index parameter space */ XLAL_CHECK_NULL ( XLALGPSDiff ( &init->startTime, &init->searchRegion.refTime ) == 0, XLAL_EINVAL, "NOTE: gridType={metric,4,spin-square,spin-age-brk} only work for refTime (%f) == startTime (%f)!\n", XLALGPSGetREAL8(&(init->searchRegion.refTime)), XLALGPSGetREAL8(&(init->startTime)) );; break; default: break; } /* which "class" of template grid to generate?: factored, or full-multidim ? */ switch ( thisScan->gridType ) { /* emulate old 'factored' grids 'sky x f0dot x f1dot x f2dot x f3dot': */ case GRID_FLAT: case GRID_ISOTROPIC: case GRID_METRIC: case GRID_FILE_SKYGRID: case GRID_METRIC_SKYFILE: /* backwards-compatibility mode */ XLAL_CHECK_NULL ( XLALInitFactoredGrid ( thisScan, init ) == XLAL_SUCCESS, XLAL_EFUNC ); break; /* ----- multi-dimensional covering of full parameter space ----- */ case GRID_FILE_FULLGRID: XLAL_CHECK_NULL ( XLALLoadFullGridFile ( thisScan, init ) == XLAL_SUCCESS, XLAL_EFUNC ); break; case GRID_SPINDOWN_SQUARE: /* square parameter space */ case GRID_SPINDOWN_AGEBRK: /* age-braking index parameter space */ { const size_t n = 2 + PULSAR_MAX_SPINS; /* Check that the reference time is the same as the start time */ XLAL_CHECK_NULL ( XLALGPSCmp ( &thisScan->spinRange.refTime, &init->startTime) == 0, XLAL_EINVAL, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option currently restricts the reference time to be the same as the start time.\n"); /* Create a vector to hold lattice tiling parameter-space points */ XLAL_CHECK_NULL ( (thisScan->spindownTilingPoint = gsl_vector_alloc(n)) != NULL, XLAL_ENOMEM, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: gsl_vector_alloc failed\n"); /* Create a lattice tiling */ XLAL_CHECK_NULL ( (thisScan->spindownTiling = XLALCreateLatticeTiling(n)) != NULL, XLAL_EFUNC ); /* Parse the sky region string and check that it consists of only one point, and set bounds on it */ SkyRegion XLAL_INIT_DECL(sky); XLAL_CHECK_NULL ( XLALParseSkyRegionString ( &sky, init->searchRegion.skyRegionString ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK_NULL ( sky.numVertices == 1, XLAL_EINVAL, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option can only handle a single sky position.\n"); XLAL_CHECK_NULL ( sky.vertices[0].system == COORDINATESYSTEM_EQUATORIAL, XLAL_EINVAL, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option only understands COORDINATESYSTEM_EQUATORIAL\n"); XLAL_CHECK_NULL ( XLALSetLatticeTilingConstantBound(thisScan->spindownTiling, 0, sky.vertices[0].longitude, sky.vertices[0].longitude) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK_NULL ( XLALSetLatticeTilingConstantBound(thisScan->spindownTiling, 1, sky.vertices[0].latitude, sky.vertices[0].latitude) == XLAL_SUCCESS, XLAL_EFUNC ); if (sky.vertices) { XLALFree (sky.vertices); } /* Set up parameter space */ if (thisScan->gridType == GRID_SPINDOWN_SQUARE) { /* square parameter space */ /* Set square bounds on the frequency and spindowns */ for (size_t i = 0; i < PULSAR_MAX_SPINS; ++i) { XLAL_CHECK_NULL ( XLALSetLatticeTilingConstantBound(thisScan->spindownTiling, 2 + i, init->searchRegion.fkdot[i], init->searchRegion.fkdot[i] + init->searchRegion.fkdotBand[i]) == XLAL_SUCCESS, XLAL_EFUNC ); } } else if (thisScan->gridType == GRID_SPINDOWN_AGEBRK) { /* age-braking index parameter space */ /* Get age and braking index from extra arguments */ const REAL8 spindownAge = init->extraArgs[0]; const REAL8 minBraking = init->extraArgs[1]; const REAL8 maxBraking = init->extraArgs[2]; /* Set age-braking index parameter space */ XLAL_CHECK_NULL ( XLAL_SUCCESS == XLALSetLatticeTilingConstantBound(thisScan->spindownTiling, 2, init->searchRegion.fkdot[0], init->searchRegion.fkdot[0] + init->searchRegion.fkdotBand[0]), XLAL_EFUNC ); XLAL_CHECK_NULL ( XLAL_SUCCESS == XLALSetLatticeTilingF1DotAgeBrakingBound(thisScan->spindownTiling, 2, 3, spindownAge, minBraking, maxBraking), XLAL_EFUNC ); XLAL_CHECK_NULL ( XLAL_SUCCESS == XLALSetLatticeTilingF2DotBrakingBound(thisScan->spindownTiling, 2, 3, 4, minBraking, maxBraking), XLAL_EFUNC ); /* This current only goes up to second spindown, so bound higher dimensions */ for (size_t i = 3; i < PULSAR_MAX_SPINS; ++i) { XLAL_CHECK_NULL ( XLAL_SUCCESS == XLALSetLatticeTilingConstantBound(thisScan->spindownTiling, 2 + i, init->searchRegion.fkdot[i], init->searchRegion.fkdot[i] + init->searchRegion.fkdotBand[i]), XLAL_EFUNC ); } } /* Create a lattice tiling with Anstar lattice and spindown metric */ gsl_matrix* metric; XLAL_CHECK_NULL ( (metric = gsl_matrix_alloc(n, n)) != NULL, XLAL_ENOMEM ); gsl_matrix_set_identity(metric); gsl_matrix_view spin_metric = gsl_matrix_submatrix(metric, 2, 2, PULSAR_MAX_SPINS, PULSAR_MAX_SPINS); XLAL_CHECK_NULL ( XLALSpindownMetric(&spin_metric.matrix, init->Tspan) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK_NULL ( XLALSetTilingLatticeAndMetric(thisScan->spindownTiling, "Ans", metric, init->metricMismatch) == XLAL_SUCCESS, XLAL_EFUNC ); /* Create iterator over flat lattice tiling */ XLAL_CHECK_NULL ( (thisScan->spindownTilingItr = XLALCreateLatticeTilingIterator(thisScan->spindownTiling, n)) != NULL, XLAL_EFUNC ); /* Cleanup */ gsl_matrix_free(metric); } break; default: XLAL_ERROR_NULL ( XLAL_EINVAL, "\nInvalid grid type '%d'\n\n", init->gridType ); break; } /* switch gridType */ /* we're ready */ thisScan->state = STATE_READY; /* return result */ return thisScan; } // XLALInitDopplerFullScan()