20 #ifndef _METRICUTILS_H
21 #define _METRICUTILS_H
23 #include <gsl/gsl_vector.h>
24 #include <gsl/gsl_matrix.h>
25 #include <lal/LALStdlib.h>
26 #include <lal/UniversalDopplerMetric.h>
49 int XLALTransformMetric( gsl_matrix **p_gpr_ij,
const gsl_matrix *transform,
const gsl_matrix *g_ij );
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 .
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 ,...
REAL8 XLALCompareMetrics(const gsl_matrix *g1_ij, const gsl_matrix *g2_ij)
Flexible comparison function between two metrics and .
int XLALCholeskyLDLTDecompMetric(gsl_matrix **p_cholesky, const gsl_matrix *g_ij)
Decompose a metric as , where is a lower-triangular matrix with unit diagonal, and is a diagonal ...
int XLALInverseTransformMetric(gsl_matrix **p_gpr_ij, const gsl_matrix *transform, const gsl_matrix *g_ij)
Apply the inverse of a transform to a metric such that , or equivalently .
int XLALDiagNormalizeMetric(gsl_matrix **p_gpr_ij, gsl_matrix **p_transform, const gsl_matrix *g_ij)
Diagonally-normalize a metric .
double XLALMetricDeterminant(const gsl_matrix *g_ij)
Compute the determinant of a metric .
int XLALTransformMetric(gsl_matrix **p_gpr_ij, const gsl_matrix *transform, const gsl_matrix *g_ij)
Apply a transform to a metric such that , or equivalently .
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 .
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.
type describing a Doppler coordinate system: lists the number of dimensions and the symbolic names of...
meta-info specifying a Doppler-metric