Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.0.1-8b258c9
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Header MetricUtils.h

Detailed Description

Various useful utility functions for working with CW parameter-space metrics.

Author
Reinhard Prix, Karl Wette

Prototypes

REAL8 XLALCompareMetrics (const gsl_matrix *g1_ij, const gsl_matrix *g2_ij)
 Flexible comparison function between two metrics \( g^1_{ij} \) and \( g^2_{ij} \) . More...
 
double XLALMetricDeterminant (const gsl_matrix *g_ij)
 Compute the determinant of a metric \( g_{ij} \) . More...
 
gsl_vector * XLALMetricEllipseBoundingBox (const gsl_matrix *g_ij, const double max_mismatch)
 Compute the extent of the bounding box of the mismatch ellipse of a metric \( g_{ij} \) . More...
 
int XLALProjectMetric (gsl_matrix **p_gpr_ij, const gsl_matrix *g_ij, const UINT4 c)
 Calculate the projected metric onto the subspace orthogonal to coordinate-axis c, namely \( g^{\prime}_{ij} = g_{ij} - ( g_{ic} g_{jc} / g_{cc} ) \) , where \( c \) is the index of the projected coordinate. More...
 
int XLALCholeskyLDLTDecompMetric (gsl_matrix **p_cholesky, const gsl_matrix *g_ij)
 Decompose a metric \( \mathbf{G} \) as \( \mathbf{G} = \mathbf{L} \mathbf{D} \mathbf{L}^{\mathrm{T}} \) , where \( \mathbf{L} \) is a lower-triangular matrix with unit diagonal, and \( \mathbf{D} \) is a diagonal matrix. More...
 
int XLALTransformMetric (gsl_matrix **p_gpr_ij, const gsl_matrix *transform, const gsl_matrix *g_ij)
 Apply a transform \( \mathbf{A} = (a_{ij}) \) to a metric \( \mathbf{G} = (g_{ij}) \) such that \( \mathbf{G} \rightarrow \mathbf{G}^{\prime} = \mathbf{A}^{\mathrm{T}} \mathbf{G} \mathbf{A} \) , or equivalently \( g_{ij} \rightarrow g^{\prime}_{kl} = g_{ij} a_{ik} a_{jl} \) . More...
 
int XLALInverseTransformMetric (gsl_matrix **p_gpr_ij, const gsl_matrix *transform, const gsl_matrix *g_ij)
 Apply the inverse of a transform \( \mathbf{A}^{-1} = \mathbf{B} = (b_{ij}) \) to a metric \( \mathbf{G} = (g_{ij}) \) such that \( \mathbf{G} \rightarrow \mathbf{G}^{\prime} = \mathbf{B}^{\mathrm{T}} \mathbf{G} \mathbf{B} \) , or equivalently \( g_{ij} \rightarrow g^{\prime}_{kl} = g_{ij} b_{ik} b_{jl} \) . More...
 
int XLALDiagNormalizeMetric (gsl_matrix **p_gpr_ij, gsl_matrix **p_transform, const gsl_matrix *g_ij)
 Diagonally-normalize a metric \( g_{ij} \) . More...
 
int XLALNaturalizeMetric (gsl_matrix **p_gpr_ij, gsl_matrix **p_transform, const gsl_matrix *g_ij, const DopplerMetricParams *metricParams)
 Return a metric in naturalized coordinates. More...
 
int XLALChangeMetricReferenceTime (gsl_matrix **p_gpr_ij, gsl_matrix **p_transform, const gsl_matrix *g_ij, const DopplerCoordinateSystem *coordSys, const double Dtau)
 Compute the transform which changes the metric reference time \( \tau_0 \rightarrow \tau_1 = \tau_0 + \Delta\tau \) . More...
 

Function Documentation

◆ XLALCompareMetrics()

REAL8 XLALCompareMetrics ( const gsl_matrix *  g1_ij,
const gsl_matrix *  g2_ij 
)

Flexible comparison function between two metrics \( g^1_{ij} \) and \( g^2_{ij} \) .

Returns maximal relative deviation \( \max_{i,j} \delta g_{ij} \) , measured in terms of diagonal entries, i.e. \( \delta g_{ij} = ( g^1_{ij} - g^2_{ij} ) / \sqrt{ g^1_{ii} g^1_{jj} } \) . This should always be well-defined, as we deal with positive-definite square matrices.

