/** Destroy a DetectorStateSeries (and set it to NULL) */ void LALDestroyDetectorStateSeries (LALStatus *status, /**< pointer to LALStatus structure */ DetectorStateSeries **detStates ) /**< pointer to vector to be destroyed */ { INITSTATUS(status); ASSERT ( detStates, status, DETECTORSTATES_ENULL, DETECTORSTATES_MSGENULL); XLALDestroyDetectorStateSeries ( (*detStates) ); (*detStates) = NULL; RETURN (status); } /* LALDestroyDetectorStateSeries() */
/** * Helper function to get rid of a multi-IFO DetectorStateSeries * Note, this is "NULL-robust" in the sense that it will not crash * on NULL-entries anywhere in this struct, so it can be used * for failure-cleanup even on incomplete structs. */ void XLALDestroyMultiDetectorStateSeries ( MultiDetectorStateSeries *mdetStates ) { UINT4 X, numDet; if ( !mdetStates ) return; numDet = mdetStates->length; if ( mdetStates->data ) { for ( X=0; X < numDet ; X ++ ) XLALDestroyDetectorStateSeries ( mdetStates->data[X] ); LALFree ( mdetStates->data ); } LALFree ( mdetStates ); return; } /* XLALDestroyMultiDetectorStateSeries() */
/** * Get the 'detector state' (ie detector-tensor, position, velocity, etc) for the given * vector of timestamps, shifted by a common time-shift \a tOffset. * * This function just calls XLALBarycenterEarth() and XLALBarycenter() for the * given vector of timestamps (shifted by tOffset) and returns the positions, * velocities and LMSTs of the detector, stored in a DetectorStateSeries. * There is also an entry containing the EarthState at each timestamp, which * can be used as input for subsequent calls to XLALBarycenter(). * * \a tOffset allows one to easily use the midpoints of SFT-timestamps, for example. * */ DetectorStateSeries * XLALGetDetectorStates ( const LIGOTimeGPSVector *timestamps, /**< array of GPS timestamps t_i */ const LALDetector *detector, /**< detector info */ const EphemerisData *edat, /**< ephemeris file data */ REAL8 tOffset /**< compute detector states at timestamps SHIFTED by tOffset */ ) { /* check input consistency */ if ( !timestamps || !detector || !edat ) { XLALPrintError ("%s: invalid NULL input, timestamps=%p, detector=%p, edat=%p\n", __func__, timestamps, detector, edat ); XLAL_ERROR_NULL ( XLAL_EINVAL ); } /* prepare return vector */ UINT4 numSteps = timestamps->length; DetectorStateSeries *ret = NULL; if ( ( ret = XLALCreateDetectorStateSeries ( numSteps )) == NULL ) { XLALPrintError ("%s: XLALCreateDetectorStateSeries(%d) failed.\n", __func__, numSteps ); XLAL_ERROR_NULL ( XLAL_EFUNC ); } /* enter detector-info into the head of the state-vector */ ret->detector = (*detector); /* set 'time-span' associated with each timestamp */ ret->deltaT = timestamps->deltaT; /* set SSB coordinate system used: EQUATORIAL for Earth-based, ECLIPTIC for LISA */ if ( detector->frDetector.prefix[0] == 'Z' ) /* LISA */ ret->system = COORDINATESYSTEM_ECLIPTIC; else /* Earth-based */ ret->system = COORDINATESYSTEM_EQUATORIAL; /* now fill all the vector-entries corresponding to different timestamps */ UINT4 i; for ( i=0; i < numSteps; i++ ) { BarycenterInput baryinput; EmissionTime emit; DetectorState *state = &(ret->data[i]); EarthState *earth = &(state->earthState); LIGOTimeGPS tgps; /* shift timestamp by tOffset */ tgps = timestamps->data[i]; XLALGPSAdd(&tgps, tOffset); /*----- first get earth-state */ if ( XLALBarycenterEarth ( earth, &tgps, edat ) != XLAL_SUCCESS ) { XLALDestroyDetectorStateSeries ( ret ); XLALPrintError("%s: XLALBarycenterEarth() failed with xlalErrno=%d\n", __func__, xlalErrno ); XLAL_ERROR_NULL ( XLAL_EFAILED ); } /*----- then get detector-specific info */ baryinput.tgps = tgps; baryinput.site = (*detector); baryinput.site.location[0] /= LAL_C_SI; baryinput.site.location[1] /= LAL_C_SI; baryinput.site.location[2] /= LAL_C_SI; baryinput.alpha = baryinput.delta = 0; /* irrelevant */ baryinput.dInv = 0; if ( XLALBarycenter ( &emit, &baryinput, earth) != XLAL_SUCCESS ) { XLALDestroyDetectorStateSeries( ret ); XLALPrintError("%s: XLALBarycenterEarth() failed with xlalErrno=%d\n", __func__, xlalErrno ); XLAL_ERROR_NULL ( XLAL_EFAILED ); } /*----- extract the output-data from this */ UINT4 j; for (j=0; j < 3; j++) /* copy detector's position and velocity */ { state->rDetector[j] = emit.rDetector[j]; state->vDetector[j] = emit.vDetector[j]; } /* for j < 3 */ /* local mean sidereal time = GMST + longitude */ state->LMST = earth->gmstRad + detector->frDetector.vertexLongitudeRadians; state->LMST = fmod (state->LMST, LAL_TWOPI ); /* normalize */ /* insert timestamp */ state->tGPS = tgps; /* compute the detector-tensor at this time-stamp in SSB-fixed Cartesian coordinates * [EQUATORIAL for Earth-based, ECLIPTIC for LISA] */ if ( XLALFillDetectorTensor ( state, detector ) != 0 ) { XLALDestroyDetectorStateSeries(ret); XLALPrintError ( "%s: XLALFillDetectorTensor() failed ... errno = %d\n\n", __func__, xlalErrno ); XLAL_ERROR_NULL ( XLAL_EFUNC ); } } /* for i < numSteps */ /* return result */ return ret; } /* XLALGetDetectorStates() */
/** * Very simple test: pick random skyposition, compute a_i, b_i using * once LALComputeAM() and once LALNewGetAMCoeffs(), and look at the errors * sum_i (a_i - a_i')^2 */ int main(int argc, char *argv[]) { LALStatus XLAL_INIT_DECL(status); int opt; /* Command-line option. */ LIGOTimeGPS startTime = {714180733, 0}; REAL8 duration = 180000; /* 50 hours */ REAL8 Tsft = 1800; /* assume 30min SFTs */ LIGOTimeGPSVector *timestamps = NULL; DetectorStateSeries *detStates = NULL; SkyPosition XLAL_INIT_DECL(skypos); EphemerisData XLAL_INIT_DECL(edat); BarycenterInput XLAL_INIT_DECL(baryinput); LALDetector *det = NULL; AMCoeffs XLAL_INIT_DECL(AMold); AMCoeffs XLAL_INIT_DECL(AMnew1); AMCoeffs XLAL_INIT_DECL(AMnew2); REAL8 alpha, delta; AMCoeffsParams XLAL_INIT_DECL(amParams); EarthState earth; UINT4 i; REAL8 maxerr01, maxerr02, maxerr12, averr01, averr02, averr12; REAL8 tolerance = 1e-2; /* be generous: allow 1% error */ struct tms buf; const CHAR *sites[] = {"H1", "L1", "V2", "G1", "T1" }; REAL8 sinzeta; /* zeta = IFO opening angle */ UINT4 pickedSite; BOOLEAN ignoreErrors = 0; /* Don't fail if tolerance exceeded */ UINT4 numChecks = 1; /* Number of times to check */ char earthEphem[] = TEST_DATA_DIR "earth00-19-DE405.dat.gz"; char sunEphem[] = TEST_DATA_DIR "sun00-19-DE405.dat.gz"; /* ----- old testing code to use 9 degree earth rotations ----- */ /* startTime.gpsSeconds = 714275242; duration = 86164; Tsft = 2154.1; */ while ((opt = LALgetopt( argc, argv, "n:qv:" )) != -1) { switch (opt) { case 'v': /* set lalDebugLevel */ break; case 'q': /* don't fail if tolerance exceeded */ ignoreErrors = 1; break; case 'n': /* number of times to check */ numChecks = atoi( LALoptarg ); break; } } /* init random-generator */ srand ( times(&buf) ); /* ----- init ephemeris ----- */ edat.ephiles.earthEphemeris = earthEphem; edat.ephiles.sunEphemeris = sunEphem; SUB ( LALInitBarycenter(&status, &edat), &status); /* ----- get timestamps ----- */ SUB ( LALMakeTimestamps ( &status, ×tamps, startTime, duration, Tsft ), &status ); /* ----- allocate memory for AM-coeffs ----- */ AMold.a = XLALCreateREAL4Vector ( timestamps->length ); AMold.b = XLALCreateREAL4Vector ( timestamps->length ); AMnew1.a = XLALCreateREAL4Vector ( timestamps->length ); AMnew1.b = XLALCreateREAL4Vector ( timestamps->length ); AMnew2.a = XLALCreateREAL4Vector ( timestamps->length ); AMnew2.b = XLALCreateREAL4Vector ( timestamps->length ); while ( numChecks-- ) { /* ----- pick detector-site at random ----- */ pickedSite = floor( 5 * (1.0 * rand() / (RAND_MAX + 1.0) ) ); /* int in [0,5) */ /* NOTE: contrary to ComputeAM() and LALGetAMCoffs(), the new function LALNewGetAMCoeffs() * computes 'a * sinzeta' and 'b * sinzeta': for the comparison we therefore need to correct * for GEO's opening-angle of 94.33degrees [JKS98]: */ if ( ! strcmp ( sites[pickedSite], "G1" ) ) sinzeta = 0.997146; else sinzeta = 1; if ( ( det = XLALGetSiteInfo ( sites[pickedSite] )) == NULL ) { XLALPrintError ("\nCall to XLALGetSiteInfo() has failed for site = '%s'... \n\n", sites[pickedSite]); return NEWGETAMCOEFFSTEST_ESUB; } /* ----- pick skyposition at random ----- */ alpha = LAL_TWOPI * (1.0 * rand() / ( RAND_MAX + 1.0 ) ); /* uniform in [0, 2pi) */ delta = LAL_PI_2 - acos ( 1 - 2.0 * rand()/RAND_MAX ); /* sin(delta) uniform in [-1,1] */ /* ----- old testing code to put source overhead ----- */ /* alpha = det->frDetector.vertexLongitudeRadians; delta = det->frDetector.vertexLatitudeRadians; */ /* ===== compute AM-coeffs the 'old way': ===== */ baryinput.site.location[0] = det->location[0]/LAL_C_SI; baryinput.site.location[1] = det->location[1]/LAL_C_SI; baryinput.site.location[2] = det->location[2]/LAL_C_SI; baryinput.alpha = alpha; baryinput.delta = delta; baryinput.dInv = 0.e0; /* amParams structure to compute a(t) and b(t) */ amParams.das = (LALDetAndSource *)LALMalloc(sizeof(LALDetAndSource)); amParams.das->pSource = (LALSource *)LALMalloc(sizeof(LALSource)); amParams.baryinput = &baryinput; amParams.earth = &earth; amParams.edat = &edat; amParams.das->pDetector = det; amParams.das->pSource->equatorialCoords.longitude = alpha; amParams.das->pSource->equatorialCoords.latitude = delta; amParams.das->pSource->orientation = 0.0; amParams.das->pSource->equatorialCoords.system = COORDINATESYSTEM_EQUATORIAL; amParams.polAngle = 0; SUB (LALComputeAM ( &status, &AMold, timestamps->data, &amParams), &status); /* ===== compute AM-coeffs the 'new way' using LALNewGetAMCoeffs() */ /* ----- get detector-state series ----- */ SUB ( LALGetDetectorStates (&status, &detStates, timestamps, det, &edat, 0 ), &status ); skypos.system = COORDINATESYSTEM_EQUATORIAL; skypos.longitude = alpha; skypos.latitude = delta; /* the 'new' and the 'newer' way ... */ SUB ( LALGetAMCoeffs ( &status, &AMnew1, detStates, skypos ), &status ); /* 'new1' */ SUB ( LALNewGetAMCoeffs ( &status, &AMnew2, detStates, skypos ), &status ); /* 'new2' */ /* ===== analyse relative errors ===== */ maxerr01 = maxerr02 = maxerr12 = 0; /* errors between 0='old', 1='new1', 2='new2' */ averr01 = averr02 = averr12 = 0; for ( i=0; i < timestamps->length; i ++ ) { /* printf("GPS time: %d s %d ns; GMST in radians: %f\n", detStates->data[i].tGPS.gpsSeconds, detStates->data[i].tGPS.gpsNanoSeconds, fmod(detStates->data[i].earthState.gmstRad,LAL_TWOPI)); printf("Old AM coeffs: a=%f, b=%f\nNew AM coeffs: a=%f, b=%f\nNEWER AM coeffs: a=%f b=%f", AMold.a->data[i], AMold.b->data[i], AMnew.a->data[i], AMnew.b->data[i], AMnewer.a->data[i], AMnewer.b->data[i]); */ REAL8 thisErr; /* compare 0-1 */ thisErr = sqrt( SQ ( AMold.a->data[i] - AMnew1.a->data[i] ) / AMold.A ); averr01 += thisErr; maxerr01 = MYMAX( thisErr, maxerr01 ); thisErr = sqrt( SQ ( AMold.b->data[i] - AMnew1.b->data[i] ) / AMold.B ); averr01 += thisErr; maxerr01 = MYMAX( thisErr, maxerr01 ); /* compare 0-2 */ thisErr = sqrt( SQ ( AMold.a->data[i] - AMnew2.a->data[i]/sinzeta ) / AMold.A ); averr02 += thisErr; maxerr02 = MYMAX( thisErr, maxerr02 ); thisErr = sqrt( SQ ( AMold.b->data[i] - AMnew2.b->data[i]/sinzeta ) / AMold.B ); averr02 += thisErr; maxerr02 = MYMAX( thisErr, maxerr02 ); /* compare 1-2 */ thisErr = sqrt( SQ ( AMnew1.a->data[i] - AMnew2.a->data[i]/sinzeta ) / AMold.A ); averr12 += thisErr; maxerr12 = MYMAX( thisErr, maxerr12 ); thisErr = sqrt( SQ ( AMnew1.b->data[i] - AMnew2.b->data[i]/sinzeta ) / AMold.B ); averr12 += thisErr; maxerr12 = MYMAX( thisErr, maxerr12 ); } averr01 /= 2.0 * timestamps->length; averr02 /= 2.0 * timestamps->length; averr12 /= 2.0 * timestamps->length; if ( lalDebugLevel ) { printf ("Parameters: IFO = %s, skypos = [%g, %g]\n", sites[pickedSite], alpha, delta ); printf ("Maximal relative errors: maxerr(0-1) = %g %%, maxerr(0-2) = %g %% maxerr(1-2) = %g %%\n", 100.0 * maxerr01, 100.0 * maxerr02, 100.0 * maxerr12 ); printf ("Average relative errors: averr(0-1) = %g %%, averr(0-2) = %g %% averr(1-2) = %g %%\n", 100.0 * averr01, 100.0 * averr02, 100.0 * averr12 ); } else printf ("%d %g %g \t %g %g %g \t %g %g %g\n", pickedSite, alpha, delta, averr01, averr02, averr12, maxerr01, maxerr02, maxerr12); if ( (averr01 > tolerance) || (averr02 > tolerance) || (averr12 > tolerance) || (maxerr01 > tolerance) ||(maxerr02 > tolerance) || (maxerr12 > tolerance) ) { XLALPrintError ("Maximal error-tolerance of %g %% was exceeded!\n", 100.0 * tolerance ); if (!ignoreErrors) return 1; } if ( lalDebugLevel ) printf("%d checks left\n", numChecks); /* ---- Clean up things that were created in this loop ---- */ XLALDestroyDetectorStateSeries ( detStates ); detStates = NULL; LALFree ( det ); LALFree ( amParams.das->pSource ); LALFree ( amParams.das ); } /* ----- free memory ----- */ XLALDestroyTimestampVector ( timestamps ); XLALDestroyREAL4Vector ( AMold.a ); XLALDestroyREAL4Vector ( AMold.b ); XLALDestroyREAL4Vector ( AMnew1.a ); XLALDestroyREAL4Vector ( AMnew1.b ); XLALDestroyREAL4Vector ( AMnew2.a ); XLALDestroyREAL4Vector ( AMnew2.b ); LALFree(edat.ephemE); LALFree(edat.ephemS); LALCheckMemoryLeaks(); return 0; /* OK */ } /* main() */
int main(int argc, char *argv[]) { LALStatus status = blank_status; ConfigVariables XLAL_INIT_DECL(config); UserVariables_t XLAL_INIT_DECL(uvar); /* register user-variables */ XLAL_CHECK ( XLALInitUserVars ( &uvar ) == XLAL_SUCCESS, XLAL_EFUNC ); /* read cmdline & cfgfile */ BOOLEAN should_exit = 0; XLAL_CHECK( XLALUserVarReadAllInput( &should_exit, argc, argv ) == XLAL_SUCCESS, XLAL_EFUNC ); if ( should_exit ) { exit(1); } if ( uvar.version ) { XLALOutputVersionString ( stdout, lalDebugLevel ); exit(0); } /* basic setup and initializations */ XLAL_CHECK ( XLALInitCode( &config, &uvar, argv[0] ) == XLAL_SUCCESS, XLAL_EFUNC ); /* ----- allocate memory for AM-coeffs ----- */ AMCoeffs AMold, AMnew1, AMnew2; /**< containers holding AM-coefs computed by 3 different AM functions */ AMold.a = XLALCreateREAL4Vector ( 1 ); AMold.b = XLALCreateREAL4Vector ( 1 ); AMnew1.a = XLALCreateREAL4Vector ( 1 ); AMnew1.b = XLALCreateREAL4Vector ( 1 ); AMnew2.a = XLALCreateREAL4Vector ( 1 ); AMnew2.b = XLALCreateREAL4Vector ( 1 ); XLAL_CHECK ( AMold.a && AMold.b && AMnew1.a && AMnew1.b && AMnew2.a && AMnew2.a, XLAL_ENOMEM, "Failed to XLALCreateREAL4Vector ( 1 )\n" ); /* ----- get detector-state series ----- */ DetectorStateSeries *detStates = NULL; XLAL_CHECK ( (detStates = XLALGetDetectorStates ( config.timestamps, config.det, config.edat, 0 )) != NULL, XLAL_EFUNC ); /* ----- compute associated SSB timing info ----- */ SSBtimes *tSSB = XLALGetSSBtimes ( detStates, config.skypos, config.timeGPS, SSBPREC_RELATIVISTIC ); XLAL_CHECK ( tSSB != NULL, XLAL_EFUNC, "XLALGetSSBtimes() failed with xlalErrno = %d\n", xlalErrno ); /* ===== 1) compute AM-coeffs the 'old way': [used in CFSv1] ===== */ BarycenterInput XLAL_INIT_DECL(baryinput); AMCoeffsParams XLAL_INIT_DECL(amParams); EarthState earth; baryinput.site.location[0] = config.det->location[0]/LAL_C_SI; baryinput.site.location[1] = config.det->location[1]/LAL_C_SI; baryinput.site.location[2] = config.det->location[2]/LAL_C_SI; baryinput.alpha = config.skypos.longitude; baryinput.delta = config.skypos.latitude; baryinput.dInv = 0.e0; /* amParams structure to compute a(t) and b(t) */ amParams.das = XLALMalloc(sizeof(*amParams.das)); amParams.das->pSource = XLALMalloc(sizeof(*amParams.das->pSource)); amParams.baryinput = &baryinput; amParams.earth = &earth; amParams.edat = config.edat; amParams.das->pDetector = config.det; amParams.das->pSource->equatorialCoords.longitude = config.skypos.longitude; amParams.das->pSource->equatorialCoords.latitude = config.skypos.latitude; amParams.das->pSource->orientation = 0.0; amParams.das->pSource->equatorialCoords.system = COORDINATESYSTEM_EQUATORIAL; amParams.polAngle = 0; LAL_CALL ( LALComputeAM ( &status, &AMold, config.timestamps->data, &amParams), &status); XLALFree ( amParams.das->pSource ); XLALFree ( amParams.das ); /* ===== 2) compute AM-coeffs the 'new way' using LALNewGetAMCoeffs() */ LALGetAMCoeffs ( &status, &AMnew1, detStates, config.skypos ); if ( status.statusCode ) { XLALPrintError ("%s: call to LALGetAMCoeffs() failed, status = %d\n\n", __func__, status.statusCode ); XLAL_ERROR ( status.statusCode & XLAL_EFUNC ); } /* ===== 3) compute AM-coeffs the 'newer way' using LALNewGetAMCoeffs() [used in CFSv2] */ LALNewGetAMCoeffs ( &status, &AMnew2, detStates, config.skypos ); if ( status.statusCode ) { XLALPrintError ("%s: call to LALNewGetAMCoeffs() failed, status = %d\n\n", __func__, status.statusCode ); XLAL_ERROR ( status.statusCode & XLAL_EFUNC ); } /* ===== 4) use standalone version of the above [used in FstatMetric_v2] */ REAL8 a0,b0; if ( XLALComputeAntennaPatternCoeffs ( &a0, &b0, &config.skypos, &config.timeGPS, config.det, config.edat ) != XLAL_SUCCESS ) { XLALPrintError ("%s: XLALComputeAntennaPatternCoeffs() failed.\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } /* ==================== output the results ==================== */ printf ("\n"); printf ("----- Input parameters:\n"); printf ("tGPS = { %d, %d }\n", config.timeGPS.gpsSeconds, config.timeGPS.gpsNanoSeconds ); printf ("Detector = %s\n", config.det->frDetector.name ); printf ("Sky position: longitude = %g rad, latitude = %g rad [equatorial coordinates]\n", config.skypos.longitude, config.skypos.latitude ); printf ("\n"); printf ("----- Antenna pattern functions (a,b):\n"); printf ("LALComputeAM: ( %-12.8g, %-12.8g) [REAL4]\n", AMold.a->data[0], AMold.b->data[0] ); printf ("LALGetAMCoeffs: ( %-12.8g, %-12.8g) [REAL4]\n", AMnew1.a->data[0], AMnew1.b->data[0] ); printf ("LALNewGetAMCoeffs: ( %-12.8g, %-12.8g) [REAL4]\n", AMnew2.a->data[0]/config.sinzeta, AMnew2.b->data[0]/config.sinzeta ); printf ("XLALComputeAntennaPatternCoeffs: ( %-12.8g, %-12.8g) [REAL8]\n", a0/config.sinzeta, b0/config.sinzeta ); printf ("\n"); printf ("----- Detector & Earth state:\n"); REAL8 *pos = detStates->data[0].rDetector; printf ("Detector position [ICRS J2000. Units=sec]: rDet = {%g, %g, %g}\n", pos[0], pos[1], pos[2] ); REAL8 *vel = detStates->data[0].vDetector; printf ("Detector velocity [ICRS J2000. Units=c]: vDet = {%g, %g, %g}\n", vel[0], vel[1], vel[2] ); printf ("Local mean sideral time: LMST = %g rad\n", detStates->data[0].LMST); printf ("\n"); printf ("----- SSB timing data:\n"); printf ("TOA difference tSSB - tDet = %g s\n", tSSB->DeltaT->data[0] ); printf ("TOA rate of change dtSSB/dtDet - 1 = %g\n", tSSB->Tdot->data[0] - 1.0 ); printf ("\n\n"); /* ----- done: free all memory */ XLAL_CHECK ( XLALDestroyConfig( &config ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyDetectorStateSeries ( detStates ); XLALDestroyREAL4Vector ( AMold.a ); XLALDestroyREAL4Vector ( AMold.b ); XLALDestroyREAL4Vector ( AMnew1.a ); XLALDestroyREAL4Vector ( AMnew1.b ); XLALDestroyREAL4Vector ( AMnew2.a ); XLALDestroyREAL4Vector ( AMnew2.b ); XLALDestroyREAL8Vector ( tSSB->DeltaT ); XLALDestroyREAL8Vector ( tSSB->Tdot ); XLALFree (tSSB); LALCheckMemoryLeaks(); return 0; } /* main */
/** * Very simple test: pick random skyposition, compute a_i, b_i using * once LALComputeAM() and once LALGetAMCoeffs(), and look at the errors * sum_i (a_i - a_i')^2 */ int main(int argc, char *argv[]) { LALStatus XLAL_INIT_DECL(status); LIGOTimeGPS startTime = {714180733, 0}; REAL8 duration = 180000; /* 50 hours */ REAL8 Tsft = 1800; /* assume 30min SFTs */ LIGOTimeGPSVector *timestamps = NULL; DetectorStateSeries *detStates = NULL; SkyPosition XLAL_INIT_DECL(skypos); EphemerisData XLAL_INIT_DECL(edat); BarycenterInput XLAL_INIT_DECL(baryinput); LALDetector *det = NULL; AMCoeffs XLAL_INIT_DECL(AMold); AMCoeffs XLAL_INIT_DECL(AMnew); REAL8 alpha, delta; AMCoeffsParams XLAL_INIT_DECL(amParams); EarthState earth; UINT4 i; REAL8 maxerr_a, maxerr_b, averr_a, averr_b; REAL8 tolerance = 1e-2; /* be generous: allow 1% error */ struct tms buf; const CHAR *sites[] = {"H1", "L1", "V2", "G1", "T1" }; UINT4 pickedSite; char earthEphem[] = TEST_DATA_DIR "earth00-19-DE405.dat.gz"; char sunEphem[] = TEST_DATA_DIR "sun00-19-DE405.dat.gz"; if ( argc == 2 && !strcmp(argv[1], "-v1") ) /* init random-generator */ srand ( times(&buf) ); /* ----- init ephemeris ----- */ edat.ephiles.earthEphemeris = earthEphem; edat.ephiles.sunEphemeris = sunEphem; SUB ( LALInitBarycenter(&status, &edat), &status); /* ----- get timestamps ----- */ SUB ( LALMakeTimestamps ( &status, ×tamps, startTime, duration, Tsft ), &status ); /* ----- allocate memory for AM-coeffs ----- */ AMold.a = XLALCreateREAL4Vector ( timestamps->length ); AMold.b = XLALCreateREAL4Vector ( timestamps->length ); AMnew.a = XLALCreateREAL4Vector ( timestamps->length ); AMnew.b = XLALCreateREAL4Vector ( timestamps->length ); /* ----- pick detector-site at random ----- */ pickedSite = floor( 5 * (1.0 * rand() / (RAND_MAX + 1.0) ) ); /* int in [0,5) */ if ( ( det = XLALGetSiteInfo ( sites[pickedSite] )) == NULL ) { XLALPrintError ("\nCall to XLALGetSiteInfo() has failed for site = '%s'... \n\n", sites[pickedSite]); return GETAMCOEFFSTEST_ESUB; } /* ----- pick skyposition at random ----- */ alpha = LAL_TWOPI * (1.0 * rand() / ( RAND_MAX + 1.0 ) ); /* uniform in [0, 2pi) */ delta = LAL_PI_2 - acos ( 1 - 2.0 * rand()/RAND_MAX ); /* sin(delta) uniform in [-1,1] */ /* ===== compute AM-coeffs the 'old way': ===== */ baryinput.site.location[0] = det->location[0]/LAL_C_SI; baryinput.site.location[1] = det->location[1]/LAL_C_SI; baryinput.site.location[2] = det->location[2]/LAL_C_SI; baryinput.alpha = alpha; baryinput.delta = delta; baryinput.dInv = 0.e0; /* amParams structure to compute a(t) and b(t) */ amParams.das = (LALDetAndSource *)LALMalloc(sizeof(LALDetAndSource)); amParams.das->pSource = (LALSource *)LALMalloc(sizeof(LALSource)); amParams.baryinput = &baryinput; amParams.earth = &earth; amParams.edat = &edat; amParams.das->pDetector = det; amParams.das->pSource->equatorialCoords.longitude = alpha; amParams.das->pSource->equatorialCoords.latitude = delta; amParams.das->pSource->orientation = 0.0; amParams.das->pSource->equatorialCoords.system = COORDINATESYSTEM_EQUATORIAL; amParams.polAngle = 0; SUB (LALComputeAM ( &status, &AMold, timestamps->data, &amParams), &status); /* ===== compute AM-coeffs the 'new way' using LALGetAMCoeffs() */ /* ----- get detector-state series ----- */ SUB ( LALGetDetectorStates (&status, &detStates, timestamps, det, &edat, 0 ), &status ); skypos.system = COORDINATESYSTEM_EQUATORIAL; skypos.longitude = alpha; skypos.latitude = delta; SUB ( LALGetAMCoeffs ( &status, &AMnew, detStates, skypos ), &status ); /* ===== analyse relative error ===== */ maxerr_a = maxerr_b = averr_a = averr_b = 0; for ( i=0; i < timestamps->length; i ++ ) { REAL8 thisErr; thisErr = sqrt( SQ ( AMold.a->data[i] - AMnew.a->data[i] ) / AMold.A ); averr_a += thisErr; maxerr_a = MYMAX( thisErr, maxerr_a ); thisErr = sqrt( SQ ( AMold.b->data[i] - AMnew.b->data[i] ) / AMold.B ); averr_b += thisErr; maxerr_b = MYMAX( thisErr, maxerr_b ); } averr_a /= timestamps->length; averr_b /= timestamps->length; if ( lalDebugLevel ) { printf ("Parameters: IFO = %s, skypos = [%g, %g]\n", sites[pickedSite], alpha, delta ); printf ("Maximal relative errors: maxerr(a) = %g %%, maxerr(b) = %g %% \n", 100.0 * maxerr_a, 100.0 * maxerr_b); printf ("Average relative errors: averr(a) = %g %%, averr(b) = %g %% \n", 100.0 * averr_a, 100.0 * averr_b ); } else printf ("%d %g %g %g %g %g %g \n", pickedSite, alpha, delta, averr_a, averr_b, maxerr_a, maxerr_b); if ( (averr_a > tolerance) || (averr_b > tolerance) || (maxerr_a > tolerance) ||(maxerr_b > tolerance)) { XLALPrintError ("Maximal error-tolerance of %g %% was exceeded!\n", 100.0 * tolerance ); return 1; } /* ----- free memory ----- */ XLALDestroyTimestampVector ( timestamps ); XLALDestroyREAL4Vector ( AMold.a ); XLALDestroyREAL4Vector ( AMold.b ); XLALDestroyREAL4Vector ( AMnew.a ); XLALDestroyREAL4Vector ( AMnew.b ); LALFree ( det ); XLALDestroyDetectorStateSeries ( detStates ); LALFree ( amParams.das->pSource ); LALFree ( amParams.das ); LALFree(edat.ephemE); LALFree(edat.ephemS); LALCheckMemoryLeaks(); return 0; /* OK */ } /* main() */
/** * Simulate a pulsar signal to best accuracy possible. * \author Reinhard Prix * \date 2005 * * The motivation for this function is to provide functions to * simulate pulsar signals <em>with the best possible accuracy</em>, * i.e. using no approximations, contrary to LALGeneratePulsarSignal(). * * Obviously this is not meant as a fast code to be used in a Monte-Carlo * simulation, but rather as a <em>reference</em> to compare other (faster) * functions agains, in order to be able to gauge the quality of a given * signal-generation routine. * * We want to calculate \f$h(t)\f$, given by * \f[ * h(t) = F_+(t)\, h_+(t) + F_\times(t) \,h_\times(t)\,, * \f] * where \f$F_+\f$ and \f$F_x\f$ are called the <em>beam-pattern</em> functions, * which depend of the wave polarization \f$\psi\f$, * the source position \f$\alpha\f$, \f$\delta\f$ and the detector position and * orientation (\f$\gamma\f$, \f$\lambda\f$, \f$L\f$ and \f$\xi\f$). The expressions for * the beam-pattern functions are given in \cite JKS98 , which we write as * \f{eqnarray}{ * F_+(t) = \sin \zeta \cos 2\psi \, a(t) + \sin \zeta \sin 2\psi \, b(t)\,,\\ * F_\times(t) = \sin\zeta \cos 2\psi \,b(t) - \sin\zeta \sin 2\psi \, a(t) \,, * \f} * where \f$\zeta\f$ is the angle between the interferometer arms, and * \f{eqnarray}{ * a(t) &=& a_1 \cos[ 2 (\alpha - T)) ] + a_2 \sin[ 2(\alpha - T)] * + a_3 \cos[ \alpha - T ] + a_4 \sin [ \alpha - T ] + a_5\,,\\ * b(t) &=& b_1 \cos[ 2(\alpha - T)] + b_2 \sin[ 2(\alpha - T) ] * + b_3 \cos[ \alpha - T ] + b_4 \sin[ \alpha - T]\,, * \f} * where \f$T\f$ is the local (mean) sidereal time of the detector, and the * time-independent coefficients \f$a_i\f$ and \f$b_i\f$ are given by * \f{eqnarray}{ * a_1 &=& \frac{1}{16} \sin 2\gamma \,(3- \cos 2\lambda)\,(3 - \cos 2\delta)\,,\\ * a_2 &=& -\frac{1}{4}\cos 2\gamma \,\sin \lambda \,(3 - \cos 2\delta) \,,\\ * a_3 &=& \frac{1}{4} \sin 2\gamma \,\sin 2\lambda \,\sin 2\delta \,\\ * a_4 &=& -\frac{1}{2} \cos 2\gamma \,\cos \lambda \,\sin 2 \delta\,,\\ * a_5 &=& \frac{3}{4} \sin 2\gamma \, \cos^2 \lambda \,\cos^2 \delta\,, * \f} * and * \f{eqnarray}{ * b_1 &=& \cos 2\gamma \,\sin \lambda \,\sin \delta\,,\\ * b_2 &=& \frac{1}{4} \sin 2\gamma \,(3-\cos 2\lambda)\, \sin \delta\,,\\ * b_3 &=& \cos 2\gamma \,\cos \lambda \,\cos\delta \,, \\ * b_4 &=& \frac{1}{2} \sin2\gamma \,\sin 2\lambda \,\cos\delta\,, * \f} * * The source model considered is a plane-wave * \f{eqnarray}{ * h_+(t) &=& A_+\, \cos \Psi(t)\,,\\ * h_\times(t) &=& A_\times \, \sin \Psi(t)\,, * \f} * where the wave-phase is \f$\Psi(t) = \Phi_0 + \Phi(t)\f$, and for an * isolated pulsar we have * \f{equation}{ * \Phi(t) = 2\pi \left[\sum_{s=0} \frac{f^{(s)}(\tau_\mathrm{ref})}{ * (s+1)!} \left( \tau(t) - \tau_\mathrm{ref} \right)^{s+1} \right]\,, * \f} * where \f$\tau_\mathrm{ref}\f$ is the "reference time" for the definition * of the pulsar-parameters \f$f^{(s)}\f$ in the solar-system barycenter * (SSB), and \f$\tau(t)\f$ is the SSB-time of the phase arriving at the * detector at UTC-time \f$t\f$, which depends on the source-position * (\f$\alpha\f$, \f$\delta\f$) and the detector-position, namely * \f{equation}{ * \tau (t) = t + \frac{ \vec{r}(t)\cdot\vec{n}}{c}\,, * \f} * where \f$\vec{r}(t)\f$ is the vector from SSB to the detector, and \f$\vec{n}\f$ * is the unit-vector pointing <em>to</em> the source. * * This is a standalone "clean-room" implementation using no other * outside-functions <em>except</em> for LALGPStoLMST1() to calculate * the local (mean) sidereal time at the detector for given GPS-time, * (which I double-checked with an independent Mathematica script), * and and XLALBarycenter() to calculate \f$\tau(t)\f$. * * NOTE: currently only isolated pulsars are supported * * NOTE2: we don't really use the highest possible accuracy right now, * as we blatently neglect all relativistic timing effects (i.e. using dT=v.n/c) * * NOTE3: no heterodyning is performed here, the time-series is generated and sampled * at the given rate, that's all! ==\> the caller needs to make sure about the * right sampling rate to use (-\>aliasing) and do the proper post-treatment... * */ REAL4TimeSeries * XLALSimulateExactPulsarSignal ( const PulsarSignalParams *params ) { XLAL_CHECK_NULL ( params != NULL, XLAL_EINVAL, "Invalid NULL input 'params'\n"); XLAL_CHECK_NULL ( params->samplingRate > 0, XLAL_EDOM, "Sampling rate must be positive, got samplingRate = %g\n", params->samplingRate ); /* don't accept heterodyning frequency */ XLAL_CHECK_NULL ( params->fHeterodyne == 0, XLAL_EINVAL, "Heterodyning frequency must be set to 0, got params->fHeterodyne = %g\n", params->fHeterodyne ); UINT4 numSpins = 3; /* get timestamps of timeseries plus detector-states */ REAL8 dt = 1.0 / params->samplingRate; LIGOTimeGPSVector *timestamps; XLAL_CHECK_NULL ( (timestamps = XLALMakeTimestamps ( params->startTimeGPS, params->duration, dt, 0 )) != NULL, XLAL_EFUNC ); UINT4 numSteps = timestamps->length; DetectorStateSeries *detStates = XLALGetDetectorStates ( timestamps, params->site, params->ephemerides, 0 ); XLAL_CHECK_NULL ( detStates != NULL, XLAL_EFUNC, "XLALGetDetectorStates() failed.\n"); XLALDestroyTimestampVector ( timestamps ); timestamps = NULL; AMCoeffs *amcoe = XLALComputeAMCoeffs ( detStates, params->pulsar.position ); XLAL_CHECK_NULL ( amcoe != NULL, XLAL_EFUNC, "XLALComputeAMCoeffs() failed.\n"); /* create output timeseries (FIXME: should really know *detector* here, not just site!!) */ const LALFrDetector *site = &(params->site->frDetector); CHAR *channel = XLALGetChannelPrefix ( site->name ); XLAL_CHECK_NULL ( channel != NULL, XLAL_EFUNC, "XLALGetChannelPrefix( %s ) failed.\n", site->name ); REAL4TimeSeries *ts = XLALCreateREAL4TimeSeries ( channel, &(detStates->data[0].tGPS), 0, dt, &emptyUnit, numSteps ); XLAL_CHECK_NULL ( ts != NULL, XLAL_EFUNC, "XLALCreateREAL4TimeSeries() failed.\n"); XLALFree ( channel ); channel = NULL; /* orientation of detector arms */ REAL8 xAzi = site->xArmAzimuthRadians; REAL8 yAzi = site->yArmAzimuthRadians; REAL8 Zeta = xAzi - yAzi; if (Zeta < 0) { Zeta = -Zeta; } if ( params->site->type == LALDETECTORTYPE_CYLBAR ) { Zeta = LAL_PI_2; } REAL8 sinZeta = sin(Zeta); /* get source skyposition */ REAL8 Alpha = params->pulsar.position.longitude; REAL8 Delta = params->pulsar.position.latitude; REAL8 vn[3]; vn[0] = cos(Delta) * cos(Alpha); vn[1] = cos(Delta) * sin(Alpha); vn[2] = sin(Delta); /* get spin-parameters (restricted to maximally 3 spindowns right now) */ REAL8 phi0 = params->pulsar.phi0; REAL8 f0 = params->pulsar.f0; REAL8 f1dot = 0, f2dot = 0, f3dot = 0; if ( params->pulsar.spindown && (params->pulsar.spindown->length > numSpins) ) { XLAL_ERROR_NULL ( XLAL_EDOM, "Currently only supports up to %d spindowns!\n", numSpins ); } if ( params->pulsar.spindown && (params->pulsar.spindown->length >= 3 ) ) { f3dot = params->pulsar.spindown->data[2]; } if ( params->pulsar.spindown && (params->pulsar.spindown->length >= 2 ) ) { f2dot = params->pulsar.spindown->data[1]; } if ( params->pulsar.spindown && (params->pulsar.spindown->length >= 1 ) ) { f1dot = params->pulsar.spindown->data[0]; } /* internally we always work with refTime = startTime->SSB, therefore * we need to translate the pulsar spin-params and initial phase to the * startTime */ REAL8 startTimeSSB = XLALGPSGetREAL8 ( &(detStates->data[0].tGPS) ) + SCALAR ( vn, detStates->data[0].rDetector ); REAL8 refTime; if ( params->pulsar.refTime.gpsSeconds != 0 ) { REAL8 refTime0 = XLALGPSGetREAL8 ( &(params->pulsar.refTime) ); REAL8 deltaRef = startTimeSSB - refTime0; LIGOTimeGPS newEpoch; PulsarSpins fkdotNew; XLALGPSSetREAL8( &newEpoch, startTimeSSB ); PulsarSpins XLAL_INIT_DECL(fkdotOld); fkdotOld[0] = f0; fkdotOld[1] = f1dot; fkdotOld[2] = f2dot; fkdotOld[3] = f3dot; REAL8 DeltaTau = XLALGPSDiff ( &newEpoch, &(params->pulsar.refTime) ); int ret = XLALExtrapolatePulsarSpins ( fkdotNew, fkdotOld, DeltaTau ); XLAL_CHECK_NULL ( ret == XLAL_SUCCESS, XLAL_EFUNC, "XLALExtrapolatePulsarSpins() failed.\n"); /* Finally, need to propagate phase */ phi0 += LAL_TWOPI * ( f0 * deltaRef + (1.0/2.0) * f1dot * deltaRef * deltaRef + (1.0/6.0) * f2dot * deltaRef * deltaRef * deltaRef + (1.0/24.0)* f3dot * deltaRef * deltaRef * deltaRef * deltaRef ); f0 = fkdotNew[0]; f1dot = fkdotNew[1]; f2dot = fkdotNew[2]; f3dot = fkdotNew[3]; refTime = startTimeSSB; } /* if refTime given */ else { /* if not given: use startTime -> SSB */ refTime = startTimeSSB; } /* get 4 amplitudes A_\mu */ REAL8 aPlus = sinZeta * params->pulsar.aPlus; REAL8 aCross = sinZeta * params->pulsar.aCross; REAL8 twopsi = 2.0 * params->pulsar.psi; REAL8 A1 = aPlus * cos(phi0) * cos(twopsi) - aCross * sin(phi0) * sin(twopsi); REAL8 A2 = aPlus * cos(phi0) * sin(twopsi) + aCross * sin(phi0) * cos(twopsi); REAL8 A3 = -aPlus * sin(phi0) * cos(twopsi) - aCross * cos(phi0) * sin(twopsi); REAL8 A4 = -aPlus * sin(phi0) * sin(twopsi) + aCross * cos(phi0) * cos(twopsi); /* main loop: generate time-series */ for ( UINT4 i = 0; i < numSteps; i++) { LIGOTimeGPS *tiGPS = &(detStates->data[i].tGPS); REAL8 ti = XLALGPSGetREAL8 ( tiGPS ); REAL8 deltati = ti - refTime; REAL8 dT = SCALAR(vn, detStates->data[i].rDetector ); REAL8 taui = deltati + dT; REAL8 phi_i = LAL_TWOPI * ( f0 * taui + (1.0/2.0) * f1dot * taui*taui + (1.0/6.0) * f2dot * taui*taui*taui + (1.0/24.0)* f3dot * taui*taui*taui*taui ); REAL8 cosphi_i = cos(phi_i); REAL8 sinphi_i = sin(phi_i); REAL8 ai = amcoe->a->data[i]; REAL8 bi = amcoe->b->data[i]; REAL8 hi = A1 * ai * cosphi_i + A2 * bi * cosphi_i + A3 * ai * sinphi_i + A4 * bi * sinphi_i; ts->data->data[i] = (REAL4)hi; } /* for i < numSteps */ XLALDestroyDetectorStateSeries( detStates ); XLALDestroyAMCoeffs ( amcoe ); return ts; } /* XLALSimulateExactPulsarSignal() */