/** * basic initializations: set-up 'ConfigVariables' */ int XLALInitCode ( ConfigVariables *cfg, const UserVariables_t *uvar, const char *app_name) { if ( !cfg || !uvar || !app_name ) { LogPrintf (LOG_CRITICAL, "%s: illegal NULL pointer input.\n\n", __func__ ); XLAL_ERROR (XLAL_EINVAL ); } /* Init ephemerides */ XLAL_CHECK ( (cfg->edat = XLALInitBarycenter ( uvar->ephemEarth, uvar->ephemSun )) != NULL, XLAL_EFUNC ); // ----- figure out which segments to use BOOLEAN manualSegments = XLALUserVarWasSet(&uvar->duration) || XLALUserVarWasSet(&uvar->startTime) || XLALUserVarWasSet(&uvar->Nseg); if ( manualSegments && uvar->segmentList ) { XLAL_ERROR ( XLAL_EDOM, "Can specify EITHER {--startTime, --duration, --Nseg} OR --segmentList\n"); } LIGOTimeGPS startTimeGPS; REAL8 duration; if ( uvar->segmentList == NULL ) { XLAL_CHECK ( uvar->Nseg >= 1, XLAL_EDOM, "Invalid input --Nseg=%d: number of segments must be >= 1\n", uvar->Nseg ); XLAL_CHECK ( uvar->duration >= 1, XLAL_EDOM, "Invalid input --duration=%f: duration must be >= 1 s\n", uvar->duration ); startTimeGPS = uvar->startTime; int ret = XLALSegListInitSimpleSegments ( &cfg->segmentList, startTimeGPS, uvar->Nseg, uvar->duration / uvar->Nseg ); XLAL_CHECK ( ret == XLAL_SUCCESS, XLAL_EFUNC, "XLALSegListInitSimpleSegments() failed with xlalErrno = %d\n", xlalErrno ); duration = uvar->duration; } else { LALSegList *segList = XLALReadSegmentsFromFile ( uvar->segmentList ); XLAL_CHECK ( segList != NULL, XLAL_EIO, "XLALReadSegmentsFromFile() failed to load segment list from file '%s', xlalErrno = %d\n", uvar->segmentList, xlalErrno ); cfg->segmentList = (*segList); // copy *contents* XLALFree ( segList ); startTimeGPS = cfg->segmentList.segs[0].start; UINT4 Nseg = cfg->segmentList.length; LIGOTimeGPS endTimeGPS = cfg->segmentList.segs[Nseg-1].end; duration = XLALGPSDiff( &endTimeGPS, &startTimeGPS ); } /* ----- figure out reference time */ LIGOTimeGPS refTimeGPS; /* treat special values first */ if ( uvar->refTime.gpsSeconds == 0 ) /* 0 = use startTime */ { refTimeGPS = uvar->startTime; } else if ( !XLALUserVarWasSet ( &uvar->refTime ) ) /* default = use mid-time of observation */ { refTimeGPS = startTimeGPS; XLALGPSAdd( &refTimeGPS, duration / 2.0 ); } else { refTimeGPS = uvar->refTime; } /* ----- get parameter-space point from user-input) */ cfg->signalParams.Amp.h0 = uvar->h0; cfg->signalParams.Amp.cosi = uvar->cosi; cfg->signalParams.Amp.psi = uvar->psi; cfg->signalParams.Amp.phi0 = uvar->phi0; { PulsarDopplerParams *dop = &(cfg->signalParams.Doppler); XLAL_INIT_MEM((*dop)); dop->refTime = refTimeGPS; dop->Alpha = uvar->Alpha; dop->Delta = uvar->Delta; dop->fkdot[0] = uvar->Freq; dop->fkdot[1] = uvar->f1dot; dop->fkdot[2] = uvar->f2dot; dop->fkdot[3] = uvar->f3dot; dop->asini = uvar->orbitasini; dop->period = uvar->orbitPeriod; dop->tp = uvar->orbitTp; dop->ecc = uvar->orbitEcc; dop->argp = uvar->orbitArgp; } /* ----- initialize IFOs and (Multi-)DetectorStateSeries ----- */ XLAL_CHECK ( XLALParseMultiLALDetector ( &cfg->multiIFO, uvar->IFOs ) == XLAL_SUCCESS, XLAL_EFUNC ); UINT4 numDet = cfg->multiIFO.length; XLAL_CHECK ( numDet >= 1, XLAL_EINVAL ); if ( uvar->sqrtSX ) { XLAL_CHECK ( XLALParseMultiNoiseFloor ( &cfg->multiNoiseFloor, uvar->sqrtSX, numDet ) == XLAL_SUCCESS, XLAL_EFUNC ); } /* ---------- translate coordinate system into internal representation ---------- */ if ( XLALDopplerCoordinateNames2System ( &cfg->coordSys, uvar->coords ) ) { LogPrintf (LOG_CRITICAL, "%s: Call to XLALDopplerCoordinateNames2System() failed. errno = %d\n\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } /* ---------- record full 'history' up to and including this application ---------- */ { CHAR *cmdline = NULL; CHAR *tmp; size_t len = strlen ( app_name ) + 1; if ( (cfg->history = XLALCalloc ( 1, sizeof(*cfg->history))) == NULL ) { LogPrintf (LOG_CRITICAL, "%s: XLALCalloc(1,%zu) failed.\n\n", __func__, sizeof(*cfg->history)); XLAL_ERROR ( XLAL_ENOMEM ); } if ( (tmp = XLALMalloc ( len )) == NULL ) { LogPrintf (LOG_CRITICAL, "%s: XLALMalloc (%zu) failed.\n\n", __func__, len ); XLAL_ERROR ( XLAL_ENOMEM ); } strcpy ( tmp, app_name ); cfg->history->app_name = tmp; /* get commandline describing search*/ if ( (cmdline = XLALUserVarGetLog ( UVAR_LOGFMT_CMDLINE )) == NULL ) { LogPrintf (LOG_CRITICAL, "%s: XLALUserVarGetLog() failed with xlalErrno = %d.\n\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } cfg->history->cmdline = cmdline; } /* record history */ return XLAL_SUCCESS; } /* XLALInitCode() */
/** * Unit test for metric functions XLALComputeDopplerPhaseMetric() * and XLALComputeDopplerFstatMetric() * * Initially modelled afer testMetricCodes.py script: * Check metric codes 'getMetric' 'FstatMetric' and 'FstatMetric_v2' by * comparing them against each other. * Given that they represent 3 very different implementations of * metric calculations, this provides a very powerful consistency test * */ static int test_XLALComputeDopplerMetrics ( void ) { int ret; const REAL8 tolPh = 0.01; // 1% tolerance on phase metrics [taken from testMetricCodes.py] // ----- load ephemeris const char earthEphem[] = TEST_DATA_DIR "earth00-19-DE200.dat.gz"; const char sunEphem[] = TEST_DATA_DIR "sun00-19-DE200.dat.gz"; EphemerisData *edat = XLALInitBarycenter ( earthEphem, sunEphem ); XLAL_CHECK ( edat != NULL, XLAL_EFUNC, "XLALInitBarycenter('%s','%s') failed with xlalErrno = %d\n", earthEphem, sunEphem, xlalErrno ); // ----- set test-parameters ---------- const LIGOTimeGPS startTimeGPS = { 792576013, 0 }; const REAL8 Tseg = 60000; const REAL8 Alpha = 1.0; const REAL8 Delta = 0.5; const REAL8 Freq = 100; const REAL8 f1dot = 0;// -1e-8; LALStringVector *detNames = XLALCreateStringVector ( "H1", "L1", "V1", NULL ); LALStringVector *sqrtSX = XLALCreateStringVector ( "1.0", "0.5", "1.5", NULL ); MultiLALDetector multiIFO; XLAL_CHECK ( XLALParseMultiLALDetector ( &multiIFO, detNames ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyStringVector ( detNames ); MultiNoiseFloor multiNoiseFloor; XLAL_CHECK ( XLALParseMultiNoiseFloor ( &multiNoiseFloor, sqrtSX, multiIFO.length ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyStringVector ( sqrtSX ); // prepare metric parameters for modern XLALComputeDopplerFstatMetric() and mid-old XLALOldDopplerFstatMetric() const DopplerCoordinateSystem coordSys = { 4, { DOPPLERCOORD_FREQ, DOPPLERCOORD_ALPHA, DOPPLERCOORD_DELTA, DOPPLERCOORD_F1DOT } }; const PulsarAmplitudeParams Amp = { 0.03, -0.3, 0.5, 0.0 }; // h0, cosi, psi, phi0 const PulsarDopplerParams dop = { .refTime = startTimeGPS, .Alpha = Alpha, .Delta = Delta, .fkdot = { Freq, f1dot }, }; LALSegList XLAL_INIT_DECL(segList); ret = XLALSegListInitSimpleSegments ( &segList, startTimeGPS, 1, Tseg ); XLAL_CHECK ( ret == XLAL_SUCCESS, XLAL_EFUNC, "XLALSegListInitSimpleSegments() failed with xlalErrno = %d\n", xlalErrno ); const DopplerMetricParams master_pars2 = { .coordSys = coordSys, .detMotionType = DETMOTION_SPIN | DETMOTION_ORBIT, .segmentList = segList, .multiIFO = multiIFO, .multiNoiseFloor = multiNoiseFloor, .signalParams = { .Amp = Amp, .Doppler = dop }, .projectCoord = - 1, // -1==no projection .approxPhase = 0, }; // ========== BEGINNING OF TEST CALLS ========== XLALPrintWarning("\n---------- ROUND 1: ephemeris-based, single-IFO phase metrics ----------\n"); { OldDopplerMetric *metric1; DopplerPhaseMetric *metric2P; REAL8 diff_2_1; DopplerMetricParams pars2 = master_pars2; pars2.multiIFO.length = 1; // truncate to first detector pars2.multiNoiseFloor.length = 1; // truncate to first detector // 1) compute metric using old FstatMetric code, now wrapped into XLALOldDopplerFstatMetric() XLAL_CHECK ( (metric1 = XLALOldDopplerFstatMetric ( OLDMETRIC_TYPE_PHASE, &pars2, edat )) != NULL, XLAL_EFUNC ); // 2) compute metric using modern UniversalDopplerMetric module: (used in lalapps_FstatMetric_v2) XLAL_CHECK ( (metric2P = XLALComputeDopplerPhaseMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); // compare metrics against each other: XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( metric2P->g_ij, metric1->g_ij )) < tolPh, XLAL_ETOL, "Error(g2,g1)= %g exceeds tolerance of %g\n", diff_2_1, tolPh ); XLALPrintWarning ("diff_2_1 = %e\n", diff_2_1 ); XLALDestroyOldDopplerMetric ( metric1 ); XLALDestroyDopplerPhaseMetric ( metric2P ); } XLALPrintWarning("\n---------- ROUND 2: Ptolemaic-based, single-IFO phase metrics ----------\n"); { OldDopplerMetric *metric1; DopplerPhaseMetric *metric2P; REAL8 diff_2_1; DopplerMetricParams pars2 = master_pars2; pars2.multiIFO.length = 1; // truncate to first detector pars2.multiNoiseFloor.length = 1; // truncate to first detector pars2.detMotionType = DETMOTION_SPIN | DETMOTION_PTOLEORBIT; // 1) compute metric using old FstatMetric code, now wrapped into XLALOldDopplerFstatMetric() XLAL_CHECK ( (metric1 = XLALOldDopplerFstatMetric ( OLDMETRIC_TYPE_PHASE, &pars2, edat )) != NULL, XLAL_EFUNC ); // 2) compute metric using modern UniversalDopplerMetric module: (used in lalapps_FstatMetric_v2) XLAL_CHECK ( (metric2P = XLALComputeDopplerPhaseMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); // compare all 3 metrics against each other: XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( metric2P->g_ij, metric1->g_ij )) < tolPh, XLAL_ETOL, "Error(g2,g1)= %g exceeds tolerance of %g\n", diff_2_1, tolPh ); XLALPrintWarning ("diff_2_1 = %e\n", diff_2_1 ); XLALDestroyOldDopplerMetric ( metric1 ); XLALDestroyDopplerPhaseMetric ( metric2P ); } XLALPrintWarning("\n---------- ROUND 3: ephemeris-based, multi-IFO F-stat metrics ----------\n"); { OldDopplerMetric *metric1; DopplerFstatMetric *metric2F; REAL8 diff_2_1; DopplerMetricParams pars2 = master_pars2; pars2.detMotionType = DETMOTION_SPIN | DETMOTION_ORBIT; pars2.multiIFO = multiIFO; // 3 IFOs pars2.multiNoiseFloor = multiNoiseFloor;// 3 IFOs // 1) compute metric using old FstatMetric code, now wrapped into XLALOldDopplerFstatMetric() XLAL_CHECK ( (metric1 = XLALOldDopplerFstatMetric ( OLDMETRIC_TYPE_FSTAT, &pars2, edat )) != NULL, XLAL_EFUNC ); // 2) compute metric using modern UniversalDopplerMetric module: (used in lalapps_FstatMetric_v2) XLAL_CHECK ( (metric2F = XLALComputeDopplerFstatMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); // compare both metrics against each other: XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( metric2F->gF_ij, metric1->gF_ij )) < tolPh, XLAL_ETOL, "Error(gF2,gF1)= %e exceeds tolerance of %e\n", diff_2_1, tolPh ); XLALPrintWarning ("gF: diff_2_1 = %e\n", diff_2_1 ); XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( metric2F->gFav_ij, metric1->gFav_ij )) < tolPh, XLAL_ETOL, "Error(gFav2,gFav1)= %e exceeds tolerance of %e\n", diff_2_1, tolPh ); XLALPrintWarning ("gFav: diff_2_1 = %e\n", diff_2_1 ); XLALDestroyOldDopplerMetric ( metric1 ); XLALDestroyDopplerFstatMetric ( metric2F ); } XLALPrintWarning("\n---------- ROUND 4: compare analytic {f,f1dot,f2dot,f3dot} phase metric vs XLALComputeDopplerPhaseMetric() ----------\n"); { DopplerPhaseMetric *metric2P; REAL8 diff_2_1; DopplerMetricParams pars2 = master_pars2; pars2.multiIFO.length = 1; // truncate to 1st detector pars2.multiNoiseFloor.length = 1; // truncate to 1st detector pars2.detMotionType = DETMOTION_SPIN | DETMOTION_ORBIT; pars2.approxPhase = 1; // use same phase-approximation as in analytic solution to improve comparison DopplerCoordinateSystem coordSys2 = { 4, { DOPPLERCOORD_FREQ, DOPPLERCOORD_F1DOT, DOPPLERCOORD_F2DOT, DOPPLERCOORD_F3DOT } }; pars2.coordSys = coordSys2; gsl_matrix* gN_ij; // a) compute metric at refTime = startTime pars2.signalParams.Doppler.refTime = startTimeGPS; XLAL_CHECK ( (metric2P = XLALComputeDopplerPhaseMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); gN_ij = NULL; XLAL_CHECK ( XLALNaturalizeMetric ( &gN_ij, NULL, metric2P->g_ij, &pars2 ) == XLAL_SUCCESS, XLAL_EFUNC ); REAL8 gStart_ij[] = { 1.0/3, 2.0/3, 6.0/5, 32.0/15, \ 2.0/3, 64.0/45, 8.0/3, 512.0/105, \ 6.0/5, 8.0/3, 36.0/7, 48.0/5, \ 32.0/15, 512.0/105, 48.0/5, 4096.0/225 }; const gsl_matrix_view gStart = gsl_matrix_view_array ( gStart_ij, 4, 4 ); // compare natural-units metric against analytic solution XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( gN_ij, &(gStart.matrix) )) < tolPh, XLAL_ETOL, "RefTime=StartTime: Error(g_ij,g_analytic)= %e exceeds tolerance of %e\n", diff_2_1, tolPh ); XLALPrintWarning ("Analytic (refTime=startTime): diff_2_1 = %e\n", diff_2_1 ); XLALDestroyDopplerPhaseMetric ( metric2P ); gsl_matrix_free ( gN_ij ); // b) compute metric at refTime = midTime pars2.signalParams.Doppler.refTime = startTimeGPS; pars2.signalParams.Doppler.refTime.gpsSeconds += Tseg / 2; XLAL_CHECK ( (metric2P = XLALComputeDopplerPhaseMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); gN_ij = NULL; XLAL_CHECK ( XLALNaturalizeMetric ( &gN_ij, NULL, metric2P->g_ij, &pars2 ) == XLAL_SUCCESS, XLAL_EFUNC ); REAL8 gMid_ij[] = { 1.0/3, 0, 1.0/5, 0, \ 0, 4.0/45, 0, 8.0/105, \ 1.0/5, 0, 1.0/7, 0, \ 0, 8.0/105, 0, 16.0/225 }; const gsl_matrix_view gMid = gsl_matrix_view_array ( gMid_ij, 4, 4 ); // compare natural-units metric against analytic solution XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( gN_ij, &(gMid.matrix) )) < tolPh, XLAL_ETOL, "RefTime=MidTime: Error(g_ij,g_analytic)= %e exceeds tolerance of %e\n", diff_2_1, tolPh ); XLALPrintWarning ("Analytic (refTime=midTime): diff_2_1 = %e\n\n", diff_2_1 ); XLALDestroyDopplerPhaseMetric ( metric2P ); gsl_matrix_free ( gN_ij ); } XLALPrintWarning("\n---------- ROUND 5: ephemeris-based, single-IFO, segment-averaged phase metrics ----------\n"); { OldDopplerMetric *metric1; DopplerPhaseMetric *metric2P; REAL8 diff_2_1; DopplerMetricParams pars2 = master_pars2; pars2.detMotionType = DETMOTION_SPIN | DETMOTION_ORBIT; pars2.multiIFO.length = 1; // truncate to first detector pars2.multiNoiseFloor.length = 1; // truncate to first detector pars2.approxPhase = 1; const UINT4 Nseg = 10; LALSegList XLAL_INIT_DECL(NsegList); ret = XLALSegListInitSimpleSegments ( &NsegList, startTimeGPS, Nseg, Tseg ); XLAL_CHECK ( ret == XLAL_SUCCESS, XLAL_EFUNC, "XLALSegListInitSimpleSegments() failed with xlalErrno = %d\n", xlalErrno ); pars2.segmentList = NsegList; LALSegList XLAL_INIT_DECL(segList_k); LALSeg segment_k; XLALSegListInit( &segList_k ); // prepare single-segment list containing segment k segList_k.arraySize = 1; segList_k.length = 1; segList_k.segs = &segment_k; // 1) compute metric using old FstatMetric code, now wrapped into XLALOldDopplerFstatMetric() metric1 = NULL; for (UINT4 k = 0; k < Nseg; ++k) { // setup 1-segment segment-list pointing k-th segment DopplerMetricParams pars2_k = pars2; pars2_k.segmentList = segList_k; pars2_k.segmentList.segs[0] = pars2.segmentList.segs[k]; // XLALOldDopplerFstatMetric() does not agree numerically with UniversalDopplerMetric when using refTime != startTime pars2_k.signalParams.Doppler.refTime = pars2_k.segmentList.segs[0].start; OldDopplerMetric *metric1_k; // per-segment coherent metric XLAL_CHECK ( (metric1_k = XLALOldDopplerFstatMetric ( OLDMETRIC_TYPE_PHASE, &pars2_k, edat )) != NULL, XLAL_EFUNC ); // manually correct reference time of metric1_k->g_ij; see Prix, "Frequency metric for CW searches" (2014-08-17), p. 4 const double dt = XLALGPSDiff( &(pars2_k.signalParams.Doppler.refTime), &(pars2.signalParams.Doppler.refTime) ); const double gFF = gsl_matrix_get( metric1_k->g_ij, 0, 0 ); const double gFA = gsl_matrix_get( metric1_k->g_ij, 0, 1 ); const double gFD = gsl_matrix_get( metric1_k->g_ij, 0, 2 ); const double gFf = gsl_matrix_get( metric1_k->g_ij, 0, 3 ); const double gAf = gsl_matrix_get( metric1_k->g_ij, 1, 3 ); const double gDf = gsl_matrix_get( metric1_k->g_ij, 2, 3 ); const double gff = gsl_matrix_get( metric1_k->g_ij, 3, 3 ); gsl_matrix_set( metric1_k->g_ij, 0, 3, gFf + gFF*dt ); gsl_matrix_set( metric1_k->g_ij, 3, 0, gsl_matrix_get( metric1_k->g_ij, 0, 3 ) ); gsl_matrix_set( metric1_k->g_ij, 1, 3, gAf + gFA*dt ); gsl_matrix_set( metric1_k->g_ij, 3, 1, gsl_matrix_get( metric1_k->g_ij, 1, 3 ) ); gsl_matrix_set( metric1_k->g_ij, 2, 3, gDf + gFD*dt ); gsl_matrix_set( metric1_k->g_ij, 3, 2, gsl_matrix_get( metric1_k->g_ij, 2, 3 ) ); gsl_matrix_set( metric1_k->g_ij, 3, 3, gff + 2*gFf*dt + gFF*dt*dt ); XLAL_CHECK ( XLALAddOldDopplerMetric ( &metric1, metric1_k ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyOldDopplerMetric ( metric1_k ); } XLAL_CHECK ( XLALScaleOldDopplerMetric ( metric1, 1.0 / Nseg ) == XLAL_SUCCESS, XLAL_EFUNC ); // 2) compute metric using modern UniversalDopplerMetric module: (used in lalapps_FstatMetric_v2) XLAL_CHECK ( (metric2P = XLALComputeDopplerPhaseMetric ( &pars2, edat )) != NULL, XLAL_EFUNC ); GPMAT( metric1->g_ij, "%0.4e" ); GPMAT( metric2P->g_ij, "%0.4e" ); // compare both metrics against each other: XLAL_CHECK ( (diff_2_1 = XLALCompareMetrics ( metric2P->g_ij, metric1->g_ij )) < tolPh, XLAL_ETOL, "Error(g2,g1)= %g exceeds tolerance of %g\n", diff_2_1, tolPh ); XLALPrintWarning ("diff_2_1 = %e\n", diff_2_1 ); XLALDestroyOldDopplerMetric ( metric1 ); XLALDestroyDopplerPhaseMetric ( metric2P ); XLALSegListClear ( &NsegList ); } XLALPrintWarning("\n---------- ROUND 6: directed binary orbital metric ----------\n"); { REAL8 Period = 68023.70496; REAL8 Omega = LAL_TWOPI / Period; REAL8 asini = 1.44; REAL8 tAsc = 897753994; REAL8 argp = 0; LIGOTimeGPS tP; XLALGPSSetREAL8 ( &tP, tAsc + argp / Omega ); const PulsarDopplerParams dopScoX1 = { .refTime = startTimeGPS, .Alpha = Alpha, .Delta = Delta, .fkdot = { Freq }, .asini = asini, .period = Period, .tp = tP }; REAL8 TspanScoX1 = 20 * 19 * 3600; // 20xPorb for long-segment regime LALSegList XLAL_INIT_DECL(segListScoX1); XLAL_CHECK ( XLALSegListInitSimpleSegments ( &segListScoX1, startTimeGPS, 1, TspanScoX1 ) == XLAL_SUCCESS, XLAL_EFUNC ); REAL8 tMid = XLALGPSGetREAL8(&startTimeGPS) + 0.5 * TspanScoX1; REAL8 DeltaMidAsc = tMid - tAsc; const DopplerCoordinateSystem coordSysScoX1 = { 6, { DOPPLERCOORD_FREQ, DOPPLERCOORD_ASINI, DOPPLERCOORD_TASC, DOPPLERCOORD_PORB, DOPPLERCOORD_KAPPA, DOPPLERCOORD_ETA } }; DopplerMetricParams pars_ScoX1 = { .coordSys = coordSysScoX1, .detMotionType = DETMOTION_SPIN | DETMOTION_ORBIT, .segmentList = segListScoX1, .multiIFO = multiIFO, .multiNoiseFloor = multiNoiseFloor, .signalParams = { .Amp = Amp, .Doppler = dopScoX1 }, .projectCoord = - 1, // -1==no projection .approxPhase = 1, }; pars_ScoX1.multiIFO.length = 1; // truncate to first detector pars_ScoX1.multiNoiseFloor.length = 1; // truncate to first detector // compute metric using modern UniversalDopplerMetric module: (used in lalapps_FstatMetric_v2) DopplerPhaseMetric *metric_ScoX1; XLAL_CHECK ( (metric_ScoX1 = XLALComputeDopplerPhaseMetric ( &pars_ScoX1, edat )) != NULL, XLAL_EFUNC ); // compute analytic metric computed from Eq.(47) in Leaci,Prix PRD91, 102003 (2015): gsl_matrix *g0_ij; XLAL_CHECK ( (g0_ij = gsl_matrix_calloc ( 6, 6 )) != NULL, XLAL_ENOMEM, "Failed to gsl_calloc a 6x6 matrix\n"); gsl_matrix_set ( g0_ij, 0, 0, pow ( LAL_PI * TspanScoX1, 2 ) / 3.0 ); gsl_matrix_set ( g0_ij, 1, 1, 2.0 * pow ( LAL_PI * Freq, 2 ) ); gsl_matrix_set ( g0_ij, 2, 2, 2.0 * pow ( LAL_PI * Freq * asini * Omega, 2 ) ); gsl_matrix_set ( g0_ij, 3, 3, 0.5 * pow ( Omega, 4 ) * pow ( Freq * asini, 2 ) * ( pow ( TspanScoX1, 2 ) / 12.0 + pow ( DeltaMidAsc, 2 ) ) ); REAL8 gPAsc = LAL_PI * pow ( Freq * asini, 2 ) * pow ( Omega, 3 ) * DeltaMidAsc; gsl_matrix_set ( g0_ij, 2, 3, gPAsc ); gsl_matrix_set ( g0_ij, 3, 2, gPAsc ); gsl_matrix_set ( g0_ij, 4, 4, 0.5 * pow ( LAL_PI * Freq * asini, 2 ) ); gsl_matrix_set ( g0_ij, 5, 5, 0.5 * pow ( LAL_PI * Freq * asini, 2 ) ); GPMAT ( metric_ScoX1->g_ij, "%0.4e" ); GPMAT ( g0_ij, "%0.4e" ); // compare metrics against each other REAL8 diff, tolScoX1 = 0.05; XLAL_CHECK ( (diff = XLALCompareMetrics ( metric_ScoX1->g_ij, g0_ij )) < tolScoX1, XLAL_ETOL, "Error(gNum,gAn)= %g exceeds tolerance of %g\n", diff, tolScoX1 ); XLALPrintWarning ("diff_Num_An = %e\n", diff ); gsl_matrix_free ( g0_ij ); XLALDestroyDopplerPhaseMetric ( metric_ScoX1 ); XLALSegListClear ( &segListScoX1 ); }