/** * Destroy the a full DopplerFullScanState structure */ void XLALDestroyDopplerFullScan ( DopplerFullScanState *scan ) { if ( scan == NULL ) { return; } if ( scan->factoredScan ) { XLALDestroyDopplerSkyScan ( &(scan->factoredScan->skyScan) ); XLALFree ( scan->factoredScan ); } if ( scan->covering ) { XLALREAL8VectorListDestroy ( scan->covering ); } if (scan->spindownTiling) { XLALDestroyLatticeTiling(scan->spindownTiling); } if (scan->spindownTilingItr) { XLALDestroyLatticeTilingIterator(scan->spindownTilingItr); } if (scan->spindownTilingPoint) { gsl_vector_free(scan->spindownTilingPoint); } if ( scan->skyRegion.vertices) { XLALFree ( scan->skyRegion.vertices); } XLALFree ( scan ); return; } // XLALDestroyDopplerFullScan()
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 MismatchTest( LatticeTiling *tiling, gsl_matrix *metric, const double max_mismatch, const UINT8 total_ref, const double mism_hist_ref[MISM_HIST_BINS] ) { const size_t n = XLALTotalLatticeTilingDimensions(tiling); // Create lattice tiling iterator and locator LatticeTilingIterator *itr = XLALCreateLatticeTilingIterator(tiling, n); XLAL_CHECK(itr != NULL, XLAL_EFUNC); LatticeTilingLocator *loc = XLALCreateLatticeTilingLocator(tiling); XLAL_CHECK(loc != NULL, XLAL_EFUNC); // Count number of points const UINT8 total = XLALTotalLatticeTilingPoints(itr); printf("Number of lattice points: %" LAL_UINT8_FORMAT "\n", total); XLAL_CHECK(imaxabs(total - total_ref) <= 1, XLAL_EFUNC, "ERROR: |total - total_ref| = |%" LAL_UINT8_FORMAT " - %" LAL_UINT8_FORMAT "| > 1", total, total_ref); // 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); // Initialise mismatch histogram counts double mism_hist[MISM_HIST_BINS] = {0}; double mism_hist_total = 0, mism_hist_out_of_range = 0; // Perform 10 injections for every template { gsl_matrix *GAMAT(injections, 3, total); gsl_matrix *GAMAT(nearest, 3, total); gsl_matrix *GAMAT(temp, 3, total); RandomParams *rng = XLALCreateRandomParams(total); XLAL_CHECK(rng != NULL, XLAL_EFUNC); for (size_t i = 0; i < 10; ++i) { // Generate random injection points XLAL_CHECK(XLALRandomLatticeTilingPoints(tiling, 0.0, rng, injections) == XLAL_SUCCESS, XLAL_EFUNC); // Find nearest lattice template points XLAL_CHECK(XLALNearestLatticeTilingPoints(loc, injections, &nearest, NULL) == XLAL_SUCCESS, XLAL_EFUNC); // Compute mismatch between injections gsl_matrix_sub(nearest, injections); gsl_blas_dsymm(CblasLeft, CblasUpper, 1.0, metric, nearest, 0.0, temp); for (size_t j = 0; j < temp->size2; ++j) { gsl_vector_view temp_j = gsl_matrix_column(temp, j); gsl_vector_view nearest_j = gsl_matrix_column(nearest, j); double mismatch = 0.0; gsl_blas_ddot(&nearest_j.vector, &temp_j.vector, &mismatch); mismatch /= max_mismatch; // Increment mismatch histogram counts ++mism_hist_total; if (mismatch < 0.0 || mismatch > 1.0) { ++mism_hist_out_of_range; } else { ++mism_hist[lround(floor(mismatch * MISM_HIST_BINS))]; } } } // Cleanup GFMAT(injections, nearest, temp); XLALDestroyRandomParams(rng); } // Normalise histogram for (size_t i = 0; i < MISM_HIST_BINS; ++i) { mism_hist[i] *= MISM_HIST_BINS / mism_hist_total; } // Print mismatch histogram and its reference printf("Mismatch histogram: "); for (size_t i = 0; i < MISM_HIST_BINS; ++i) { printf(" %0.3f", mism_hist[i]); } printf("\n"); printf("Reference histogram:"); for (size_t i = 0; i < MISM_HIST_BINS; ++i) { printf(" %0.3f", mism_hist_ref[i]); } printf("\n"); // Determine error between mismatch histogram and its reference double mism_hist_error = 0.0; for (size_t i = 0; i < MISM_HIST_BINS; ++i) { mism_hist_error += fabs(mism_hist[i] - mism_hist_ref[i]); } mism_hist_error /= MISM_HIST_BINS; printf("Mismatch histogram error: %0.3e\n", mism_hist_error); const double mism_hist_error_tol = 5e-2; if (mism_hist_error >= mism_hist_error_tol) { XLAL_ERROR(XLAL_EFAILED, "ERROR: mismatch histogram error exceeds %0.3e\n", mism_hist_error_tol); } // Check fraction of injections out of histogram range const double mism_out_of_range = mism_hist_out_of_range / mism_hist_total; printf("Fraction of points out of histogram range: %0.3e\n", mism_out_of_range); const double mism_out_of_range_tol = 2e-3; if (mism_out_of_range > mism_out_of_range_tol) { XLAL_ERROR(XLAL_EFAILED, "ERROR: fraction of points out of histogram range exceeds %0.3e\n", mism_out_of_range_tol); } // Perform 10 injections outside parameter space { gsl_matrix *GAMAT(injections, 3, 10); gsl_matrix *GAMAT(nearest, n, total); RandomParams *rng = XLALCreateRandomParams(total); XLAL_CHECK(rng != NULL, XLAL_EFUNC); // Generate random injection points outside parameter space XLAL_CHECK(XLALRandomLatticeTilingPoints(tiling, 5.0, rng, injections) == XLAL_SUCCESS, XLAL_EFUNC); // Find nearest lattice template points XLAL_CHECK(XLALNearestLatticeTilingPoints(loc, injections, &nearest, NULL) == XLAL_SUCCESS, XLAL_EFUNC); // Cleanup GFMAT(injections, nearest); XLALDestroyRandomParams(rng); } // Cleanup XLALDestroyLatticeTiling(tiling); XLALDestroyLatticeTilingIterator(itr); XLALDestroyLatticeTilingLocator(loc); GFMAT(metric, points); LALCheckMemoryLeaks(); printf("\n"); fflush(stdout); return XLAL_SUCCESS; }