Parameters
[in]g1_ijMetric to compare \( g^1_{ij} \) .
[in]g2_ijMetric to compare \( g^2_{ij} \) .

Definition at line 52 of file MetricUtils.c.

◆ XLALMetricDeterminant()

double XLALMetricDeterminant ( const gsl_matrix *  g_ij)

Compute the determinant of a metric \( g_{ij} \) .

Parameters
[in]g_ijParameter-space metric \( g_{ij} \)

Definition at line 98 of file MetricUtils.c.

◆ XLALMetricEllipseBoundingBox()

gsl_vector * XLALMetricEllipseBoundingBox ( const gsl_matrix *  g_ij,
const double  max_mismatch 
)

Compute the extent of the bounding box of the mismatch ellipse of a metric \( g_{ij} \) .

Parameters
[in]g_ijParameter-space metric \( g_{ij} \)
[in]max_mismatchMaximum prescribed mismatch

Definition at line 131 of file MetricUtils.c.

◆ XLALProjectMetric()

int XLALProjectMetric ( gsl_matrix **  p_gpr_ij,
const gsl_matrix *  g_ij,
const UINT4  c 
)

Calculate the projected metric onto the subspace orthogonal to coordinate-axis c, namely \( g^{\prime}_{ij} = g_{ij} - ( g_{ic} g_{jc} / g_{cc} ) \) , where \( c \) is the index of the projected coordinate.

Note
*p_gpr_ij will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
[in,out]p_gpr_ijPointer to projected matrix \( g^{\prime}_{ij} \)
[in]g_ijMatrix to project \( g_{ij} \)
[in]cIndex of projected coordinate

Definition at line 182 of file MetricUtils.c.

◆ XLALCholeskyLDLTDecompMetric()

int XLALCholeskyLDLTDecompMetric ( gsl_matrix **  p_cholesky,
const gsl_matrix *  g_ij 
)

Decompose a metric \( \mathbf{G} \) as \( \mathbf{G} = \mathbf{L} \mathbf{D} \mathbf{L}^{\mathrm{T}} \) , where \( \mathbf{L} \) is a lower-triangular matrix with unit diagonal, and \( \mathbf{D} \) is a diagonal matrix.

This decomposition may be useful if the metric cannot yet be guaranteed to be positive definite.

Parameters
[in,out]p_choleskyPointer to decomposition; stores \( \mathbf{L} \) in lower triangular part \( \mathbf{D} \) on diagonal
[in]g_ijMatrix to decompose \( \mathbf{G} \)

Definition at line 231 of file MetricUtils.c.

◆ XLALTransformMetric()

int XLALTransformMetric ( gsl_matrix **  p_gpr_ij,
const gsl_matrix *  transform,
const gsl_matrix *  g_ij 
)

Apply a transform \( \mathbf{A} = (a_{ij}) \) to a metric \( \mathbf{G} = (g_{ij}) \) such that \( \mathbf{G} \rightarrow \mathbf{G}^{\prime} = \mathbf{A}^{\mathrm{T}} \mathbf{G} \mathbf{A} \) , or equivalently \( g_{ij} \rightarrow g^{\prime}_{kl} = g_{ij} a_{ik} a_{jl} \) .

Note
*p_gpr_ij will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
[in,out]p_gpr_ijPointer to transformed matrix \( \mathbf{G}^{\prime} \)
[in]transformTransform to apply \( \mathbf{A} \)
[in]g_ijMatrix to transform \( \mathbf{G} \)

Definition at line 286 of file MetricUtils.c.

◆ XLALInverseTransformMetric()

int XLALInverseTransformMetric ( gsl_matrix **  p_gpr_ij,
const gsl_matrix *  transform,
const gsl_matrix *  g_ij 
)

Apply the inverse of a transform \( \mathbf{A}^{-1} = \mathbf{B} = (b_{ij}) \) to a metric \( \mathbf{G} = (g_{ij}) \) such that \( \mathbf{G} \rightarrow \mathbf{G}^{\prime} = \mathbf{B}^{\mathrm{T}} \mathbf{G} \mathbf{B} \) , or equivalently \( g_{ij} \rightarrow g^{\prime}_{kl} = g_{ij} b_{ik} b_{jl} \) .

