Compute the supersky metrics and coordinate transforms of [29] and [32] .
Prototypes | |
| SuperskyMetrics * | XLALComputeSuperskyMetrics (const SuperskyMetricType metric_type, const size_t spindowns, const LIGOTimeGPS *ref_time, const LALSegList *segments, const double fiducial_freq, const MultiLALDetector *detectors, const MultiNoiseFloor *detector_weights, const DetectorMotionType detector_motion, const EphemerisData *ephemerides) |
Compute the supersky metrics, which are returned in a SuperskyMetrics struct. More... | |
| SuperskyMetrics * | XLALCopySuperskyMetrics (const SuperskyMetrics *metrics) |
Copy a SuperskyMetrics struct. More... | |
| void | XLALDestroySuperskyMetrics (SuperskyMetrics *metrics) |
Destroy a SuperskyMetrics struct. More... | |
| int | XLALFITSWriteSuperskyMetrics (FITSFile *file, const SuperskyMetrics *metrics) |
Write a SuperskyMetrics struct to a FITS file. More... | |
| int | XLALFITSReadSuperskyMetrics (FITSFile *file, SuperskyMetrics **metrics) |
Read a SuperskyMetrics struct from a FITS file. More... | |
| int | XLALSuperskyMetricsDimensions (const SuperskyMetrics *metrics, size_t *spindowns) |
| Return dimensions of the supersky metrics. More... | |
| int | XLALScaleSuperskyMetricsFiducialFreq (SuperskyMetrics *metrics, const double new_fiducial_freq) |
| Scale all supersky metrics and their coordinate transform data to a new fiducial frequency. More... | |
| 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 supersky metrics have the same frequency spacing for the given maximum mismatches. More... | |
| int | XLALSetPhysicalPointSuperskyRefTime (PulsarDopplerParams *out_phys, const SuperskyTransformData *rssky_transf) |
| Set the reference time of a physical point to that of the reduced supersky coordinates. More... | |
| int | XLALConvertPhysicalToSuperskyPoint (gsl_vector *out_rssky, const PulsarDopplerParams *in_phys, const SuperskyTransformData *rssky_transf) |
| Convert a point from physical to supersky coordinates. More... | |
| 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. More... | |
| int | XLALConvertSuperskyToSuperskyPoint (gsl_vector *out_rssky, const SuperskyTransformData *out_rssky_transf, const gsl_vector *in_rssky, const gsl_vector *ref_rssky, const SuperskyTransformData *in_rssky_transf) |
| Convert a point between supersky coordinates. More... | |
| int | XLALConvertPhysicalToSuperskyPoints (gsl_matrix **out_rssky, const gsl_matrix *in_phys, const SuperskyTransformData *rssky_transf) |
| Convert a set of points from physical to supersky coordinates. More... | |
| int | XLALConvertSuperskyToPhysicalPoints (gsl_matrix **out_phys, const gsl_matrix *in_rssky, const SuperskyTransformData *rssky_transf) |
| Convert a set of points from supersky to physical coordinates. More... | |
| 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 \( (\alpha, \delta) \) for a lattice tiling using the reduced supersky metric. More... | |
| 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-space bounds on the reduced supersky coordinates \( (n_a,n_b) \) for the patch indexed by patch_index. More... | |
| 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 \( f^{(s)} \) for a lattice tiling using the reduced supersky metric. More... | |
| 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 \( f^{(s)} \) for a lattice tiling using the reduced supersky metric. More... | |
| 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 supersky lattice tiling. More... | |
| 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 lattice tiling in another set of reduced supersky coordinates. More... | |
| 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 supersky coordinates. More... | |
Data Structures | |
| struct | SuperskyMetrics |
| Computed supersky metrics, returned by XLALComputeSuperskyMetrics(). More... | |
Enumerations | |
| enum | SuperskyMetricType { SUPERSKY_METRIC_TYPE , SUPERSKY_DIRECTED_METRIC_TYPE , MAX_METRIC_TYPE } |
| Type of supersky metric to compute. More... | |
| SuperskyMetrics * XLALComputeSuperskyMetrics | ( | const SuperskyMetricType | metric_type, |
| const size_t | spindowns, | ||
| const LIGOTimeGPS * | ref_time, | ||
| const LALSegList * | segments, | ||
| const double | fiducial_freq, | ||
| const MultiLALDetector * | detectors, | ||
| const MultiNoiseFloor * | detector_weights, | ||
| const DetectorMotionType | detector_motion, | ||
| const EphemerisData * | ephemerides | ||
| ) |
Compute the supersky metrics, which are returned in a SuperskyMetrics struct.
| [in] | metric_type | Type of supersky metric to compute |
| [in] | spindowns | Number of frequency+spindown coordinates |
| [in] | ref_time | Reference time for the metrics |
| [in] | segments | List of segments to compute metrics over |
| [in] | fiducial_freq | Fiducial frequency for sky-position coordinates |
| [in] | detectors | List of detectors to average metrics over |
| [in] | detector_weights | Weights used to combine single-detector metrics (default: unit weights) |
| [in] | detector_motion | Which detector motion to use |
| [in] | ephemerides | Earth/Sun ephemerides |
Definition at line 603 of file SuperskyMetrics.c.
| SuperskyMetrics * XLALCopySuperskyMetrics | ( | const SuperskyMetrics * | metrics | ) |
Copy a SuperskyMetrics struct.
| [in] | metrics | Supersky metrics struct |
Definition at line 719 of file SuperskyMetrics.c.
| void XLALDestroySuperskyMetrics | ( | SuperskyMetrics * | metrics | ) |
Destroy a SuperskyMetrics struct.
| [in] | metrics | Supersky metrics struct |
Definition at line 758 of file SuperskyMetrics.c.
| int XLALFITSWriteSuperskyMetrics | ( | FITSFile * | file, |
| const SuperskyMetrics * | metrics | ||
| ) |
Write a SuperskyMetrics struct to a FITS file.
| [in] | file | FITS file pointer |
| [in] | metrics | Supersky metrics struct |
Definition at line 794 of file SuperskyMetrics.c.
| int XLALFITSReadSuperskyMetrics | ( | FITSFile * | file, |
| SuperskyMetrics ** | metrics | ||
| ) |
Read a SuperskyMetrics struct from a FITS file.
| [in] | file | FITS file pointer |
| [out] | metrics | Supersky metrics struct |
Definition at line 844 of file SuperskyMetrics.c.
| int XLALSuperskyMetricsDimensions | ( | const SuperskyMetrics * | metrics, |
| size_t * | spindowns | ||
| ) |
Return dimensions of the supersky metrics.
| [in] | metrics | Supersky metrics struct |
| [out] | spindowns | Number of spindown dimensions |
Definition at line 908 of file SuperskyMetrics.c.
| int XLALScaleSuperskyMetricsFiducialFreq | ( | SuperskyMetrics * | metrics, |
| const double | new_fiducial_freq | ||
| ) |
Scale all supersky metrics and their coordinate transform data to a new fiducial frequency.
| [in] | metrics | Supersky metrics struct |
| [in] | new_fiducial_freq | New fiducial frequency |
Definition at line 959 of file SuperskyMetrics.c.
| 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 supersky metrics have the same frequency spacing for the given maximum mismatches.
| [in] | metrics | Supersky metrics struct |
| [in] | coh_max_mismatch | Maximum coherent mismatch |
| [in] | semi_max_mismatch | Maximum semicoherent mismatch |
Definition at line 984 of file SuperskyMetrics.c.
| int XLALSetPhysicalPointSuperskyRefTime | ( | PulsarDopplerParams * | out_phys, |
| const SuperskyTransformData * | rssky_transf | ||
| ) |
Set the reference time of a physical point to that of the reduced supersky coordinates.
| [out] | out_phys | Output point in physical coordinates |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
Definition at line 1045 of file SuperskyMetrics.c.
| int XLALConvertPhysicalToSuperskyPoint | ( | gsl_vector * | out_rssky, |
| const PulsarDopplerParams * | in_phys, | ||
| const SuperskyTransformData * | rssky_transf | ||
| ) |
Convert a point from physical to supersky coordinates.
| [out] | out_rssky | Output point in supersky coordinates |
| [in] | in_phys | Input point in physical coordinates |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
Definition at line 1122 of file SuperskyMetrics.c.
| 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.
| [out] | out_phys | Output point in physical coordinates |
| [in] | in_rssky | Input point in supersky coordinates |
| ref_rssky | [in,optional] Reference point in supersky coordinates | |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
Definition at line 1190 of file SuperskyMetrics.c.
| int XLALConvertSuperskyToSuperskyPoint | ( | gsl_vector * | out_rssky, |
| const SuperskyTransformData * | out_rssky_transf, | ||
| const gsl_vector * | in_rssky, | ||
| const gsl_vector * | ref_rssky, | ||
| const SuperskyTransformData * | in_rssky_transf | ||
| ) |
Convert a point between supersky coordinates.
The vectors out_rssky and in_rssky may be the same.
| [out] | out_rssky | Output point in supersky coordinates |
| [in] | out_rssky_transf | Output reduced supersky coordinate transform data |
| [in] | in_rssky | Input point in supersky coordinates |
| ref_rssky | [in,optional] Reference point in supersky coordinates | |
| [in] | in_rssky_transf | Input reduced supersky coordinate transform data |
Definition at line 1254 of file SuperskyMetrics.c.
| int XLALConvertPhysicalToSuperskyPoints | ( | gsl_matrix ** | out_rssky, |
| const gsl_matrix * | in_phys, | ||
| const SuperskyTransformData * | rssky_transf | ||
| ) |
Convert a set of points from physical to supersky coordinates.
| [out] | out_rssky | Columns are output point in supersky coordinates |
| [in] | in_phys | Columns are input point in physical coordinates |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
Definition at line 1280 of file SuperskyMetrics.c.
| int XLALConvertSuperskyToPhysicalPoints | ( | gsl_matrix ** | out_phys, |
| const gsl_matrix * | in_rssky, | ||
| const SuperskyTransformData * | rssky_transf | ||
| ) |
Convert a set of points from supersky to physical coordinates.
| [out] | out_phys | Columns are output point in physical coordinates |
| [in] | in_rssky | Columns are input point in supersky coordinates |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
Definition at line 1326 of file SuperskyMetrics.c.
| 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 \( (\alpha, \delta) \) for a lattice tiling using the reduced supersky metric.
The metric and coordinate transform data must be supplied, since they will be transformed such that the physical sky region maps to a region in the reduced supersky coordinates \( (n_a,n_b) \) which may be covered by the lattice tiling.
| [in] | tiling | Lattice tiling |
| [in,out] | rssky_metric | Reduced supersky metric |
| [in,out] | rssky_transf | Reduced supersky coordinate transform data |
| [in] | alpha1 | First bound on sky position right ascension |
| [in] | alpha2 | Second bound on sky position right ascension |
| [in] | delta1 | First bound on sky position declination |
| [in] | delta2 | Second bound on sky position declination |
Definition at line 1475 of file SuperskyMetrics.c.
| 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-space bounds on the reduced supersky coordinates \( (n_a,n_b) \) for the patch indexed by patch_index.
| [in] | tiling | Lattice tiling |
| [in] | rssky_metric | Reduced supersky metric |
| [in] | max_mismatch | Maximum prescribed mismatch |
| [in] | patch_count | Number of equal-area patches to divide sky into |
| [in] | patch_index | Index of the patch for which to compute bounds |
Definition at line 2125 of file SuperskyMetrics.c.
| 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 \( f^{(s)} \) for a lattice tiling using the reduced supersky metric.
| [in] | tiling | Lattice tiling |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
| [in] | s | Spindown order; 0=frequency, 1=first spindown, etc. |
| [in] | bound1 | First bound on frequency/spindown |
| [in] | bound2 | Second bound on frequency/spindown |
Definition at line 2388 of file SuperskyMetrics.c.
| 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 \( f^{(s)} \) for a lattice tiling using the reduced supersky metric.
| [in] | tiling | Lattice tiling |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
| [in] | s | Spindown order; 0=frequency, 1=first spindown, etc. |
| [in] | padding | Whether bounds are padded |
Definition at line 2434 of file SuperskyMetrics.c.
| 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 supersky lattice tiling.
| [in] | tiling | Lattice tiling |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
| [out] | min_phys | Minimum physical range |
| [out] | max_phys | Maximum physical range |
Definition at line 2524 of file SuperskyMetrics.c.
| 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 lattice tiling in another set of reduced supersky coordinates.
| [in] | tiling | Lattice tiling |
| [in] | rssky_transf | Reduced supersky coordinate transform data |
| [in] | rssky2_transf | Other reduced supersky coordinate transform data |
| [out] | min_rssky2 | Minimum range of other reduced supersky coordinates |
| [out] | max_rssky2 | Maximum range of other reduced supersky coordinates |
Definition at line 2615 of file SuperskyMetrics.c.
| 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 supersky coordinates.
| [in] | tiling | Lattice tiling |
| [in] | min_rssky | Minimum range of reduced supersky coordinates |
| [in] | max_rssky | Maximum range of reduced supersky coordinates |
Definition at line 2649 of file SuperskyMetrics.c.
| enum SuperskyMetricType |
Type of supersky metric to compute.
| Enumerator | |
|---|---|
| SUPERSKY_METRIC_TYPE | Metric for all-sky searches. |
| SUPERSKY_DIRECTED_METRIC_TYPE | Metric for directed searches. |
| MAX_METRIC_TYPE | |
Definition at line 51 of file SuperskyMetrics.h.