LALPulsar 7.1.1.1-eeff03c
SuperskyMetrics.h
Go to the documentation of this file.
1//
2// Copyright (C) 2014--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#ifndef _SUPERSKYMETRICS_H
21#define _SUPERSKYMETRICS_H
22
23#include <stdbool.h>
24#include <gsl/gsl_matrix.h>
25#include <lal/LALStdlib.h>
26#include <lal/UniversalDopplerMetric.h>
27#include <lal/LatticeTiling.h>
28#include <lal/FITSFileIO.h>
29
30#ifdef __cplusplus
31extern "C" {
32#endif
33
34///
35/// \defgroup SuperskyMetrics_h Header SuperskyMetrics.h
36/// \ingroup lalpulsar_metric
37/// \author Karl Wette
38/// \brief Compute the supersky metrics and coordinate transforms of \cite WettePrix2013a and \cite Wette2015a .
39///
40/// @{
41///
42
43///
44/// Supersky metric coordinate transform data
45///
46typedef struct tagSuperskyTransformData SuperskyTransformData;
47
48///
49/// Type of supersky metric to compute.
50///
51typedef enum tagSuperskyMetricType {
52 SUPERSKY_METRIC_TYPE, ///< Metric for all-sky searches
53 SUPERSKY_DIRECTED_METRIC_TYPE, ///< Metric for directed searches
56
57///
58/// Computed supersky metrics, returned by XLALComputeSuperskyMetrics().
59///
60#ifdef SWIG // SWIG interface directives
61SWIGLAL( ARRAY_MULTIPLE_LENGTHS( tagSuperskyMetrics, num_segments ) );
62#endif // SWIG
63typedef struct tagSuperskyMetrics {
64 size_t num_segments; ///< Number of segments
65
66#ifdef SWIG // SWIG interface directives
67 SWIGLAL( ARRAY_1D( SuperskyMetrics, gsl_matrix *, coh_rssky_metric, size_t, num_segments ) );
68#endif // SWIG
69 gsl_matrix **coh_rssky_metric; ///< Coherent reduced supersky metric (2-dimensional sky) for each segment
70#ifdef SWIG // SWIG interface directives
71 SWIGLAL( ARRAY_1D( SuperskyMetrics, SuperskyTransformData *, coh_rssky_transf, size_t, num_segments ) );
72#endif // SWIG
73 SuperskyTransformData **coh_rssky_transf; ///< Coherent reduced supersky metric coordinate transform data for each segment
74
75 gsl_matrix *semi_rssky_metric; ///< Semicoherent reduced supersky metric (2-dimensional sky)
76 SuperskyTransformData *semi_rssky_transf; ///< Semicoherent reduced supersky metric coordinate transform data
77
79
80///
81/// Compute the supersky metrics, which are returned in a \c SuperskyMetrics struct.
82///
84 const SuperskyMetricType metric_type, ///< [in] Type of supersky metric to compute
85 const size_t spindowns, ///< [in] Number of frequency+spindown coordinates
86 const LIGOTimeGPS *ref_time, ///< [in] Reference time for the metrics
87 const LALSegList *segments, ///< [in] List of segments to compute metrics over
88 const double fiducial_freq, ///< [in] Fiducial frequency for sky-position coordinates
89 const MultiLALDetector *detectors, ///< [in] List of detectors to average metrics over
90 const MultiNoiseFloor *detector_weights, ///< [in] Weights used to combine single-detector metrics (default: unit weights)
91 const DetectorMotionType detector_motion, ///< [in] Which detector motion to use
92 const EphemerisData *ephemerides ///< [in] Earth/Sun ephemerides
93);
94
95///
96/// Copy a \c SuperskyMetrics struct.
97///
99 const SuperskyMetrics *metrics ///< [in] Supersky metrics struct
100);
101
102///
103/// Destroy a \c SuperskyMetrics struct.
104///
106 SuperskyMetrics *metrics ///< [in] Supersky metrics struct
107);
108
109///
110/// Write a \c SuperskyMetrics struct to a FITS file.
111///
113 FITSFile *file, ///< [in] FITS file pointer
114 const SuperskyMetrics *metrics ///< [in] Supersky metrics struct
115);
116
117///
118/// Read a \c SuperskyMetrics struct from a FITS file.
119///
121 FITSFile *file, ///< [in] FITS file pointer
122 SuperskyMetrics **metrics ///< [out] Supersky metrics struct
123);
124
125///
126/// Return dimensions of the supersky metrics.
127///
129 const SuperskyMetrics *metrics, ///< [in] Supersky metrics struct
130 size_t *spindowns ///< [out] Number of spindown dimensions
131);
132
133///
134/// Scale all supersky metrics and their coordinate transform data to a new fiducial frequency.
135///
137 SuperskyMetrics *metrics, ///< [in] Supersky metrics struct
138 const double new_fiducial_freq ///< [in] New fiducial frequency
139);
140
141///
142/// Project and rescale the reduced supersky metrics in the frequency dimension, such that all
143/// reduced supersky metrics have the same frequency spacing for the given maximum mismatches.
144///
146 SuperskyMetrics *metrics, ///< [in] Supersky metrics struct
147 const double coh_max_mismatch, ///< [in] Maximum coherent mismatch
148 const double semi_max_mismatch ///< [in] Maximum semicoherent mismatch
149);
150
151///
152/// Set the reference time of a physical point to that of the reduced supersky coordinates.
153///
155 PulsarDopplerParams *out_phys, ///< [out] Output point in physical coordinates
156 const SuperskyTransformData *rssky_transf ///< [in] Reduced supersky coordinate transform data
157);
158
159///
160/// Convert a point from physical to supersky coordinates.
161///
163 gsl_vector *out_rssky, ///< [out] Output point in supersky coordinates
164 const PulsarDopplerParams *in_phys, ///< [in] Input point in physical coordinates
165 const SuperskyTransformData *rssky_transf ///< [in] Reduced supersky coordinate transform data
166);
167
168///
169/// Convert a point from supersky to physical coordinates.
170///
172 PulsarDopplerParams *out_phys, ///< [out] Output point in physical coordinates
173 const gsl_vector *in_rssky, ///< [in] Input point in supersky coordinates
174 const gsl_vector *ref_rssky, ///< [in,optional] Reference point in supersky coordinates
175 const SuperskyTransformData *rssky_transf ///< [in] Reduced supersky coordinate transform data
176);
177
178///
179/// Convert a point between supersky coordinates. The vectors \c out_rssky and \c in_rssky may be the same.
180///
182 gsl_vector *out_rssky, ///< [out] Output point in supersky coordinates
183 const SuperskyTransformData *out_rssky_transf,///< [in] Output reduced supersky coordinate transform data
184 const gsl_vector *in_rssky, ///< [in] Input point in supersky coordinates
185 const gsl_vector *ref_rssky, ///< [in,optional] Reference point in supersky coordinates
186 const SuperskyTransformData *in_rssky_transf ///< [in] Input reduced supersky coordinate transform data
187);
188
189///
190/// Convert a set of points from physical to supersky coordinates.
191///
192#ifdef SWIG // SWIG interface directives
193SWIGLAL( INOUT_STRUCTS( gsl_matrix **, out_rssky ) );
194#endif
196 gsl_matrix **out_rssky, ///< [out] Columns are output point in supersky coordinates
197 const gsl_matrix *in_phys, ///< [in] Columns are input point in physical coordinates
198 const SuperskyTransformData *rssky_transf ///< [in] Reduced supersky coordinate transform data
199);
200
201///
202/// Convert a set of points from supersky to physical coordinates.
203///
204#ifdef SWIG // SWIG interface directives
205SWIGLAL( INOUT_STRUCTS( gsl_matrix **, out_phys ) );
206#endif
208 gsl_matrix **out_phys, ///< [out] Columns are output point in physical coordinates
209 const gsl_matrix *in_rssky, ///< [in] Columns are input point in supersky coordinates
210 const SuperskyTransformData *rssky_transf ///< [in] Reduced supersky coordinate transform data
211);
212
213#ifdef SWIG // SWIG interface directives
214SWIGLAL( COPYINOUT_ARRAYS( gsl_matrix, rssky_metric, rssky_transf ) );
215#endif // SWIG
216
217///
218/// Set parameter-space bounds on the physical sky position \f$ (\alpha, \delta) \f$ for a lattice
219/// tiling using the reduced supersky metric. The metric and coordinate transform data must be supplied,
220/// since they will be transformed such that the physical sky region maps to a region in the reduced
221/// supersky coordinates \f$ (n_a,n_b) \f$ which may be covered by the lattice tiling.
222///
224 LatticeTiling *tiling, ///< [in] Lattice tiling
225 gsl_matrix *rssky_metric, ///< [in,out] Reduced supersky metric
226 SuperskyTransformData *rssky_transf, ///< [in,out] Reduced supersky coordinate transform data
227 const double alpha1, ///< [in] First bound on sky position right ascension
228 const double alpha2, ///< [in] Second bound on sky position right ascension
229 const double delta1, ///< [in] First bound on sky position declination
230 const double delta2 ///< [in] Second bound on sky position declination
231);
232
233///
234/// Divide the reduced supersky parameter space into \p patch_count equal-area patches, and
235/// set parameter-space bounds on the reduced supersky coordinates \f$ (n_a,n_b) \f$ for the
236/// patch indexed by \p patch_index.
237///
239 LatticeTiling *tiling, ///< [in] Lattice tiling
240 const gsl_matrix *rssky_metric, ///< [in] Reduced supersky metric
241 const double max_mismatch, ///< [in] Maximum prescribed mismatch
242 const UINT4 patch_count, ///< [in] Number of equal-area patches to divide sky into
243 const UINT4 patch_index ///< [in] Index of the patch for which to compute bounds
244);
245
246#ifdef SWIG // SWIG interface directives
247SWIGLAL_CLEAR( COPYINOUT_ARRAYS( gsl_matrix, rssky_metric, rssky_transf ) );
248#endif // SWIG
249
250///
251/// Set parameter-space bounds on the physical frequency/spindowns \f$ f^{(s)} \f$ for a lattice
252/// tiling using the reduced supersky metric.
253///
255 LatticeTiling *tiling, ///< [in] Lattice tiling
256 const SuperskyTransformData *rssky_transf, ///< [in] Reduced supersky coordinate transform data
257 const size_t s, ///< [in] Spindown order; 0=frequency, 1=first spindown, etc.
258 const double bound1, ///< [in] First bound on frequency/spindown
259 const double bound2 ///< [in] Second bound on frequency/spindown
260);
261
262///
263/// Set parameter-space bound padding on the physical frequency/spindowns \f$ f^{(s)} \f$ for a lattice
264/// tiling using the reduced supersky metric.
265///
267 LatticeTiling *tiling, ///< [in] Lattice tiling
268 const SuperskyTransformData *rssky_transf, ///< [in] Reduced supersky coordinate transform data
269 const size_t s, ///< [in] Spindown order; 0=frequency, 1=first spindown, etc.
270 const bool padding ///< [in] Whether bounds are padded
271);
272
273///
274/// Register a lattice tiling callback function which computes the physical range covered
275/// by a reduced supersky lattice tiling.
276///
277#ifdef SWIG // SWIG interface directives
278SWIGLAL( OUTPUT_OWNED_BY_1ST_ARG( PulsarDopplerParams **, min_phys ) );
279SWIGLAL( OUTPUT_OWNED_BY_1ST_ARG( PulsarDopplerParams **, max_phys ) );
280#endif
282 LatticeTiling *tiling, ///< [in] Lattice tiling
283 const SuperskyTransformData *rssky_transf, ///< [in] Reduced supersky coordinate transform data
284 const PulsarDopplerParams **min_phys, ///< [out] Minimum physical range
285 const PulsarDopplerParams **max_phys ///< [out] Maximum physical range
286);
287
288///
289/// Register a lattice tiling callback function which computes the range covered by a
290/// reduced supersky lattice tiling in another set of reduced supersky coordinates.
291///
292#ifdef SWIG // SWIG interface directives
293SWIGLAL( OUTPUT_OWNED_BY_1ST_ARG( gsl_vector **, min_rssky2 ) );
294SWIGLAL( OUTPUT_OWNED_BY_1ST_ARG( gsl_vector **, max_rssky2 ) );
295#endif
297 LatticeTiling *tiling, ///< [in] Lattice tiling
298 const SuperskyTransformData *rssky_transf, ///< [in] Reduced supersky coordinate transform data
299 const SuperskyTransformData *rssky2_transf, ///< [in] Other reduced supersky coordinate transform data
300 const gsl_vector **min_rssky2, ///< [out] Minimum range of other reduced supersky coordinates
301 const gsl_vector **max_rssky2 ///< [out] Maximum range of other reduced supersky coordinates
302);
303
304///
305/// Set parameter-space bounds on an entire lattice tiling given minimum and maximum
306/// ranges in reduced supersky coordinates.
307///
309 LatticeTiling *tiling, ///< [in] Lattice tiling
310 const gsl_vector *min_rssky, ///< [in] Minimum range of reduced supersky coordinates
311 const gsl_vector *max_rssky ///< [in] Maximum range of reduced supersky coordinates
312);
313
314/// @}
315
316#ifdef __cplusplus
317}
318#endif
319
320#endif // _SUPERSKYMETRICS_H
321
322// Local Variables:
323// c-file-style: "linux"
324// c-basic-offset: 2
325// End:
struct tagFITSFile FITSFile
Representation of a FITS file.
Definition: FITSFileIO.h:54
uint32_t UINT4
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.
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...
void XLALDestroySuperskyMetrics(SuperskyMetrics *metrics)
Destroy a SuperskyMetrics struct.
SuperskyMetrics * XLALCopySuperskyMetrics(const SuperskyMetrics *metrics)
Copy a SuperskyMetrics struct.
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 XLALConvertPhysicalToSuperskyPoints(gsl_matrix **out_rssky, const gsl_matrix *in_phys, const SuperskyTransformData *rssky_transf)
Convert a set of points from physical to supersky coordinates.
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...
SuperskyMetricType
Type of supersky metric to compute.
int XLALScaleSuperskyMetricsFiducialFreq(SuperskyMetrics *metrics, const double new_fiducial_freq)
Scale all supersky metrics and their coordinate transform data to a new fiducial frequency.
int XLALSetPhysicalPointSuperskyRefTime(PulsarDopplerParams *out_phys, const SuperskyTransformData *rssky_transf)
Set the reference time of a physical point to that of the reduced supersky coordinates.
int XLALConvertPhysicalToSuperskyPoint(gsl_vector *out_rssky, const PulsarDopplerParams *in_phys, const SuperskyTransformData *rssky_transf)
Convert a point from physical to supersky coordinates.
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 XLALFITSWriteSuperskyMetrics(FITSFile *file, const SuperskyMetrics *metrics)
Write a SuperskyMetrics struct to a FITS file.
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 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.
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...
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.
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...
int XLALFITSReadSuperskyMetrics(FITSFile *file, SuperskyMetrics **metrics)
Read a SuperskyMetrics struct from a FITS file.
@ SUPERSKY_METRIC_TYPE
Metric for all-sky searches.
@ SUPERSKY_DIRECTED_METRIC_TYPE
Metric for directed searches.
@ MAX_METRIC_TYPE
DetectorMotionType
Bitfield of different types of detector-motion to use in order to compute the Doppler-metric.
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
array of detectors definitions 'LALDetector'
array of detector-specific 'noise floors' (ie PSD values), assumed constant over the frequency-band o...
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
Computed supersky metrics, returned by XLALComputeSuperskyMetrics().
gsl_matrix * semi_rssky_metric
Semicoherent reduced supersky metric (2-dimensional sky)
gsl_matrix ** coh_rssky_metric
Coherent reduced supersky metric (2-dimensional sky) for each segment.
size_t num_segments
Number of segments.
SuperskyTransformData * semi_rssky_transf
Semicoherent reduced supersky metric coordinate transform data.
SuperskyTransformData ** coh_rssky_transf
Coherent reduced supersky metric coordinate transform data for each segment.
Supersky metric coordinate transform data.