LALPulsar  6.1.0.1-c9a8ef6
Header UniversalDopplerMetric.h

Detailed Description

Author
Reinhard Prix, Karl Wette

Function to compute the full F-statistic metric, including antenna-pattern functions from multi-detector, derived in [20] .

Prototypes

static double XLALAverage_am1_am2_Phi_i_Phi_j (const intparams_t *params, double *relerr_max)
 Integrate a general quadruple product CW_am1_am2_Phi_i_Phi_j() from 0 to 1. More...
 
static double CW_am1_am2_Phi_i_Phi_j (double tt, void *params)
 For gsl-integration: general quadruple product between two antenna-pattern functions am1, am2 in {a(t),b(t)} and two phase-derivatives phi_i * phi_j, i.e. More...
 
REAL8 XLALComputePhaseDerivative (REAL8 t, const PulsarDopplerParams *dopplerPoint, DopplerCoordinateID coordID, const EphemerisData *edat, const LALDetector *site, BOOLEAN includeRoemer)
 
static double CW_Phi_i (double tt, void *params)
 Partial derivative of continuous-wave (CW) phase, with respect to Doppler coordinate 'i' := intparams_t->phderiv. More...
 
int XLALDetectorPosVel (PosVel3D_t *spin_posvel, PosVel3D_t *orbit_posvel, const LIGOTimeGPS *tGPS, const LALDetector *site, const EphemerisData *edat, DetectorMotionType detMotionType)
 Given a GPS time and detector, return the current position (and velocity) of the detector. More...
 
