/* ----- function definitions ---------- */ int main ( int argc, char *argv[] ) { LALStatus status; UserInput_t uvar_s; UserInput_t *uvar = &uvar_s; INIT_MEM ( status ); INIT_MEM ( uvar_s ); struct tms buf; uvar->randSeed = times(&buf); // ---------- register all our user-variable ---------- XLALregBOOLUserStruct ( help, 'h', UVAR_HELP , "Print this help/usage message"); XLALregINTUserStruct ( randSeed, 's', UVAR_OPTIONAL, "Specify random-number seed for reproducible noise."); /* read cmdline & cfgfile */ XLAL_CHECK ( XLALUserVarReadAllInput ( argc, argv ) == XLAL_SUCCESS, XLAL_EFUNC ); if ( uvar->help ) { /* if help was requested, we're done */ exit (0); } srand ( uvar->randSeed ); REAL8 startTimeREAL8 = 714180733; REAL8 duration = 180000; /* 50 hours */ REAL8 Tsft = 1800; /* assume 30min SFTs */ char earthEphem[] = TEST_DATA_DIR "earth00-19-DE200.dat.gz"; char sunEphem[] = TEST_DATA_DIR "sun00-19-DE200.dat.gz"; //REAL8 tolerance = 2e-10; /* same algorithm, should be basically identical results */ LIGOTimeGPS startTime, refTime; XLALGPSSetREAL8 ( &startTime, startTimeREAL8 ); refTime = startTime; // pick skyposition at random ----- */ SkyPosition skypos; skypos.longitude = LAL_TWOPI * (1.0 * rand() / ( RAND_MAX + 1.0 ) ); // alpha uniform in [0, 2pi) skypos.latitude = LAL_PI_2 - acos ( 1 - 2.0 * rand()/RAND_MAX ); // sin(delta) uniform in [-1,1] skypos.system = COORDINATESYSTEM_EQUATORIAL; // pick binary orbital parameters: // somewhat inspired by Sco-X1 parameters from S2-paper (PRD76, 082001 (2007), gr-qc/0605028) // but with a more extreme eccentricity, and random argp REAL8 argp = LAL_TWOPI * (1.0 * rand() / ( RAND_MAX + 1.0 ) ); // uniform in [0, 2pi) BinaryOrbitParams orbit; XLALGPSSetREAL8 ( &orbit.tp, 731163327 ); // time of observed periapsis passage (in SSB) orbit.argp = argp; // argument of periapsis (radians) orbit.asini = 1.44; // projected, normalized orbital semi-major axis (s) */ orbit.ecc = 1e-2; // relatively large value, for better testing orbit.period = 68023; // period (s) : about ~18.9h // ----- step 0: prepare test-case input for calling the BinarySSB-functions // setup detectors const char *sites[3] = { "H1", "L1", "V1" }; UINT4 numDetectors = sizeof( sites ) / sizeof ( sites[0] ); MultiLALDetector multiIFO; multiIFO.length = numDetectors; for ( UINT4 X = 0; X < numDetectors; X ++ ) { LALDetector *det = XLALGetSiteInfo ( sites[X] ); XLAL_CHECK ( det != NULL, XLAL_EFUNC, "XLALGetSiteInfo ('%s') failed for detector X=%d\n", sites[X], X ); multiIFO.sites[X] = (*det); // struct copy XLALFree ( det ); } // load ephemeris EphemerisData *edat = XLALInitBarycenter ( earthEphem, sunEphem ); XLAL_CHECK ( edat != NULL, XLAL_EFUNC, "XLALInitBarycenter('%s','%s') failed\n", earthEphem, sunEphem ); // setup multi-timeseries MultiLIGOTimeGPSVector *multiTS; XLAL_CHECK ( (multiTS = XLALCalloc ( 1, sizeof(*multiTS))) != NULL, XLAL_ENOMEM ); XLAL_CHECK ( (multiTS->data = XLALCalloc (numDetectors, sizeof(*multiTS->data))) != NULL, XLAL_ENOMEM ); multiTS->length = numDetectors; for ( UINT4 X = 0; X < numDetectors; X ++ ) { multiTS->data[X] = XLALMakeTimestamps ( startTime, duration, Tsft, 0 ); XLAL_CHECK ( multiTS->data[X] != NULL, XLAL_EFUNC, "XLALMakeTimestamps() failed.\n"); } /* for X < numIFOs */ // generate detector-states MultiDetectorStateSeries *multiDetStates = XLALGetMultiDetectorStates ( multiTS, &multiIFO, edat, 0 ); XLAL_CHECK ( multiDetStates != NULL, XLAL_EFUNC, "XLALGetMultiDetectorStates() failed.\n"); // generate isolated-NS SSB times MultiSSBtimes *multiSSBIn = XLALGetMultiSSBtimes ( multiDetStates, skypos, refTime, SSBPREC_RELATIVISTICOPT ); XLAL_CHECK ( multiSSBIn != NULL, XLAL_EFUNC, "XLALGetMultiSSBtimes() failed.\n"); // ----- step 1: compute reference-result using old LALGetMultiBinarytimes() MultiSSBtimes *multiBinary_ref = NULL; LALGetMultiBinarytimes (&status, &(multiBinary_ref), multiSSBIn, multiDetStates, &orbit, refTime ); XLAL_CHECK ( status.statusCode == 0, XLAL_EFAILED, "LALGetMultiBinarytimes() failed with status = %d : '%s'\n", status.statusCode, status.statusDescription ); // ----- step 2: compute test-result using new XLALAddMultiBinaryTimes() MultiSSBtimes *multiBinary_test = NULL; PulsarDopplerParams doppler; memset(&doppler, 0, sizeof(doppler)); doppler.tp = orbit.tp; doppler.argp = orbit.argp; doppler.asini = orbit.asini; doppler.ecc = orbit.ecc; doppler.period = orbit.period; XLAL_CHECK ( XLALAddMultiBinaryTimes ( &multiBinary_test, multiSSBIn, &doppler ) == XLAL_SUCCESS, XLAL_EFUNC ); // ----- step 3: compare results REAL8 err_DeltaT, err_Tdot; REAL8 tolerance = 1e-10; int ret = XLALCompareMultiSSBtimes ( &err_DeltaT, &err_Tdot, multiBinary_ref, multiBinary_test ); XLAL_CHECK ( ret == XLAL_SUCCESS, XLAL_EFUNC, "XLALCompareMultiSSBtimes() failed.\n"); XLALPrintWarning ( "INFO: err(DeltaT) = %g, err(Tdot) = %g\n", err_DeltaT, err_Tdot ); XLAL_CHECK ( err_DeltaT < tolerance, XLAL_ETOL, "error(DeltaT) = %g exceeds tolerance of %g\n", err_DeltaT, tolerance ); XLAL_CHECK ( err_Tdot < tolerance, XLAL_ETOL, "error(Tdot) = %g exceeds tolerance of %g\n", err_Tdot, tolerance ); // ---- step 4: clean-up memory XLALDestroyUserVars(); XLALDestroyEphemerisData ( edat ); XLALDestroyMultiSSBtimes ( multiBinary_test ); XLALDestroyMultiSSBtimes ( multiBinary_ref ); XLALDestroyMultiSSBtimes ( multiSSBIn ); XLALDestroyMultiTimestamps ( multiTS ); XLALDestroyMultiDetectorStateSeries ( multiDetStates ); // check for memory-leaks LALCheckMemoryLeaks(); return XLAL_SUCCESS; } // main()
/** * Handle user-input and set up shop accordingly, and do all * consistency-checks on user-input. */ int XLALInitMakefakedata ( ConfigVars_t *cfg, UserVariables_t *uvar ) { XLAL_CHECK ( cfg != NULL, XLAL_EINVAL, "Invalid NULL input 'cfg'\n" ); XLAL_CHECK ( uvar != NULL, XLAL_EINVAL, "Invalid NULL input 'uvar'\n"); cfg->VCSInfoString = XLALGetVersionString(0); XLAL_CHECK ( cfg->VCSInfoString != NULL, XLAL_EFUNC, "XLALGetVersionString(0) failed.\n" ); // version info was requested: output then exit if ( uvar->version ) { printf ("%s\n", cfg->VCSInfoString ); exit (0); } /* if requested, log all user-input and code-versions */ if ( uvar->logfile ) { XLAL_CHECK ( XLALWriteMFDlog ( uvar->logfile, cfg ) == XLAL_SUCCESS, XLAL_EFUNC, "XLALWriteMFDlog() failed with xlalErrno = %d\n", xlalErrno ); } /* Init ephemerides */ XLAL_CHECK ( (cfg->edat = XLALInitBarycenter ( uvar->ephemEarth, uvar->ephemSun )) != NULL, XLAL_EFUNC ); /* check for negative fMin and Band, which would break the fMin_eff, fBand_eff calculation below */ XLAL_CHECK ( uvar->fmin >= 0, XLAL_EDOM, "Invalid negative frequency fMin=%f!\n\n", uvar->fmin ); XLAL_CHECK ( uvar->Band > 0, XLAL_EDOM, "Invalid non-positive frequency band Band=%f!\n\n", uvar->Band ); // ---------- check user-input consistency ---------- // ----- check if frames + frame channels given BOOLEAN have_frames = (uvar->inFrames != NULL); BOOLEAN have_channels= (uvar->inFrChannels != NULL); XLAL_CHECK ( !(have_frames || have_channels) || (have_frames && have_channels), XLAL_EINVAL, "Need both --inFrames and --inFrChannels, or NONE\n"); // ----- IFOs : only from one of {--IFOs, --noiseSFTs, --inFrChannels}: mutually exclusive BOOLEAN have_IFOs = (uvar->IFOs != NULL); BOOLEAN have_noiseSFTs = (uvar->noiseSFTs != NULL); XLAL_CHECK ( have_frames || have_IFOs || have_noiseSFTs, XLAL_EINVAL, "Need one of --IFOs, --noiseSFTs or --inFrChannels to determine detectors\n"); if ( have_frames ) { XLAL_CHECK ( !have_IFOs && !have_noiseSFTs, XLAL_EINVAL, "If --inFrames given, cannot handle --IFOs or --noiseSFTs input\n"); XLAL_CHECK ( XLALParseMultiLALDetector ( &(cfg->multiIFO), uvar->inFrChannels ) == XLAL_SUCCESS, XLAL_EFUNC ); } else { // !have_frames XLAL_CHECK ( !(have_IFOs && have_noiseSFTs), XLAL_EINVAL, "Cannot handle both --IFOs and --noiseSFTs input\n"); } if ( have_IFOs ) { XLAL_CHECK ( XLALParseMultiLALDetector ( &(cfg->multiIFO), uvar->IFOs ) == XLAL_SUCCESS, XLAL_EFUNC ); } // ----- TIMESTAMPS: either from --timestampsFiles, --startTime+duration, or --noiseSFTs BOOLEAN have_startTime = XLALUserVarWasSet ( &uvar->startTime ); BOOLEAN have_duration = XLALUserVarWasSet ( &uvar->duration ); BOOLEAN have_timestampsFiles = ( uvar->timestampsFiles != NULL ); // need BOTH startTime+duration or none XLAL_CHECK ( ( have_duration && have_startTime) || !( have_duration || have_startTime ), XLAL_EINVAL, "Need BOTH {--startTime,--duration} or NONE\n"); // at least one of {startTime,timestamps,noiseSFTs,inFrames} required XLAL_CHECK ( have_timestampsFiles || have_startTime || have_noiseSFTs || have_frames, XLAL_EINVAL, "Need at least one of {--timestampsFiles, --startTime+duration, --noiseSFTs, --inFrames}\n" ); // don't allow timestamps + {startTime+duration OR noiseSFTs} XLAL_CHECK ( !have_timestampsFiles || !(have_startTime||have_noiseSFTs), XLAL_EINVAL, "--timestampsFiles incompatible with {--noiseSFTs or --startTime+duration}\n"); // note, however, that we DO allow --noiseSFTs and --startTime+duration, which will act as a constraint // on the noise-SFTs to load in // don't allow --SFToverlap with either --noiseSFTs OR --timestampsFiles XLAL_CHECK ( uvar->SFToverlap >= 0, XLAL_EDOM ); BOOLEAN haveOverlap = ( uvar->SFToverlap > 0 ); XLAL_CHECK ( !haveOverlap || !( have_noiseSFTs || have_timestampsFiles ), XLAL_EINVAL, "--SFToverlap incompatible with {--noiseSFTs or --timestampsFiles}\n" ); // now handle the 3 mutually-exclusive cases: have_noiseSFTs || have_timestampsFiles || have_startTime (only) if ( have_noiseSFTs ) { SFTConstraints XLAL_INIT_DECL(constraints); if ( have_startTime && have_duration ) // use optional (startTime+duration) as constraints, { LIGOTimeGPS minStartTime, maxStartTime; minStartTime = uvar->startTime; maxStartTime = uvar->startTime; XLALGPSAdd ( &maxStartTime, uvar->duration ); constraints.minStartTime = &minStartTime; constraints.maxStartTime = &maxStartTime; char bufGPS1[32], bufGPS2[32]; XLALPrintWarning ( "Only noise-SFTs between GPS [%s, %s] will be used!\n", XLALGPSToStr(bufGPS1, &minStartTime), XLALGPSToStr(bufGPS2, &maxStartTime) ); } /* if start+duration given */ XLAL_CHECK ( (cfg->noiseCatalog = XLALSFTdataFind ( uvar->noiseSFTs, &constraints )) != NULL, XLAL_EFUNC ); XLAL_CHECK ( cfg->noiseCatalog->length > 0, XLAL_EINVAL, "No noise-SFTs matching (start+duration, timestamps) were found!\n" ); XLAL_CHECK ( (cfg->multiNoiseCatalogView = XLALGetMultiSFTCatalogView ( cfg->noiseCatalog )) != NULL, XLAL_EFUNC ); // extract multi-timestamps from the multi-SFT-catalog view XLAL_CHECK ( (cfg->multiTimestamps = XLALTimestampsFromMultiSFTCatalogView ( cfg->multiNoiseCatalogView )) != NULL, XLAL_EFUNC ); // extract IFOs from multi-SFT catalog XLAL_CHECK ( XLALMultiLALDetectorFromMultiSFTCatalogView ( &(cfg->multiIFO), cfg->multiNoiseCatalogView ) == XLAL_SUCCESS, XLAL_EFUNC ); } // endif have_noiseSFTs else if ( have_timestampsFiles ) { XLAL_CHECK ( (cfg->multiTimestamps = XLALReadMultiTimestampsFiles ( uvar->timestampsFiles )) != NULL, XLAL_EFUNC ); XLAL_CHECK ( (cfg->multiTimestamps->length > 0) && (cfg->multiTimestamps->data != NULL), XLAL_EINVAL, "Got empty timestamps-list from XLALReadMultiTimestampsFiles()\n" ); for ( UINT4 X=0; X < cfg->multiTimestamps->length; X ++ ) { cfg->multiTimestamps->data[X]->deltaT = uvar->Tsft; // Tsft information not given by timestamps-file } } // endif have_timestampsFiles else if ( have_startTime && have_duration ) { XLAL_CHECK ( ( cfg->multiTimestamps = XLALMakeMultiTimestamps ( uvar->startTime, uvar->duration, uvar->Tsft, uvar->SFToverlap, cfg->multiIFO.length )) != NULL, XLAL_EFUNC ); } // endif have_startTime // check if the user asked for Gaussian white noise to be produced (sqrtSn[X]!=0), otherwise leave noise-floors at 0 if ( uvar->sqrtSX != NULL ) { XLAL_CHECK ( XLALParseMultiNoiseFloor ( &(cfg->multiNoiseFloor), uvar->sqrtSX, cfg->multiIFO.length ) == XLAL_SUCCESS, XLAL_EFUNC ); } else { cfg->multiNoiseFloor.length = cfg->multiIFO.length; // values remain at their default sqrtSn[X] = 0; } #ifdef HAVE_LIBLALFRAME // if user requested time-series data from frames to be added: try to read the frames now if ( have_frames ) { UINT4 numDetectors = uvar->inFrChannels->length; XLAL_CHECK ( uvar->inFrames->length == numDetectors, XLAL_EINVAL, "Need equal number of channel names (%d) as frame specifications (%d)\n", uvar->inFrChannels->length, numDetectors ); XLAL_CHECK ( (cfg->inputMultiTS = XLALCalloc ( 1, sizeof(*cfg->inputMultiTS))) != NULL, XLAL_ENOMEM ); cfg->inputMultiTS->length = numDetectors; XLAL_CHECK ( (cfg->inputMultiTS->data = XLALCalloc ( numDetectors, sizeof(cfg->inputMultiTS->data[0]) )) != NULL, XLAL_ENOMEM ); if ( cfg->multiTimestamps == NULL ) { XLAL_CHECK ( (cfg->multiTimestamps = XLALCalloc ( 1, sizeof(*cfg->multiTimestamps) )) != NULL, XLAL_ENOMEM ); XLAL_CHECK ( (cfg->multiTimestamps->data = XLALCalloc ( numDetectors, sizeof(cfg->multiTimestamps->data[0]))) != NULL, XLAL_ENOMEM ); cfg->multiTimestamps->length = numDetectors; } for ( UINT4 X = 0; X < numDetectors; X ++ ) { LALCache *cache; XLAL_CHECK ( (cache = XLALCacheImport ( uvar->inFrames->data[X] )) != NULL, XLAL_EFUNC, "Failed to import cache file '%s'\n", uvar->inFrames->data[X] ); // this is a sorted cache, so extract its time-range: REAL8 cache_tStart = cache->list[0].t0; REAL8 cache_tEnd = cache->list[cache->length-1].t0 + cache->list[cache->length-1].dt; REAL8 cache_duration = (cache_tEnd - cache_tStart); LIGOTimeGPS ts_start; REAL8 ts_duration; // check that it's consistent with timestamps, if given, otherwise create timestamps from this if ( cfg->multiTimestamps->data[X] != NULL ) // FIXME: implicitly assumes timestamps are sorted, which is not guaranteed by timestamps-reading from file { const LIGOTimeGPSVector *timestampsX = cfg->multiTimestamps->data[X]; REAL8 tStart = XLALGPSGetREAL8( ×tampsX->data[0] ); REAL8 tEnd = XLALGPSGetREAL8( ×tampsX->data[timestampsX->length-1]) + timestampsX->deltaT; XLAL_CHECK ( tStart >= cache_tStart && tEnd <= cache_tEnd, XLAL_EINVAL, "Detector X=%d: Requested timestamps-range [%.0f, %.0f]s outside of cache range [%.0f,%.0f]s\n", X, tStart, tEnd, cache_tStart, cache_tEnd ); XLALGPSSetREAL8 ( &ts_start, tStart ); ts_duration = (tEnd - tStart); } else { XLALGPSSetREAL8 ( &ts_start, (REAL8)cache_tStart + 1); // cache times can apparently be by rounded up or down by 1s, so shift by 1s to be safe ts_duration = cache_duration - 1; XLAL_CHECK ( (cfg->multiTimestamps->data[X] = XLALMakeTimestamps ( ts_start, ts_duration, uvar->Tsft, uvar->SFToverlap ) ) != NULL, XLAL_EFUNC ); } // ----- now open frame stream and read *all* the data within this time-range [FIXME] ---------- LALFrStream *stream; XLAL_CHECK ( (stream = XLALFrStreamCacheOpen ( cache )) != NULL, XLAL_EFUNC, "Failed to open stream from cache file '%s'\n", uvar->inFrames->data[X] ); XLALDestroyCache ( cache ); const char *channel = uvar->inFrChannels->data[X]; size_t limit = 0; // unlimited read REAL8TimeSeries *ts; XLAL_CHECK ( (ts = XLALFrStreamInputREAL8TimeSeries ( stream, channel, &ts_start, ts_duration, limit )) != NULL, XLAL_EFUNC, "Frame reading failed for stream created for '%s': ts_start = {%d,%d}, duration=%.0f\n", uvar->inFrames->data[X], ts_start.gpsSeconds, ts_start.gpsNanoSeconds, ts_duration ); cfg->inputMultiTS->data[X] = ts; XLAL_CHECK ( XLALFrStreamClose ( stream ) == XLAL_SUCCESS, XLAL_EFUNC, "Stream closing failed for cache file '%s'\n", uvar->inFrames->data[X] ); } // for X < numDetectors } // if inFrames // if user requested timeseries *output* to frame files, handle deprecated options XLAL_CHECK ( !(uvar->TDDframedir && uvar->outFrameDir), XLAL_EINVAL, "Specify only ONE of {--TDDframedir or --outFrameDir} or NONE\n"); if ( uvar->TDDframedir ) { cfg->outFrameDir = uvar->TDDframedir; } else if ( uvar->outFrameDir ) { cfg->outFrameDir = uvar->outFrameDir; } #endif return XLAL_SUCCESS; } /* XLALInitMakefakedata() */
/** * Handle user-input and check its validity. * Load ephemeris and calculate AM-coefficients (stored globally) */ void Initialize (LALStatus *status, struct CommandLineArgsTag *CLA) { EphemerisData *edat=NULL; /* Stores earth/sun ephemeris data for barycentering */ BarycenterInput baryinput; /* Stores detector location and other barycentering data */ EarthState earth; AMCoeffsParams *amParams; LIGOTimeGPS *midTS=NULL; /* Time stamps for amplitude modulation coefficients */ LALDetector *Detector; /* Our detector*/ INT4 k; INITSTATUS(status); ATTATCHSTATUSPTR (status); if ( XLALUserVarWasSet ( &(CLA->nTsft) ) ) CLA->duration = 1.0 * CLA->nTsft * CLA->Tsft; /* read or generate SFT timestamps */ if ( XLALUserVarWasSet(&(CLA->timestamps)) ) { XLAL_CHECK_LAL ( status, ( timestamps = XLALReadTimestampsFile ( CLA->timestamps ) ) != NULL, XLAL_EFUNC ); if ( (CLA->nTsft > 0) && ( (UINT4)CLA->nTsft < timestamps->length ) ) /* truncate if required */ timestamps->length = CLA->nTsft; CLA->nTsft = timestamps->length; } /* if have_timestamps */ else { LIGOTimeGPS tStart; tStart.gpsSeconds = CLA->gpsStart; tStart.gpsNanoSeconds = 0; XLAL_CHECK_LAL ( status, ( timestamps = XLALMakeTimestamps( tStart, CLA->duration, CLA->Tsft, 0 ) ) != NULL, XLAL_EFUNC ); CLA->nTsft = timestamps->length; } /* no timestamps */ /*---------- initialize detector ---------- */ { BOOLEAN have_IFO = XLALUserVarWasSet ( &CLA->IFO ); BOOLEAN have_detector = XLALUserVarWasSet ( &CLA->detector ); CHAR *IFO; if ( !have_IFO && !have_detector ) { fprintf (stderr, "\nNeed to specify the detector (--IFO) !\n\n"); ABORT (status, SEMIANALYTIC_EINPUT, SEMIANALYTIC_MSGEINPUT); } if ( have_IFO ) IFO = CLA->IFO; else IFO = CLA->detector; if ( ( Detector = XLALGetSiteInfo ( IFO ) ) == NULL ) { ABORT (status, SEMIANALYTIC_EINPUT, SEMIANALYTIC_MSGEINPUT); } } /* ---------- load ephemeris-files ---------- */ { edat = XLALInitBarycenter( CLA->ephemEarth, CLA->ephemSun ); if ( !edat ) { XLALPrintError("XLALInitBarycenter failed: could not load Earth ephemeris '%s' and Sun ephemeris '%s'\n", CLA->ephemEarth, CLA->ephemSun); ABORT (status, SEMIANALYTIC_EINPUT, SEMIANALYTIC_MSGEINPUT); } } /* ephemeris-reading */ /* ---------- calculate AM-coefficients ---------- */ /* prepare call to barycentering routing */ baryinput.site.location[0] = Detector->location[0]/LAL_C_SI; baryinput.site.location[1] = Detector->location[1]/LAL_C_SI; baryinput.site.location[2] = Detector->location[2]/LAL_C_SI; baryinput.alpha = CLA->Alpha; baryinput.delta = CLA->Delta; baryinput.dInv = 0.e0; /* amParams structure to compute a(t) and b(t) */ /* Allocate space for amParams stucture */ /* Here, amParams->das is the Detector and Source info */ amParams = (AMCoeffsParams *)LALMalloc(sizeof(AMCoeffsParams)); amParams->das = (LALDetAndSource *)LALMalloc(sizeof(LALDetAndSource)); amParams->das->pSource = (LALSource *)LALMalloc(sizeof(LALSource)); /* Fill up AMCoeffsParams structure */ amParams->baryinput = &baryinput; amParams->earth = &earth; amParams->edat = edat; amParams->das->pDetector = Detector; amParams->das->pSource->equatorialCoords.system = COORDINATESYSTEM_EQUATORIAL; amParams->das->pSource->equatorialCoords.longitude = CLA->Alpha; amParams->das->pSource->equatorialCoords.latitude = CLA->Delta; amParams->das->pSource->orientation = 0.0; amParams->polAngle = amParams->das->pSource->orientation ; /* These two have to be the same!!!!!!!!!*/ /* Allocate space for AMCoeffs */ XLAL_INIT_MEM(amc); TRY ( LALSCreateVector(status->statusPtr, &(amc.a), (UINT4) CLA->nTsft), status); TRY ( LALSCreateVector(status->statusPtr, &(amc.b), (UINT4) CLA->nTsft), status); /* Mid point of each SFT */ midTS = (LIGOTimeGPS *)LALCalloc(CLA->nTsft,sizeof(LIGOTimeGPS)); for(k=0; k < CLA->nTsft; k++) { /* FIXME: loss of precision; consider midTS[k] = timestamps->data[k]; XLALGPSAdd(&midTS[k], 0.5*CLA->Tsft); */ REAL8 teemp=0.0; teemp = XLALGPSGetREAL8(&(timestamps->data[k])); teemp += 0.5*CLA->Tsft; XLALGPSSetREAL8(&(midTS[k]), teemp); } TRY ( LALComputeAM(status->statusPtr, &amc, midTS, amParams), status); /* Free memory */ XLALDestroyTimestampVector ( timestamps); LALFree(midTS); LALFree(Detector); XLALDestroyEphemerisData(edat); LALFree(amParams->das->pSource); LALFree(amParams->das); LALFree(amParams); DETATCHSTATUSPTR (status); RETURN(status); } /* ParseUserInput() */
// ---------- main ---------- int main ( int argc, char *argv[] ) { XLAL_CHECK ( argc == 1, XLAL_EINVAL, "No input arguments allowed.\n" ); XLAL_CHECK ( argv != NULL, XLAL_EINVAL ); // ----- load ephemeris EphemerisData *ephem; XLAL_CHECK ( (ephem = XLALInitBarycenter ( TEST_DATA_DIR "earth00-19-DE405.dat.gz", TEST_DATA_DIR "sun00-19-DE405.dat.gz" )) != NULL, XLAL_EFUNC ); // ----- setup injection and data parameters LALStringVector *detNames = NULL; XLAL_CHECK ( (detNames = XLALCreateStringVector ( "H1", "L1", NULL )) != NULL, XLAL_EFUNC ); UINT4 numDetectors = detNames->length; // generate and assume some gaussian noise floors MultiNoiseFloor XLAL_INIT_DECL(injectSqrtSX); MultiNoiseFloor XLAL_INIT_DECL(assumeSqrtSX); injectSqrtSX.length = numDetectors; assumeSqrtSX.length = numDetectors; for ( UINT4 X = 0; X < numDetectors; X ++ ) { injectSqrtSX.sqrtSn[X] = 0; // don't inject random noise to keep errors deterministic and informative (resampling differs much more on noise) assumeSqrtSX.sqrtSn[X] = 1.0 + 2.0*X; } LIGOTimeGPS startTime = {711595934, 0}; REAL8 Tspan = 20 * 3600; LIGOTimeGPS endTime = startTime; XLALGPSAdd( &endTime, Tspan ); REAL8 Tsft = 1800; LIGOTimeGPS refTime = { startTime.gpsSeconds - 2.3 * Tspan, 0 }; MultiLIGOTimeGPSVector *multiTimestamps; XLAL_CHECK ( ( multiTimestamps = XLALCalloc ( 1, sizeof(*multiTimestamps))) != NULL, XLAL_ENOMEM ); XLAL_CHECK ( ( multiTimestamps->data = XLALCalloc ( numDetectors, sizeof(multiTimestamps->data[0]) )) != NULL, XLAL_ENOMEM ); multiTimestamps->length = numDetectors; LIGOTimeGPS startTimeX = startTime; for ( UINT4 X=0; X < numDetectors; X ++ ) { XLAL_CHECK ( (multiTimestamps->data[X] = XLALMakeTimestamps ( startTimeX, Tspan, Tsft, 0 ) ) != NULL, XLAL_EFUNC ); XLALGPSAdd ( &startTimeX, 0.5 * Tspan ); // shift start-times by 1/2 Tspan for each detector Tspan *= 2.0; } // for X < numDetectors // shift a few timestamps around to create gaps UINT4 numSFTsPerDet = multiTimestamps->data[0]->length; multiTimestamps->data[0]->data[numSFTsPerDet-1].gpsSeconds += 10000; multiTimestamps->data[0]->data[numSFTsPerDet-2].gpsSeconds += 5000; multiTimestamps->data[1]->data[0].gpsSeconds -= 10000; multiTimestamps->data[1]->data[1].gpsSeconds -= 5000; SFTCatalog *catalog; XLAL_CHECK ( (catalog = XLALMultiAddToFakeSFTCatalog ( NULL, detNames, multiTimestamps )) != NULL, XLAL_EFUNC ); // ----- CW sources to injet ---------- REAL8 Freq = 100.0; PulsarParamsVector *injectSources; XLAL_CHECK ( (injectSources = XLALCreatePulsarParamsVector(1)) != NULL, XLAL_EFUNC ); injectSources->data[0].Amp.h0 = 1; injectSources->data[0].Amp.cosi = 0.5; injectSources->data[0].Amp.psi = 0.1; injectSources->data[0].Amp.phi0 = 1.2; REAL8 asini = 0; // 1.4; // sco-X1 like REAL8 Period = 0; // 19 * 3600;// sco-X1 like REAL8 ecc = 0; // 0.1; // much larger than ScoX1 PulsarDopplerParams XLAL_INIT_DECL(Doppler); Doppler.Alpha = 0.5; Doppler.Delta = -0.5; Doppler.fkdot[0] = Freq; Doppler.fkdot[1] = -1e-9; Doppler.refTime = refTime; Doppler.asini = asini; Doppler.ecc = ecc; Doppler.tp = startTime; Doppler.period = Period; Doppler.argp = 0.5; injectSources->data[0].Doppler = Doppler; REAL8 dFreq = 0.1 / Tspan; // 10x finer than native FFT resolution REAL8 mis = 0.5; REAL8 df1dot = sqrt( 720.0 * mis ) / (LAL_PI * Tspan * Tspan); // metric (f-projected) stepsize for given mismatch mis REAL8 dSky = 1e4 / (Freq * Tspan); // rough estimate of a 'metric' sky step, eg. Eq.(118) in \cite Prix07 REAL8 dPeriod = 3600; UINT4 numFreqBins = 1000; UINT4 numf1dotPoints = 2; UINT4 numSkyPoints = 2; UINT4 numPeriodPoints = 2; PulsarSpinRange XLAL_INIT_DECL(spinRange); spinRange.refTime = refTime; memcpy ( &spinRange.fkdot, &injectSources->data[0].Doppler.fkdot, sizeof(spinRange.fkdot) ); spinRange.fkdotBand[0] = (numFreqBins - 1)*dFreq - 10*LAL_REAL8_EPS; spinRange.fkdotBand[1] = (numf1dotPoints - 1)*df1dot - 10*LAL_REAL8_EPS; Doppler.fkdot[0] -= 0.4 * spinRange.fkdotBand[0]; REAL8 minCoverFreq, maxCoverFreq; XLAL_CHECK ( XLALCWSignalCoveringBand ( &minCoverFreq, &maxCoverFreq, &startTime, &endTime, &spinRange, asini, Period, ecc ) == XLAL_SUCCESS, XLAL_EFUNC ); // ----- setup optional Fstat arguments FstatOptionalArgs optionalArgs = FstatOptionalArgsDefaults; optionalArgs.injectSources = injectSources; optionalArgs.injectSqrtSX = &injectSqrtSX; optionalArgs.assumeSqrtSX = &assumeSqrtSX; // ----- prepare input data with injection for all available methods FstatInput *input[FMETHOD_END]; FstatResults *results[FMETHOD_END]; for ( UINT4 iMethod = FMETHOD_START; iMethod < FMETHOD_END; iMethod ++ ) { results[iMethod] = NULL; if ( !XLALFstatMethodIsAvailable(iMethod) ) { continue; } optionalArgs.FstatMethod = iMethod; XLAL_CHECK ( (input[iMethod] = XLALCreateFstatInput ( catalog, minCoverFreq, maxCoverFreq, dFreq, ephem, &optionalArgs )) != NULL, XLAL_EFUNC ); } FstatQuantities whatToCompute = (FSTATQ_2F | FSTATQ_FAFB); // ----- loop over all templates {sky, f1dot, period} for ( UINT4 iSky = 0; iSky < numSkyPoints; iSky ++ ) { for ( UINT4 if1dot = 0; if1dot < numf1dotPoints; if1dot ++ ) { for ( UINT4 iPeriod = 0; iPeriod < numPeriodPoints; iPeriod ++ ) { // ----- loop over all available methods and compare Fstat results FstatMethodType firstMethod = FMETHOD_START; for ( UINT4 iMethod = FMETHOD_START; iMethod < FMETHOD_END; iMethod ++ ) { if ( !XLALFstatMethodIsAvailable(iMethod) ) { continue; } if ( firstMethod == FMETHOD_START ) { // keep track of first available method found firstMethod = iMethod; } XLAL_CHECK ( XLALComputeFstat ( &results[iMethod], input[iMethod], &Doppler, numFreqBins, whatToCompute ) == XLAL_SUCCESS, XLAL_EFUNC ); if ( lalDebugLevel & LALINFOBIT ) { FILE *fp; char fname[1024]; XLAL_INIT_MEM ( fname ); snprintf ( fname, sizeof(fname)-1, "twoF%s-iSky%02d-if1dot%02d-iPeriod%02d.dat", XLALGetFstatMethodName(iMethod), iSky, if1dot, iPeriod ); XLAL_CHECK ( (fp = fopen ( fname, "wb" )) != NULL, XLAL_EFUNC ); for ( UINT4 k = 0; k < results[iMethod]->numFreqBins; k ++ ) { REAL8 Freq0 = results[iMethod]->doppler.fkdot[0]; REAL8 Freq_k = Freq0 + k * results[iMethod]->dFreq; if ( whatToCompute & FSTATQ_FAFB ) { fprintf ( fp, "%20.16g %10.4g %10.4g %10.4g %10.4g %10.4g\n", Freq_k, results[iMethod]->twoF[k], crealf(results[iMethod]->Fa[k]), cimagf(results[iMethod]->Fa[k]), crealf(results[iMethod]->Fb[k]), cimagf(results[iMethod]->Fb[k]) ); } else { fprintf ( fp, "%20.16g %10.4g\n", Freq_k, results[iMethod]->twoF[k] ); } } // for k < numFreqBins fclose(fp); } // if info // compare to first result if ( iMethod != firstMethod ) { XLALPrintInfo ("Comparing results between method '%s' and '%s'\n", XLALGetFstatMethodName(firstMethod), XLALGetFstatMethodName(iMethod) ); if ( compareFstatResults ( results[firstMethod], results[iMethod] ) != XLAL_SUCCESS ) { XLALPrintError ("Comparison between method '%s' and '%s' failed\n", XLALGetFstatMethodName(firstMethod), XLALGetFstatMethodName(iMethod) ); XLAL_ERROR ( XLAL_EFUNC ); } } } // for i < FMETHOD_END Doppler.period += dPeriod; } // for iPeriod < numPeriodPoints Doppler.fkdot[1] += df1dot; } // for if1dot < numf1dotPoints Doppler.Alpha += dSky; } // for iSky < numSkyPoints // free remaining memory for ( UINT4 iMethod=FMETHOD_START; iMethod < FMETHOD_END; iMethod ++ ) { if ( !XLALFstatMethodIsAvailable(iMethod) ) { continue; } XLALDestroyFstatInput ( input[iMethod] ); XLALDestroyFstatResults ( results[iMethod] ); } // for i < FMETHOD_END XLALDestroyPulsarParamsVector ( injectSources ); XLALDestroySFTCatalog ( catalog ); XLALDestroyMultiTimestamps ( multiTimestamps ); XLALDestroyStringVector ( detNames ); XLALDestroyEphemerisData ( ephem ); LALCheckMemoryLeaks(); return XLAL_SUCCESS; } // 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() */