Note
*p_gpr_ij will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
[in,out]p_gpr_ijPointer to transformed matrix \( \mathbf{G}^{\prime} \)
[in]transformTransform \( \mathbf{A} \) , the inverse of which to apply
[in]g_ijMatrix to transform \( \mathbf{G} \)

Definition at line 341 of file MetricUtils.c.

◆ XLALDiagNormalizeMetric()

int XLALDiagNormalizeMetric ( gsl_matrix **  p_gpr_ij,
gsl_matrix **  p_transform,
const gsl_matrix *  g_ij 
)

Diagonally-normalize a metric \( g_{ij} \) .

Diagonally-normalize means normalize metric by its diagonal, namely apply the transformation \( g_{ij} \rightarrow g^{\prime}_{ij} = g_{ij} / \sqrt{g_{ii} g_{jj}} \) , resulting in a lower condition number and unit diagonal elements.

If p_transform is non-NULL, return the diagonal-normalization transform in *p_transform. If p_gpr_ij is non-NULL, apply the transform to the metric *p_gpr_ij.

Note
*p_gpr_ij and *p_transform will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
p_gpr_ij[in,out,optional] Pointer to transformed matrix \( g^{\prime}_{ij} \)
p_transform[in,out,optional] Pointer to diagonal-normalization transform
[in]g_ijMatrix to transform \( g_{ij} \)

Definition at line 391 of file MetricUtils.c.

◆ XLALNaturalizeMetric()

int XLALNaturalizeMetric ( gsl_matrix **  p_gpr_ij,
gsl_matrix **  p_transform,
const gsl_matrix *  g_ij,
const DopplerMetricParams metricParams 
)

Return a metric in naturalized coordinates.

Frequency coordinates of spindown order \( s \) are scaled by

\[ \frac{2\pi}{(s+1)!} \left(\frac{\overline{\Delta T}}{2}\right)^{s+1} \]

where \( \overline{\Delta T}\equiv\sum_{k}^{N} \Delta T_k \) is the average segment-length over all \( N \) segments.

Sky coordinates are scaled by Holgers' units, see Eq.(44) in PRD82,042002(2010), without the equatorial rotation in alpha:

\[ \frac{2\pi f R_{E} \cos(\delta_{D})}{c} \]

where \( f \) is the frequency and \( R_{E} \) the Earth radius, and \( \delta_{D} \) is the detectors latitude.

If p_transform is non-NULL, return the naturalization transform in *p_transform. If p_gpr_ij is non-NULL, apply the transform to the metric *p_gpr_ij.

Note
*p_gpr_ij and *p_transform will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
p_gpr_ij[in,out,optional] Pointer to transformed matrix \( g^{\prime}_{ij} \)
p_transform[in,out,optional] Pointer to naturalization transform
[in]g_ijMatrix to transform \( g_{ij} \)
[in]metricParamsInput parameters used to calculate naturalization transform

Definition at line 456 of file MetricUtils.c.

◆ XLALChangeMetricReferenceTime()

int XLALChangeMetricReferenceTime ( gsl_matrix **  p_gpr_ij,
gsl_matrix **  p_transform,
const gsl_matrix *  g_ij,
const DopplerCoordinateSystem coordSys,
const double  Dtau 
)

Compute the transform which changes the metric reference time \( \tau_0 \rightarrow \tau_1 = \tau_0 + \Delta\tau \) .

If p_transform is non-NULL, return the reference-time transform in *p_transform. If p_gpr_ij is non-NULL, apply the transform to the metric *p_gpr_ij.

Note
*p_gpr_ij and *p_transform will be allocated if NULL. *p_gpr_ij and g_ij may point to the same matrix.
Parameters
p_gpr_ij[in,out,optional] Pointer to transformed matrix \( g^{\prime}_{ij} \)
p_transform[in,out,optional] Pointer to reference time transform
[in]g_ijMatrix to transform \( g_{ij} \)
[in]coordSysCoordinate system of metric
[in]DtauDifference between new and old reference times \( \Delta\tau \)

Definition at line 578 of file MetricUtils.c.