int XLALPtolemaicPosVel (PosVel3D_t *posvel, const LIGOTimeGPS *tGPS)
 Compute position and velocity assuming a purely "Ptolemaic" orbital motion (i.e. More...
 
static double XLALCovariance_Phi_ij (const MultiLALDetector *multiIFO, const MultiNoiseFloor *multiNoiseFloor, const LALSegList *segList, const intparams_t *params, double *relerr_max)
 Compute a pure phase-deriv covariance \( [\phi_i, \phi_j] = \langle phi_i phi_j\rangle - \langle phi_i\rangle\langle phi_j\rangle \) which gives a component of the "phase metric". More...
 
DopplerPhaseMetricXLALComputeDopplerPhaseMetric (const DopplerMetricParams *metricParams, const EphemerisData *edat)
 Calculate an approximate "phase-metric" with the specified parameters. More...
 
void XLALDestroyDopplerPhaseMetric (DopplerPhaseMetric *metric)
 Free a DopplerPhaseMetric structure. More...
 
DopplerFstatMetricXLALComputeDopplerFstatMetric (const DopplerMetricParams *metricParams, const EphemerisData *edat)
 Calculate the general (single-segment coherent, or multi-segment semi-coherent) full (multi-IFO) Fstat-metrix and the Fisher-matrix derived in [20] . More...
 
void XLALDestroyDopplerFstatMetric (DopplerFstatMetric *metric)
 Free a DopplerFstatMetric structure. More...
 
FmetricAtoms_tXLALComputeAtomsForFmetric (const DopplerMetricParams *metricParams, const EphemerisData *edat)
 Function to the compute the FmetricAtoms_t, from which the F-metric and Fisher-matrix can be computed. More...
 
int XLALParseDetectorMotionString (const CHAR *detMotionString)
 Parse a detector-motion type string into the corresponding enum-number,. More...
 
const CHARXLALDetectorMotionName (DetectorMotionType detMotionType)
 Provide a pointer to a static string containing the DopplerCoordinate-name cooresponding to the enum DopplerCoordinateID. More...
 
int XLALParseDopplerCoordinateString (const CHAR *coordName)
 Parse a DopplerCoordinate-name into the corresponding DopplerCoordinateID. More...
 
int XLALDopplerCoordinateNames2System (DopplerCoordinateSystem *coordSys, const LALStringVector *coordNames)
 Given a LALStringVector of coordinate-names, parse them into a 'DopplerCoordinateSystem', namely a list of coordinate-IDs. More...
 
int XLALFindDopplerCoordinateInSystem (const DopplerCoordinateSystem *coordSys, const DopplerCoordinateID coordID)
 Given a coordinate ID 'coordID', return its dimension within the given coordinate system 'coordSys', or return -1 if 'coordID' is not found. More...
 
const CHARXLALDopplerCoordinateName (DopplerCoordinateID coordID)
 Provide a pointer to a static string containing the DopplerCoordinate-name cooresponding to the enum DopplerCoordinateID. More...
 
const CHARXLALDopplerCoordinateHelp (DopplerCoordinateID coordID)
 Provide a pointer to a static string containing the a descriptive 'help-string' describing the coordinate DopplerCoordinateID. More...
 
CHARXLALDopplerCoordinateHelpAll (void)
 Return a string (allocated here) containing a full name - helpstring listing for all doppler-coordinates DopplerCoordinateID allowed by UniversalDopplerMetric.c. More...
 
void XLALDestroyFmetricAtoms (FmetricAtoms_t *atoms)
 Free a FmetricAtoms_t structure, allowing any pointers to be NULL. More...
 
static DopplerFstatMetricXLALComputeFmetricFromAtoms (const FmetricAtoms_t *atoms, REAL8 cosi, REAL8 psi)
 Compute the 'F-metric' gF_ij (and also gFav_ij, m1_ij, m2_ij, m3_ij) from the given FmetricAtoms and the signal amplitude parameters. More...
 
static gsl_matrix * XLALComputeFisherFromAtoms (const FmetricAtoms_t *atoms, PulsarAmplitudeParams Amp)
 Function to compute full 4+n dimensional Fisher matric for the full CW parameter-space of Amplitude + Doppler parameters ! More...
 
void XLALequatorialVect2ecliptic (vect3D_t out, const vect3D_t in)
 Convert 3-D vector from equatorial into ecliptic coordinates. More...
 
void XLALeclipticVect2equatorial (vect3D_t out, const vect3D_t in)
 Convert 3-D vector from ecliptic into equatorial coordinates. More...
 
void XLALmatrix33_in_vect3 (vect3D_t out, mat33_t mat, const vect3D_t in)
 compute matrix product mat . More...
 
vect3Dlist_tXLALComputeOrbitalDerivatives (UINT4 maxorder, const LIGOTimeGPS *tGPS, const EphemerisData *edat)
 Compute time-derivatives up to 'maxorder' of the Earths' orbital position vector \( r_{\mathrm{orb}}(t) \) . More...
 
void XLALDestroyVect3Dlist (vect3Dlist_t *list)
 
static UINT4 findHighestGCSpinOrder (const DopplerCoordinateSystem *coordSys)
 Return the highest 'global-correlation' spindown order found in this coordinate system. More...
 

Data Structures

struct  vect2Dlist_t
 variable-length list of 2D-vectors More...
 
struct  vect3Dlist_t
 variable-length list of 3D vectors More...
 
struct  PosVel3D_t
 Small Container to hold two 3D vectors: position and velocity. More...
 
struct  DopplerCoordinateSystem
 type describing a Doppler coordinate system: lists the number of dimensions and the symbolic names of the coordinates. More...
 
struct  DopplerMetricParams
 meta-info specifying a Doppler-metric More...
 
struct  FmetricAtoms_t
 Struct to hold the 'atoms', ie weighted phase-derivative averages like \( \langle a^2 \partial_i \phi \partial_j \phi\rangle> \) from which the F-metric is computed, but also the full Fisher-matrix. More...
 
struct  DopplerFstatMetric
 Struct to hold the output of XLALComputeDopplerFstatMetric(), including meta-info on the number of dimensions, the coordinate-system and type of metric. More...
 
struct  DopplerPhaseMetric
 Struct to hold the output of XLALComputeDopplerPhaseMetric(), including meta-info on the number of dimensions, the coordinate-system and type of metric. More...
 

Typedefs

typedef REAL8 vect2D_t[2]
 2D vector More...
 
typedef REAL8 vect3D_t[3]
 3D vector More...
 
typedef REAL8 mat33_t[3][3]
 3x3 matrix, useful for spatial 3D vector operations More...
 

Enumerations

enum  DetectorMotionType {
  DETMOTION_SPIN = 0x01 , DETMOTION_SPINZ = 0x02 , DETMOTION_SPINXY = 0x03 , DETMOTION_MASKSPIN = 0x0F ,
  DETMOTION_ORBIT = 0x10 , DETMOTION_PTOLEORBIT = 0x20 , DETMOTION_MASKORBIT = 0xF0
}
 Bitfield of different types of detector-motion to use in order to compute the Doppler-metric. More...
 
enum  DopplerCoordinateID {
  DOPPLERCOORD_NONE = -1 , DOPPLERCOORD_FREQ , DOPPLERCOORD_F1DOT , DOPPLERCOORD_F2DOT ,
  DOPPLERCOORD_F3DOT , DOPPLERCOORD_F4DOT , DOPPLERCOORD_LASTFDOT = DOPPLERCOORD_F4DOT , DOPPLERCOORD_GC_NU0 ,
  DOPPLERCOORD_GC_NU1 , DOPPLERCOORD_GC_NU2 , DOPPLERCOORD_GC_NU3 , DOPPLERCOORD_GC_NU4 ,
  DOPPLERCOORD_ALPHA , DOPPLERCOORD_DELTA , DOPPLERCOORD_N2X_EQU , DOPPLERCOORD_N2Y_EQU ,
  DOPPLERCOORD_N2X_ECL , DOPPLERCOORD_N2Y_ECL , DOPPLERCOORD_N3X_EQU , DOPPLERCOORD_N3Y_EQU ,
  DOPPLERCOORD_N3Z_EQU , DOPPLERCOORD_N3X_ECL , DOPPLERCOORD_N3Y_ECL , DOPPLERCOORD_N3Z_ECL ,
  DOPPLERCOORD_N3SX_EQU , DOPPLERCOORD_N3SY_EQU , DOPPLERCOORD_N3OX_ECL , DOPPLERCOORD_N3OY_ECL ,
  DOPPLERCOORD_N3OZ_ECL , DOPPLERCOORD_ASINI , DOPPLERCOORD_TASC , DOPPLERCOORD_PORB ,
  DOPPLERCOORD_KAPPA , DOPPLERCOORD_ETA , DOPPLERCOORD_VP , DOPPLERCOORD_DASC ,
  DOPPLERCOORD_KAPPAP , DOPPLERCOORD_ETAP , DOPPLERCOORD_LAST
}
 enum listing symbolic 'names' for all Doppler Coordinates supported by the metric codes in FstatMetric More...
 

Macros

#define POW2(a)   ( (a) * (a) )
 Function to compute the full F-statistic metric, including antenna-pattern functions from multi-detector, as derived in [20] . More...
 
#define DOPPLERMETRIC_MAX_DIM   60
 should be large enough for a long time ... More...
 

Files

file  old-FstatMetric.c
 The only purpose of this file is to serve as a backwards-comparison check for XLALDopplerFstatMetric(). This used to be a standalone-code 'lalapps_FstatMetric', and was XLALified and moved into the test-directory, main() was wrapped into the forwards-compatible function XLALOldDopplerFstatMetric() and called in UniversalDopplerMetricTest for comparison.
 
file  UniversalDopplerMetricTest.c
 Tests for exported functions in UniversalDopplerMetric.
 

Function Documentation

◆ XLALAverage_am1_am2_Phi_i_Phi_j()

static double XLALAverage_am1_am2_Phi_i_Phi_j ( const intparams_t params,
double *  relerr_max 
)
static

Integrate a general quadruple product CW_am1_am2_Phi_i_Phi_j() from 0 to 1.

This implements the expression \( \langle<q_1 q_2 \phi_i \phi_j\rangle \) for single-IFO average over the observation time.

The input parameters correspond to CW_am1_am2_Phi_i_Phi_j()

Definition at line 285 of file UniversalDopplerMetric.c.

◆ CW_am1_am2_Phi_i_Phi_j()

static double CW_am1_am2_Phi_i_Phi_j ( double  tt,
void *  params 
)
static

For gsl-integration: general quadruple product between two antenna-pattern functions am1, am2 in {a(t),b(t)} and two phase-derivatives phi_i * phi_j, i.e.

compute an expression of the form \( q_1(t) q_2(t) \phi_i(t) \phi_j(t) \) , where \( q_i = \{a(t), b(t)\} \) .

NOTE: this can be 'truncated' to any sub-expression by using AMCOMP_NONE for antenna-pattern component and DOPPLERCOORD_NONE for DopplerCoordinate, eg in this way this function can be used to compute \( a^2(t), b^2(t), a(t) b(t) \) , or \( phi_i(t) phi_j(t) \) .

Definition at line 359 of file UniversalDopplerMetric.c.

◆ XLALComputePhaseDerivative()

REAL8 XLALComputePhaseDerivative ( REAL8  t,
const PulsarDopplerParams dopplerPoint,
DopplerCoordinateID  coordID,
const EphemerisData edat,
const LALDetector site,
BOOLEAN  includeRoemer 
)
Parameters
ttime 't' to compute derivative at: (effectively interpreted as SSB time if approxPhase given)
dopplerPointphase-evolution ('doppler') parameters to compute phase-derivative at
coordIDcoord-ID of coordinate wrt which to compute phase derivative
edatephemeris data
siteoptional detector (if NULL: use H1 as default)
includeRoemerwhether to include Roemer-delay correction 'tSSB = t + dRoemer(t)' or not

Definition at line 437 of file UniversalDopplerMetric.c.

◆ CW_Phi_i()

static double CW_Phi_i ( double  tt,
void *  params 
)
static

Partial derivative of continuous-wave (CW) phase, with respect to Doppler coordinate 'i' := intparams_t->phderiv.

Time is in 'natural units' of Tspan, i.e. tt is in [0, 1] corresponding to GPS-times in [startTime, startTime + Tspan ]

< Frequency [Units: Hz].

< Global correlation frequency [Units: Hz]. Activates 'reduced' detector position.

< First spindown [Units: Hz/s].

< Global correlation first spindown [Units: Hz/s]. Activates 'reduced' detector position.

< Second spindown [Units: Hz/s^2].

< Global correlation second spindown [Units: Hz/s^2]. Activates 'reduced' detector position.

< Third spindown [Units: Hz/s^3].

< Global correlation third spindown [Units: Hz/s^3]. Activates 'reduced' detector position.

< Fourth spindown [Units: Hz/s^4].

< Global correlation fourth spindown [Units: Hz/s^4]. Activates 'reduced' detector position.

< Right ascension [Units: radians]. Uses 'reduced' detector position.

< Declination [Units: radians]. Uses 'reduced' detector position.

< X component of constrained sky position in equatorial coordinates [Units: none]. Uses 'reduced' detector position.

< Y component of constrained sky position in equatorial coordinates [Units: none]. Uses 'reduced' detector position.

< X component of constrained sky position in ecliptic coordinates [Units: none]. Uses 'reduced' detector position.

< Y component of constrained sky position in ecliptic coordinates [Units: none]. Uses 'reduced' detector position.

< X component of unconstrained super-sky position in equatorial coordinates [Units: none].

< Y component of unconstrained super-sky position in equatorial coordinates [Units: none].

< Z component of unconstrained super-sky position in equatorial coordinates [Units: none].

< X component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< Y component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< Z component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< X spin-component of unconstrained super-sky position in equatorial coordinates [Units: none].

< Y spin-component of unconstrained super-sky position in equatorial coordinates [Units: none].

< X orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< Y orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< Z orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

< Projected semimajor axis of binary orbit in small-eccentricy limit (ELL1 model) [Units: light seconds].

< Time of ascension (neutron star crosses line of nodes moving away from observer) for binary orbit (ELL1 model) [Units: GPS s].

< Period of binary orbit (ELL1 model) [Units: s].

< Lagrange parameter 'kappa = ecc * cos(argp)', ('ecc' = eccentricity, 'argp' = argument of periapse) of binary orbit (ELL1 model) [Units: none]

< Lagrange parameter 'eta = ecc * sin(argp) of binary orbit (ELL1 model) [Units: none]

< Distance traversed on the arc of binary orbit (ELL1 model) 'dasc = 2 * pi * (ap/porb) * tasc' [Units: light second]."

< Rescaled (by asini) differential-coordinate 'dvp = asini * dOMEGA', ('OMEGA' = 2 * pi/'porb') of binary orbit (ELL1 model) [Units: (light second)/(GPS second)].

< Rescaled (by asini) differential-coordinate 'dkappap = asini * dkappa' [Units: light seconds].

< Rescaled (by asini) differential-coordinate 'detap = asini * deta' [Units: light seconds].

Definition at line 483 of file UniversalDopplerMetric.c.

◆ XLALDetectorPosVel()

int XLALDetectorPosVel ( PosVel3D_t spin_posvel,
PosVel3D_t orbit_posvel,
const LIGOTimeGPS tGPS,
const LALDetector site,
const EphemerisData edat,
DetectorMotionType  detMotionType 
)

Given a GPS time and detector, return the current position (and velocity) of the detector.

< No spin motion

< Full spin motion

< Ecliptic-Z component of spin motion only

< Ecliptic-X+Y components of spin motion only

< No orbital motion

< Ephemeris-based orbital motion

< Ptolemaic (circular) orbital motion

Parameters
[out]spin_posvelinstantaneous sidereal position and velocity vector
[out]orbit_posvelinstantaneous orbital position and velocity vector
[in]tGPSGPS time
[in]sitedetector info
[in]edatephemeris data
[in]detMotionTypedetector motion type

Definition at line 761 of file UniversalDopplerMetric.c.

◆ XLALPtolemaicPosVel()

int XLALPtolemaicPosVel ( PosVel3D_t posvel,
const LIGOTimeGPS tGPS 
)

Compute position and velocity assuming a purely "Ptolemaic" orbital motion (i.e.

on a circle) around the sun, approximating Earth's orbit

Parameters
[out]posvelinstantaneous position and velocity vector
[in]tGPSGPS time

Definition at line 886 of file UniversalDopplerMetric.c.

◆ XLALCovariance_Phi_ij()

static double XLALCovariance_Phi_ij ( const MultiLALDetector multiIFO,
const MultiNoiseFloor multiNoiseFloor,
const LALSegList segList,
const intparams_t params,
double *  relerr_max 
)
static

Compute a pure phase-deriv covariance \( [\phi_i, \phi_j] = \langle phi_i phi_j\rangle - \langle phi_i\rangle\langle phi_j\rangle \) which gives a component of the "phase metric".

NOTE: for passing unit noise-weights, set MultiNoiseFloor->length=0 (but multiNoiseFloor==NULL is invalid)

Parameters
[in]multiIFOdetectors to use
[in]multiNoiseFloorcorresponding noise floors for weights, NULL means unit-weights
[in]segListsegment list
[in]paramsintegration parameters
[in]relerr_maxmaximal error for integration

Definition at line 928 of file UniversalDopplerMetric.c.

◆ XLALComputeDopplerPhaseMetric()

DopplerPhaseMetric * XLALComputeDopplerPhaseMetric ( const DopplerMetricParams metricParams,
const EphemerisData edat 
)

Calculate an approximate "phase-metric" with the specified parameters.

Note: if this function is called with multiple detectors, the phase components are averaged over detectors as well as time. This is a somewhat ad-hoc approach; if you want a more rigorous multi-detector metric you need to use the full Fstat-metric, as computed by XLALComputeDopplerFstatMetric().

Return NULL on error.

Parameters
metricParamsinput parameters determining the metric calculation
edatephemeris data

Definition at line 1202 of file UniversalDopplerMetric.c.

◆ XLALDestroyDopplerPhaseMetric()

void XLALDestroyDopplerPhaseMetric ( DopplerPhaseMetric metric)

Free a DopplerPhaseMetric structure.

Definition at line 1427 of file UniversalDopplerMetric.c.

◆ XLALComputeDopplerFstatMetric()

DopplerFstatMetric * XLALComputeDopplerFstatMetric ( const DopplerMetricParams metricParams,
const EphemerisData edat 
)

Calculate the general (single-segment coherent, or multi-segment semi-coherent) full (multi-IFO) Fstat-metrix and the Fisher-matrix derived in [20] .

The semi-coherent metrics \( g_{ij} \) over \( N \) segments are computed according to

\[ \overline{g}_{ij} \equiv \frac{1}{N} \sum_{k=1}^{N} g_{ij,k} \]

where \( g_{ij,k} \) is the coherent single-segment metric of segment k

Note: The returned DopplerFstatMetric struct contains the matrices g_ij (the phase metric), gF_ij (the F-metric), gFav_ij (the average F-metric), m1_ij, m2_ij, m3_ij (auxiliary matrices) and Fisher_ab (the full 4+n dimensional Fisher matrix).

The returned metric struct also carries the meta-info about the metrics in the field 'DopplerMetricParams meta'.

Note2: for backwards-compatibility, we treat params->Nseg==0 equivalent to Nseg==1, ie compute the coherent, single-segment metrics

Return NULL on error.

Parameters
metricParamsinput parameters determining the metric calculation
edatephemeris data

Definition at line 1466 of file UniversalDopplerMetric.c.

◆ XLALDestroyDopplerFstatMetric()

void XLALDestroyDopplerFstatMetric ( DopplerFstatMetric metric)

Free a DopplerFstatMetric structure.

Definition at line 1606 of file UniversalDopplerMetric.c.

◆ XLALComputeAtomsForFmetric()

FmetricAtoms_t * XLALComputeAtomsForFmetric ( const DopplerMetricParams metricParams,
const EphemerisData edat 
)

Function to the compute the FmetricAtoms_t, from which the F-metric and Fisher-matrix can be computed.

NOTE: if MultiNoiseFloor.length=0, unit noise weights are assumed.

Parameters
metricParamsinput parameters determining the metric calculation
edatephemeris data

Definition at line 1632 of file UniversalDopplerMetric.c.

◆ XLALParseDetectorMotionString()

int XLALParseDetectorMotionString ( const CHAR detMotionString)

Parse a detector-motion type string into the corresponding enum-number,.

Definition at line 1901 of file UniversalDopplerMetric.c.

◆ XLALDetectorMotionName()

const CHAR * XLALDetectorMotionName ( DetectorMotionType  detMotionType)

Provide a pointer to a static string containing the DopplerCoordinate-name cooresponding to the enum DopplerCoordinateID.

Definition at line 1923 of file UniversalDopplerMetric.c.

◆ XLALParseDopplerCoordinateString()

int XLALParseDopplerCoordinateString ( const CHAR coordName)

Parse a DopplerCoordinate-name into the corresponding DopplerCoordinateID.

Definition at line 1943 of file UniversalDopplerMetric.c.

◆ XLALDopplerCoordinateNames2System()

int XLALDopplerCoordinateNames2System ( DopplerCoordinateSystem coordSys,
const LALStringVector coordNames 
)

Given a LALStringVector of coordinate-names, parse them into a 'DopplerCoordinateSystem', namely a list of coordinate-IDs.

Parameters
[out]coordSyspointer to coord-system struct
[in]coordNameslist of coordinate names

Definition at line 1968 of file UniversalDopplerMetric.c.

◆ XLALFindDopplerCoordinateInSystem()

int XLALFindDopplerCoordinateInSystem ( const DopplerCoordinateSystem coordSys,
const DopplerCoordinateID  coordID 
)

Given a coordinate ID 'coordID', return its dimension within the given coordinate system 'coordSys', or return -1 if 'coordID' is not found.

Definition at line 1996 of file UniversalDopplerMetric.c.

◆ XLALDopplerCoordinateName()

const CHAR * XLALDopplerCoordinateName ( DopplerCoordinateID  coordID)

Provide a pointer to a static string containing the DopplerCoordinate-name cooresponding to the enum DopplerCoordinateID.

Definition at line 2013 of file UniversalDopplerMetric.c.

◆ XLALDopplerCoordinateHelp()

const CHAR * XLALDopplerCoordinateHelp ( DopplerCoordinateID  coordID)

Provide a pointer to a static string containing the a descriptive 'help-string' describing the coordinate DopplerCoordinateID.

Definition at line 2033 of file UniversalDopplerMetric.c.

◆ XLALDopplerCoordinateHelpAll()

CHAR * XLALDopplerCoordinateHelpAll ( void  )

Return a string (allocated here) containing a full name - helpstring listing for all doppler-coordinates DopplerCoordinateID allowed by UniversalDopplerMetric.c.

Definition at line 2052 of file UniversalDopplerMetric.c.

◆ XLALDestroyFmetricAtoms()

void XLALDestroyFmetricAtoms ( FmetricAtoms_t atoms)

Free a FmetricAtoms_t structure, allowing any pointers to be NULL.

Definition at line 2104 of file UniversalDopplerMetric.c.

◆ XLALComputeFmetricFromAtoms()

static DopplerFstatMetric * XLALComputeFmetricFromAtoms ( const FmetricAtoms_t atoms,
REAL8  cosi,
REAL8  psi 
)
static

Compute the 'F-metric' gF_ij (and also gFav_ij, m1_ij, m2_ij, m3_ij) from the given FmetricAtoms and the signal amplitude parameters.

Definition at line 2143 of file UniversalDopplerMetric.c.

◆ XLALComputeFisherFromAtoms()

static gsl_matrix * XLALComputeFisherFromAtoms ( const FmetricAtoms_t atoms,
PulsarAmplitudeParams  Amp 
)
static

Function to compute full 4+n dimensional Fisher matric for the full CW parameter-space of Amplitude + Doppler parameters !

Definition at line 2283 of file UniversalDopplerMetric.c.

◆ XLALequatorialVect2ecliptic()

void XLALequatorialVect2ecliptic ( vect3D_t  out,
const vect3D_t  in 
)

Convert 3-D vector from equatorial into ecliptic coordinates.

Definition at line 2395 of file UniversalDopplerMetric.c.

◆ XLALeclipticVect2equatorial()

void XLALeclipticVect2equatorial ( vect3D_t  out,
const vect3D_t  in 
)

Convert 3-D vector from ecliptic into equatorial coordinates.

Definition at line 2408 of file UniversalDopplerMetric.c.

◆ XLALmatrix33_in_vect3()

void XLALmatrix33_in_vect3 ( vect3D_t  out,
mat33_t  mat,
const vect3D_t  in 
)

compute matrix product mat .

vect

Definition at line 2421 of file UniversalDopplerMetric.c.

◆ XLALComputeOrbitalDerivatives()

vect3Dlist_t * XLALComputeOrbitalDerivatives ( UINT4  maxorder,
const LIGOTimeGPS tGPS,
const EphemerisData edat 
)

Compute time-derivatives up to 'maxorder' of the Earths' orbital position vector \( r_{\mathrm{orb}}(t) \) .

Algorithm: using 5-point differentiation expressions on r_orb(t) returned from XLALBarycenterEarth().

Returns a vector of derivatives \( \frac{d^n\,r_{\mathrm{orb}}}{d\,t^n} \) at the given GPS time. Note, the return vector includes the zeroth-order derivative, so we return (maxorder + 1) derivatives: n = 0 ... maxorder

Parameters
[in]maxorderhighest derivative-order to compute
[in]tGPSGPS time at which to compute the derivatives
[in]edatephemeris data

Definition at line 2446 of file UniversalDopplerMetric.c.

◆ XLALDestroyVect3Dlist()

void XLALDestroyVect3Dlist ( vect3Dlist_t list)

Definition at line 2548 of file UniversalDopplerMetric.c.

◆ findHighestGCSpinOrder()

static UINT4 findHighestGCSpinOrder ( const DopplerCoordinateSystem coordSys)
static

Return the highest 'global-correlation' spindown order found in this coordinate system.

Counting nu0 = order1, nu1 = order2, nu2 = order3, ..., order = 0 therefore means there are no GC spin coordinates at all

Definition at line 2570 of file UniversalDopplerMetric.c.

Typedef Documentation

◆ vect2D_t

typedef REAL8 vect2D_t[2]

2D vector

Definition at line 61 of file UniversalDopplerMetric.h.

◆ vect3D_t

typedef REAL8 vect3D_t[3]

3D vector

Definition at line 63 of file UniversalDopplerMetric.h.

◆ mat33_t

typedef REAL8 mat33_t[3][3]

3x3 matrix, useful for spatial 3D vector operations

Definition at line 65 of file UniversalDopplerMetric.h.

Enumeration Type Documentation

◆ DetectorMotionType

Bitfield of different types of detector-motion to use in order to compute the Doppler-metric.

Enumerator
DETMOTION_SPIN 

Full spin motion.

DETMOTION_SPINZ 

Ecliptic-Z component of spin motion only.

DETMOTION_SPINXY 

Ecliptic-X+Y components of spin motion only.

DETMOTION_MASKSPIN 

Mask for spin motion bits.

DETMOTION_ORBIT 

Ephemeris-based orbital motion.

DETMOTION_PTOLEORBIT 

Ptolemaic (circular) orbital motion.

DETMOTION_MASKORBIT 

Mask for orbital motion bits.

Definition at line 94 of file UniversalDopplerMetric.h.

◆ DopplerCoordinateID

enum listing symbolic 'names' for all Doppler Coordinates supported by the metric codes in FstatMetric

Enumerator
DOPPLERCOORD_NONE 

No Doppler component.

DOPPLERCOORD_FREQ 

Frequency [Units: Hz].

DOPPLERCOORD_F1DOT 

First spindown [Units: Hz/s].

DOPPLERCOORD_F2DOT 

Second spindown [Units: Hz/s^2].

DOPPLERCOORD_F3DOT 

Third spindown [Units: Hz/s^3].

DOPPLERCOORD_F4DOT 

Fourth spindown [Units: Hz/s^4].

DOPPLERCOORD_LASTFDOT 
DOPPLERCOORD_GC_NU0 

Global correlation frequency [Units: Hz].

Activates 'reduced' detector position.

DOPPLERCOORD_GC_NU1 

Global correlation first spindown [Units: Hz/s].

Activates 'reduced' detector position.

DOPPLERCOORD_GC_NU2 

Global correlation second spindown [Units: Hz/s^2].

Activates 'reduced' detector position.

DOPPLERCOORD_GC_NU3 

Global correlation third spindown [Units: Hz/s^3].

Activates 'reduced' detector position.

DOPPLERCOORD_GC_NU4 

Global correlation fourth spindown [Units: Hz/s^4].

Activates 'reduced' detector position.

DOPPLERCOORD_ALPHA 

Right ascension [Units: radians].

Uses 'reduced' detector position.

DOPPLERCOORD_DELTA 

Declination [Units: radians].

Uses 'reduced' detector position.

DOPPLERCOORD_N2X_EQU 

X component of contrained sky position in equatorial coordinates [Units: none].

Uses 'reduced' detector position.

DOPPLERCOORD_N2Y_EQU 

Y component of contrained sky position in equatorial coordinates [Units: none].

Uses 'reduced' detector position.

DOPPLERCOORD_N2X_ECL 

X component of contrained sky position in ecliptic coordinates [Units: none].

Uses 'reduced' detector position.

DOPPLERCOORD_N2Y_ECL 

Y component of contrained sky position in ecliptic coordinates [Units: none].

Uses 'reduced' detector position.

DOPPLERCOORD_N3X_EQU 

X component of unconstrained super-sky position in equatorial coordinates [Units: none].

DOPPLERCOORD_N3Y_EQU 

Y component of unconstrained super-sky position in equatorial coordinates [Units: none].

DOPPLERCOORD_N3Z_EQU 

Z component of unconstrained super-sky position in equatorial coordinates [Units: none].

DOPPLERCOORD_N3X_ECL 

X component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_N3Y_ECL 

Y component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_N3Z_ECL 

Z component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_N3SX_EQU 

X spin-component of unconstrained super-sky position in equatorial coordinates [Units: none].

DOPPLERCOORD_N3SY_EQU 

Y spin-component of unconstrained super-sky position in equatorial coordinates [Units: none].

DOPPLERCOORD_N3OX_ECL 

X orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_N3OY_ECL 

Y orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_N3OZ_ECL 

Z orbit-component of unconstrained super-sky position in ecliptic coordinates [Units: none].

DOPPLERCOORD_ASINI 

Projected semimajor axis of binary orbit in small-eccentricy limit (ELL1 model) [Units: light seconds].

DOPPLERCOORD_TASC 

Time of ascension (neutron star crosses line of nodes moving away from observer) for binary orbit (ELL1 model) [Units: GPS seconds].

DOPPLERCOORD_PORB 

Period of binary orbit (ELL1 model) [Units: s].

DOPPLERCOORD_KAPPA 

Lagrange parameter 'kappa = ecc * cos(argp)', ('ecc' = eccentricity, 'argp' = argument of periapse) of binary orbit (ELL1 model) [Units: none].

DOPPLERCOORD_ETA 

Lagrange parameter 'eta = ecc * sin(argp) of binary orbit (ELL1 model) [Units: none].

DOPPLERCOORD_VP 

Rescaled (by asini) differential-coordinate 'dvp = asini * dOMEGA', ('OMEGA' = 2 * pi/'porb') of binary orbit (ELL1 model) [Units: (light second)/(GPS second)].

DOPPLERCOORD_DASC 

Distance traversed on the arc of binary orbit (ELL1 model) 'dasc = 2 * pi * (ap/porb) * tasc' [Units: light second].

DOPPLERCOORD_KAPPAP 

Rescaled (by asini) differential-coordinate 'dkappap = asini * dkappa' [Units: light seconds].

DOPPLERCOORD_ETAP 

Rescaled (by asini) differential-coordinate 'detap = asini * deta' [Units: light seconds].

DOPPLERCOORD_LAST 

Definition at line 110 of file UniversalDopplerMetric.h.

Macro Definition Documentation

◆ POW2

#define POW2 (   a)    ( (a) * (a) )

Function to compute the full F-statistic metric, including antenna-pattern functions from multi-detector, as derived in [20] .

Author
Reinhard Prix, Karl Wette shortcuts for integer powers

Definition at line 60 of file UniversalDopplerMetric.c.

◆ DOPPLERMETRIC_MAX_DIM

#define DOPPLERMETRIC_MAX_DIM   60

should be large enough for a long time ...

Definition at line 166 of file UniversalDopplerMetric.h.