/* void RunGeneratePulsarSignalTest(LALStatus *status, int argc, char **argv) */ /* 02/02/05 gam */ void RunGeneratePulsarSignalTest(LALStatus *status) { INT4 i, j, k; /* all purpose indices */ INT4 iSky = 0; /* for loop over sky positions */ INT4 jDeriv = 0; /* for loop over spindowns */ INT4 testNumber = 0; /* which test is being run */ /* The input and output structs used by the functions being tested */ PulsarSignalParams *pPulsarSignalParams = NULL; REAL4TimeSeries *signalvec = NULL; SFTParams *pSFTParams = NULL; SkyConstAndZeroPsiAMResponse *pSkyConstAndZeroPsiAMResponse; SFTandSignalParams *pSFTandSignalParams; SFTVector *outputSFTs = NULL; SFTVector *fastOutputSFTs = NULL; REAL4 renorm; /* to renormalize SFTs to account for different effective sample rates */ /* containers for detector and ephemeris data */ LALDetector cachedDetector; CHAR IFO[6] = "LHO"; EphemerisData *edat = NULL; char earthFile[] = TEST_DATA_DIR "earth00-19-DE405.dat.gz"; char sunFile[] = TEST_DATA_DIR "sun00-19-DE405.dat.gz"; /* containers for sky position and spindown data */ REAL8 **skyPosData; REAL8 **freqDerivData; INT4 numSkyPosTotal = 8; INT4 numSpinDown = 1; INT4 numFreqDerivTotal = 2; INT4 numFreqDerivIncludingNoSpinDown; /* set below */ INT4 iFreq = 0; REAL8 cosTmpDEC; REAL8 tmpDeltaRA; REAL8 DeltaRA = 0.01; REAL8 DeltaDec = 0.01; REAL8 DeltaFDeriv1 = -1.0e-10; REAL8 DeltaFDeriv2 = 0.0; REAL8 DeltaFDeriv3 = 0.0; REAL8 DeltaFDeriv4 = 0.0; REAL8 DeltaFDeriv5 = 0.0; /* containers for number of SFTs, times to generate SFTs, reference time, and gap */ INT4 numSFTs = 4; REAL8 f0SFT = 943.12; REAL8 bandSFT = 0.5; UINT4 gpsStartTimeSec = 765432109; /* GPS start time of data requested seconds; example = Apr 08 2004 04:01:36 UTC */ UINT4 gpsStartTimeNan = 0; /* GPS start time of data requested nanoseconds */ LIGOTimeGPSVector *timeStamps; LIGOTimeGPS GPSin; /* holder for reference-time for pulsar parameters at the detector; will convert to SSB! */ REAL8 sftGap = 73.0; /* extra artificial gap between SFTs */ REAL8 tSFT = 1800.0; REAL8 duration = tSFT + (numSFTs - 1)*(tSFT + sftGap); INT4 nBinsSFT = (INT4)(bandSFT*tSFT + 0.5); /* additional parameters that determine what signals to test; note f0SGNL and bandSGNL must be compatible with f0SFT and bandSFT */ REAL8 f0SGNL = f0SFT + bandSFT/2.0; REAL8 dfSGNL = 1.0/tSFT; REAL8 bandSGNL = 0.005; INT4 nBinsSGNL = (INT4)(bandSGNL*tSFT + 0.5); INT4 nsamples = 18000; /* nsamples from SFT header; 2 x this would be the effective number of time samples used to create an SFT from raw data */ INT4 Dterms = 3; /* 09/07/05 gam; use Dterms to fill in SFT bins with fake data as per LALDemod else fill in bin with zero */ REAL8 h_0 = 7.0e-22; /* Source amplitude; use arbitrary small number for default */ REAL8 cosIota; /* cosine of inclination angle iota of the source */ /* variables for comparing differences */ REAL4 maxDiffSFTMod, diffAtMaxPower; REAL4 overallMaxDiffSFTMod; REAL4 maxMod, fastMaxMod; INT4 jMaxMod, jFastMaxMod; REAL4 tmpDiffSFTMod, sftMod, fastSFTMod; REAL4 smallMod = 1.e-30; REAL4 epsDiffMod; REAL4 epsBinErrorRate; /* 10/12/04 gam; Allowed bin error rate */ INT4 binErrorCount = 0; /* 10/12/04 gam; Count number of bin errors */ /* randval is always set to a default value or given a random value to generate certain signal parameters or mismatch */ REAL4 randval; #ifdef INCLUDE_RANDVAL_MISMATCH INT4 seed=0; INT4 rndCount; RandomParams *randPar=NULL; FILE *fpRandom; #endif INITSTATUS(status); ATTATCHSTATUSPTR(status); /* generate timeStamps */ timeStamps = (LIGOTimeGPSVector *)LALMalloc(sizeof(LIGOTimeGPSVector)); timeStamps->data =(LIGOTimeGPS *)LALMalloc (numSFTs*sizeof(LIGOTimeGPS)); timeStamps->length = numSFTs; for (i = 0; i < numSFTs; i++) { timeStamps->data[i].gpsSeconds = gpsStartTimeSec + (UINT4)(i*(tSFT + sftGap)); timeStamps->data[i].gpsNanoSeconds = 0; } /* for i < numSFTs */ /* generate skyPosData */ skyPosData=(REAL8 **)LALMalloc(numSkyPosTotal*sizeof(REAL8 *)); for(iSky=0;iSky<numSkyPosTotal;iSky++) { skyPosData[iSky] = (REAL8 *)LALMalloc(2*sizeof(REAL8)); /* Try fairly random sky positions skyPosData[iSky][0] = RA, skyPosData[iSky][1] = DEC */ if (iSky == 0) { skyPosData[iSky][0] = 0.02*LAL_TWOPI; skyPosData[iSky][1] = 0.03*LAL_PI/2.0; } else if (iSky == 1) { skyPosData[iSky][0] = 0.23*LAL_TWOPI; skyPosData[iSky][1] = -0.57*LAL_PI/2.0; } else if (iSky == 2) { skyPosData[iSky][0] = 0.47*LAL_TWOPI; skyPosData[iSky][1] = 0.86*LAL_PI/2.0; } else if (iSky == 3) { skyPosData[iSky][0] = 0.38*LAL_TWOPI; skyPosData[iSky][1] = -0.07*LAL_PI/2.0; } else if (iSky == 4) { skyPosData[iSky][0] = 0.65*LAL_TWOPI; skyPosData[iSky][1] = 0.99*LAL_PI/2.0; } else if (iSky == 5) { skyPosData[iSky][0] = 0.72*LAL_TWOPI; skyPosData[iSky][1] = -0.99*LAL_PI/2.0; } else if (iSky == 6) { skyPosData[iSky][0] = 0.81*LAL_TWOPI; skyPosData[iSky][1] = 0.19*LAL_PI/2.0; } else if (iSky == 7) { skyPosData[iSky][0] = 0.99*LAL_TWOPI; skyPosData[iSky][1] = 0.01*LAL_PI/2.0; } else { skyPosData[iSky][0] = 0.0; skyPosData[iSky][1] = 0.0; } /* END if (k == 0) ELSE ... */ } /* END for(iSky=0;iSky<numSkyPosTotal;iSky++) */ freqDerivData = NULL; /* 02/02/05 gam */ if (numSpinDown > 0) { freqDerivData=(REAL8 **)LALMalloc(numFreqDerivTotal*sizeof(REAL8 *)); for(jDeriv=0;jDeriv<numFreqDerivTotal;jDeriv++) { freqDerivData[jDeriv] = (REAL8 *)LALMalloc(numSpinDown*sizeof(REAL8)); if (jDeriv == 0) { for(k=0;k<numSpinDown;k++) { freqDerivData[jDeriv][k] = 0.0; } } else if (jDeriv == 1) { for(k=0;k<numSpinDown;k++) { freqDerivData[jDeriv][k] = 1.e-9; /* REALLY THIS IS ONLY GOOD FOR numSpinDown = 1; */ } } else { for(k=0;k<numSpinDown;k++) { freqDerivData[jDeriv][k] = 0.0; } } /* END if (k == 0) ELSE ... */ } /* END for(jDeriv=0;jDeriv<numFreqDerivTotal;jDeriv++) */ numFreqDerivIncludingNoSpinDown = numFreqDerivTotal; } else { numFreqDerivIncludingNoSpinDown = 1; /* Even if numSpinDown = 0 still need to count case of zero spindown. */ } /* END if (numSpinDown > 0) ELSE ... */ /* Initialize ephemeris data */ XLAL_CHECK_LAL (status, ( edat = XLALInitBarycenter( earthFile, sunFile ) ) != NULL, XLAL_EFUNC); /* Allocate memory for PulsarSignalParams and initialize */ pPulsarSignalParams = (PulsarSignalParams *)LALMalloc(sizeof(PulsarSignalParams)); pPulsarSignalParams->pulsar.position.system = COORDINATESYSTEM_EQUATORIAL; pPulsarSignalParams->pulsar.spindown = NULL; if (numSpinDown > 0) { LALDCreateVector(status->statusPtr, &(pPulsarSignalParams->pulsar.spindown),((UINT4)numSpinDown)); CHECKSTATUSPTR (status); } pPulsarSignalParams->orbit.asini = 0 /* isolated pulsar */; pPulsarSignalParams->transfer = NULL; pPulsarSignalParams->dtDelayBy2 = 0; pPulsarSignalParams->dtPolBy2 = 0; /* Set up pulsar site */ if (strstr(IFO, "LHO")) { cachedDetector = lalCachedDetectors[LALDetectorIndexLHODIFF]; } else if (strstr(IFO, "LLO")) { cachedDetector = lalCachedDetectors[LALDetectorIndexLLODIFF]; } else if (strstr(IFO, "GEO")) { cachedDetector = lalCachedDetectors[LALDetectorIndexGEO600DIFF]; } else if (strstr(IFO, "VIRGO")) { cachedDetector = lalCachedDetectors[LALDetectorIndexVIRGODIFF]; } else if (strstr(IFO, "TAMA")) { cachedDetector = lalCachedDetectors[LALDetectorIndexTAMA300DIFF]; } else { /* "Invalid or null IFO" */ ABORT( status, GENERATEPULSARSIGNALTESTC_EIFO, GENERATEPULSARSIGNALTESTC_MSGEIFO); } pPulsarSignalParams->site = &cachedDetector; pPulsarSignalParams->ephemerides = edat; pPulsarSignalParams->startTimeGPS.gpsSeconds = (INT4)gpsStartTimeSec; pPulsarSignalParams->startTimeGPS.gpsNanoSeconds = (INT4)gpsStartTimeNan; pPulsarSignalParams->duration = (UINT4)duration; pPulsarSignalParams->samplingRate = (REAL8)ceil(2.0*bandSFT); /* Make sampleRate an integer so that T*samplingRate = integer for integer T */ pPulsarSignalParams->fHeterodyne = f0SFT; GPSin.gpsSeconds = timeStamps->data[0].gpsSeconds; GPSin.gpsNanoSeconds = timeStamps->data[0].gpsNanoSeconds; /* Allocate memory for SFTParams and initialize */ pSFTParams = (SFTParams *)LALCalloc(1, sizeof(SFTParams)); pSFTParams->Tsft = tSFT; pSFTParams->timestamps = timeStamps; pSFTParams->noiseSFTs = NULL; #ifdef INCLUDE_RANDVAL_MISMATCH /* Initial seed and randPar to use LALUniformDeviate to generate random mismatch during Monte Carlo. */ fpRandom = fopen("/dev/urandom","r"); rndCount = fread(&seed, sizeof(INT4),1, fpRandom); fclose(fpRandom); /* seed = 1234; */ /* Test value */ LALCreateRandomParams(status->statusPtr, &randPar, seed); CHECKSTATUSPTR (status); #endif /* allocate memory for structs needed by LALComputeSkyAndZeroPsiAMResponse and LALFastGeneratePulsarSFTs */ pSkyConstAndZeroPsiAMResponse = (SkyConstAndZeroPsiAMResponse *)LALMalloc(sizeof(SkyConstAndZeroPsiAMResponse)); pSkyConstAndZeroPsiAMResponse->skyConst = (REAL8 *)LALMalloc((2*numSpinDown*(numSFTs+1)+2*numSFTs+3)*sizeof(REAL8)); pSkyConstAndZeroPsiAMResponse->fPlusZeroPsi = (REAL4 *)LALMalloc(numSFTs*sizeof(REAL4)); pSkyConstAndZeroPsiAMResponse->fCrossZeroPsi = (REAL4 *)LALMalloc(numSFTs*sizeof(REAL4)); pSFTandSignalParams = (SFTandSignalParams *)LALMalloc(sizeof(SFTandSignalParams)); /* create lookup table (LUT) values for doing trig */ /* pSFTandSignalParams->resTrig = 64; */ /* length sinVal and cosVal; resolution of trig functions = 2pi/resTrig */ /* pSFTandSignalParams->resTrig = 128; */ /* 10/08/04 gam; length sinVal and cosVal; domain = -2pi to 2pi inclusive; resolution = 4pi/resTrig */ pSFTandSignalParams->resTrig = 0; /* 10/12/04 gam; turn off using LUTs since this is more typical. */ /* 02/02/05 gam; if NOT pSFTandSignalParams->resTrig > 0 should not create trigArg etc... */ if (pSFTandSignalParams->resTrig > 0) { pSFTandSignalParams->trigArg = (REAL8 *)LALMalloc((pSFTandSignalParams->resTrig+1)*sizeof(REAL8)); pSFTandSignalParams->sinVal = (REAL8 *)LALMalloc((pSFTandSignalParams->resTrig+1)*sizeof(REAL8)); pSFTandSignalParams->cosVal = (REAL8 *)LALMalloc((pSFTandSignalParams->resTrig+1)*sizeof(REAL8)); for (k=0; k<=pSFTandSignalParams->resTrig; k++) { /* pSFTandSignalParams->trigArg[k]= ((REAL8)LAL_TWOPI) * ((REAL8)k) / ((REAL8)pSFTandSignalParams->resTrig); */ /* 10/08/04 gam */ pSFTandSignalParams->trigArg[k]= -1.0*((REAL8)LAL_TWOPI) + 2.0 * ((REAL8)LAL_TWOPI) * ((REAL8)k) / ((REAL8)pSFTandSignalParams->resTrig); pSFTandSignalParams->sinVal[k]=sin( pSFTandSignalParams->trigArg[k] ); pSFTandSignalParams->cosVal[k]=cos( pSFTandSignalParams->trigArg[k] ); } } pSFTandSignalParams->pSigParams = pPulsarSignalParams; pSFTandSignalParams->pSFTParams = pSFTParams; pSFTandSignalParams->nSamples = nsamples; pSFTandSignalParams->Dterms = Dterms; /* 09/07/05 gam */ /* ********************************************************/ /* */ /* START SECTION: LOOP OVER SKY POSITIONS */ /* */ /* ********************************************************/ for(iSky=0;iSky<numSkyPosTotal;iSky++) { /* set source sky position declination (DEC) */ randval = 0.5; /* Gives default value */ #ifdef INCLUDE_SEQUENTIAL_MISMATCH randval = ( (REAL4)(iSky) )/( (REAL4)(numSkyPosTotal) ); #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif pPulsarSignalParams->pulsar.position.latitude = skyPosData[iSky][1] + (((REAL8)randval) - 0.5)*DeltaDec; cosTmpDEC = cos(skyPosData[iSky][1]); if (cosTmpDEC != 0.0) { tmpDeltaRA = DeltaRA/cosTmpDEC; } else { tmpDeltaRA = 0.0; /* We are at a celestial pole */ } /* set source sky position right ascension (RA) */ randval = 0.5; /* Gives default value */ #ifdef INCLUDE_SEQUENTIAL_MISMATCH randval = ( (REAL4)(iSky) )/( (REAL4)(numSkyPosTotal) ); #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif pPulsarSignalParams->pulsar.position.longitude = skyPosData[iSky][0] + (((REAL8)randval) - 0.5)*tmpDeltaRA; /* Find reference time in SSB for this sky positions */ int ret = XLALConvertGPS2SSB ( &(pPulsarSignalParams->pulsar.refTime), GPSin, pPulsarSignalParams ); if ( ret != XLAL_SUCCESS ) { XLALPrintError ("XLALConvertGPS2SSB() failed with xlalErrno = %d\n", xlalErrno ); ABORTXLAL (status); } /* one per sky position fill in SkyConstAndZeroPsiAMResponse for use with LALFastGeneratePulsarSFTs */ LALComputeSkyAndZeroPsiAMResponse (status->statusPtr, pSkyConstAndZeroPsiAMResponse, pSFTandSignalParams); CHECKSTATUSPTR (status); /* ********************************************************/ /* */ /* START SECTION: LOOP OVER SPINDOWN */ /* */ /* ********************************************************/ for(jDeriv=0;jDeriv<numFreqDerivIncludingNoSpinDown;jDeriv++) { /* source spindown parameters */ if (numSpinDown > 0) { for(k=0;k<numSpinDown;k++) { randval = 0.5; /* Gives default value */ #ifdef INCLUDE_SEQUENTIAL_MISMATCH if (freqDerivData[jDeriv][k] < 0.0) { randval = ( (REAL4)(iSky) )/( (REAL4)(numSkyPosTotal) ); } else { randval = 0.5; /* If derivative is not negative (i.e., it is zero) then do not add in mismatch; keep it zero. */ } #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif if (k == 0) { pPulsarSignalParams->pulsar.spindown->data[k] = freqDerivData[jDeriv][k] + (((REAL8)randval) - 0.5)*DeltaFDeriv1; } else if (k == 1) { pPulsarSignalParams->pulsar.spindown->data[k] = freqDerivData[jDeriv][k] + (((REAL8)randval) - 0.5)*DeltaFDeriv2; } else if (k == 2) { pPulsarSignalParams->pulsar.spindown->data[k] = freqDerivData[jDeriv][k] + (((REAL8)randval) - 0.5)*DeltaFDeriv3; } else if (k == 3) { pPulsarSignalParams->pulsar.spindown->data[k] = freqDerivData[jDeriv][k] + (((REAL8)randval) - 0.5)*DeltaFDeriv4; } else if (k == 4) { pPulsarSignalParams->pulsar.spindown->data[k] = freqDerivData[jDeriv][k] + (((REAL8)randval) - 0.5)*DeltaFDeriv5; } /* END if (k == 0) ELSE ... */ } } /* ***************************************************/ /* */ /* START SECTION: LOOP OVER FREQUENCIES */ /* */ /* ***************************************************/ for(iFreq=0;iFreq<nBinsSGNL;iFreq++) { /* set source orientation psi */ randval = 0.5; /* Gives default value */ #ifdef INCLUDE_SEQUENTIAL_MISMATCH randval = ( (REAL4)(iFreq) )/( (REAL4)(nBinsSGNL) ); #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif pPulsarSignalParams->pulsar.psi = (randval - 0.5) * ((REAL4)LAL_PI_2); /* set angle between source spin axis and direction from source to SSB, cosIota */ randval = 1.0; /* Gives default value */ #ifdef INCLUDE_SEQUENTIAL_MISMATCH randval = ( (REAL4)(iFreq) )/( (REAL4)(nBinsSGNL) ); #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif cosIota = 2.0*((REAL8)randval) - 1.0; /* h_0 is fixed above; get A_+ and A_x from h_0 and cosIota */ pPulsarSignalParams->pulsar.aPlus = (REAL4)(0.5*h_0*(1.0 + cosIota*cosIota)); pPulsarSignalParams->pulsar.aCross = (REAL4)(h_0*cosIota); /* get random value for phi0 */ randval = 0.125; /* Gives default value pi/4*/ #ifdef INCLUDE_SEQUENTIAL_MISMATCH randval = ( (REAL4)(iFreq) )/( (REAL4)(nBinsSGNL) ); #endif #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif pPulsarSignalParams->pulsar.phi0 = ((REAL8)randval) * ((REAL8)LAL_TWOPI); /* randval steps through various mismatches to frequencies centered on a bin */ /* Note that iFreq = nBinsSGNL/2 gives randval = 0.5 gives no mismatch from frequency centered on a bin */ randval = ( (REAL4)(iFreq) )/( (REAL4)(nBinsSGNL) ); #ifdef INCLUDE_RANDVAL_MISMATCH LALUniformDeviate(status->statusPtr, &randval, randPar); CHECKSTATUSPTR (status); #endif pPulsarSignalParams->pulsar.f0 = f0SGNL + iFreq*dfSGNL + (((REAL8)randval) - 0.5)*dfSGNL; testNumber++; /* Update count of which test we about to do. */ /* FIRST: Use LALGeneratePulsarSignal and LALSignalToSFTs to generate outputSFTs */ signalvec = NULL; LALGeneratePulsarSignal(status->statusPtr, &signalvec, pPulsarSignalParams); CHECKSTATUSPTR (status); outputSFTs = NULL; LALSignalToSFTs(status->statusPtr, &outputSFTs, signalvec, pSFTParams); CHECKSTATUSPTR (status); #ifdef PRINT_OUTPUTSFT if (testNumber == TESTNUMBER_TO_PRINT) { i=SFTINDEX_TO_PRINT; /* index of which outputSFT to output */ fprintf(stdout,"iFreq = %i, inject h_0 = %23.10e \n",iFreq,h_0); fprintf(stdout,"iFreq = %i, inject cosIota = %23.10e, A_+ = %23.10e, A_x = %23.10e \n",iFreq,cosIota,pPulsarSignalParams->pulsar.aPlus,pPulsarSignalParams->pulsar.aCross); fprintf(stdout,"iFreq = %i, inject psi = %23.10e \n",iFreq,pPulsarSignalParams->pulsar.psi); fprintf(stdout,"iFreq = %i, inject phi0 = %23.10e \n",iFreq,pPulsarSignalParams->pulsar.phi0); fprintf(stdout,"iFreq = %i, search f0 = %23.10e, inject f0 = %23.10e \n",iFreq,f0SGNL + iFreq*dfSGNL,pPulsarSignalParams->pulsar.f0); fprintf(stdout,"outputSFTs->data[%i].data->length = %i \n",i,outputSFTs->data[i].data->length); renorm = ((REAL4)nsamples)/((REAL4)(outputSFTs->data[i].data->length - 1)); for(j=0;j<nBinsSFT;j++) { /* fprintf(stdout,"%i %g %g \n", j, renorm*outputSFTs->data[i].data->data[j].re, renorm*outputSFTs->data[i].data->data[j].im); */ fprintf(stdout,"%i %g \n",j,renorm*renorm*outputSFTs->data[i].data->data[j].re*outputSFTs->data[i].data->data[j].re + renorm*renorm*outputSFTs->data[i].data->data[j].im*outputSFTs->data[i].data->data[j].im); fflush(stdout); } } #endif /* SECOND: Use LALComputeSkyAndZeroPsiAMResponse and LALFastGeneratePulsarSFTs to generate outputSFTs */ LALFastGeneratePulsarSFTs (status->statusPtr, &fastOutputSFTs, pSkyConstAndZeroPsiAMResponse, pSFTandSignalParams); CHECKSTATUSPTR (status); #ifdef PRINT_FASTOUTPUTSFT if (testNumber == TESTNUMBER_TO_PRINT) { REAL4 fPlus; REAL4 fCross; i=SFTINDEX_TO_PRINT; /* index of which outputSFT to output */ fPlus = pSkyConstAndZeroPsiAMResponse->fPlusZeroPsi[i]*cos(2.0*pPulsarSignalParams->pulsar.psi) + pSkyConstAndZeroPsiAMResponse->fCrossZeroPsi[i]*sin(2.0*pPulsarSignalParams->pulsar.psi); fCross = pSkyConstAndZeroPsiAMResponse->fCrossZeroPsi[i]*cos(2.0*pPulsarSignalParams->pulsar.psi) - pSkyConstAndZeroPsiAMResponse->fPlusZeroPsi[i]*sin(2.0*pPulsarSignalParams->pulsar.psi); fprintf(stdout,"iFreq = %i, inject h_0 = %23.10e \n",iFreq,h_0); fprintf(stdout,"iFreq = %i, inject cosIota = %23.10e, A_+ = %23.10e, A_x = %23.10e \n",iFreq,cosIota,pPulsarSignalParams->pulsar.aPlus,pPulsarSignalParams->pulsar.aCross); fprintf(stdout,"iFreq = %i, inject psi = %23.10e \n",iFreq,pPulsarSignalParams->pulsar.psi); fprintf(stdout,"iFreq = %i, fPlus, fCross = %23.10e, %23.10e \n",iFreq,fPlus,fCross); fprintf(stdout,"iFreq = %i, inject phi0 = %23.10e \n",iFreq,pPulsarSignalParams->pulsar.phi0); fprintf(stdout,"iFreq = %i, search f0 = %23.10e, inject f0 = %23.10e \n",iFreq,f0SGNL + iFreq*dfSGNL,pPulsarSignalParams->pulsar.f0); fprintf(stdout,"fastOutputSFTs->data[%i].data->length = %i \n",i,fastOutputSFTs->data[i].data->length); fflush(stdout); for(j=0;j<nBinsSFT;j++) { /* fprintf(stdout,"%i %g %g \n",j,fastOutputSFTs->data[i].data->data[j].re,fastOutputSFTs->data[i].data->data[j].im); */ fprintf(stdout,"%i %g \n",j,fastOutputSFTs->data[i].data->data[j].re*fastOutputSFTs->data[i].data->data[j].re + fastOutputSFTs->data[i].data->data[j].im*fastOutputSFTs->data[i].data->data[j].im); fflush(stdout); } } #endif /* find maximum difference in power */ epsDiffMod = 0.20; /* maximum allowed percent difference */ /* 10/12/04 gam */ overallMaxDiffSFTMod = 0.0; for (i = 0; i < numSFTs; i++) { renorm = ((REAL4)nsamples)/((REAL4)(outputSFTs->data[i].data->length - 1)); maxDiffSFTMod = 0.0; diffAtMaxPower = 0.0; maxMod = 0.0; fastMaxMod = 0.0; jMaxMod = -1; jFastMaxMod = -1; /* Since doppler shifts can move the signal by an unknown number of bins search the whole band for max modulus: */ for(j=0;j<nBinsSFT;j++) { sftMod = renorm*renorm*crealf(outputSFTs->data[i].data->data[j])*crealf(outputSFTs->data[i].data->data[j]) + renorm*renorm*cimagf(outputSFTs->data[i].data->data[j])*cimagf(outputSFTs->data[i].data->data[j]); sftMod = sqrt(sftMod); fastSFTMod = crealf(fastOutputSFTs->data[i].data->data[j])*crealf(fastOutputSFTs->data[i].data->data[j]) + cimagf(fastOutputSFTs->data[i].data->data[j])*cimagf(fastOutputSFTs->data[i].data->data[j]); fastSFTMod = sqrt(fastSFTMod); if (fabs(sftMod) > smallMod) { tmpDiffSFTMod = fabs((sftMod - fastSFTMod)/sftMod); if (tmpDiffSFTMod > maxDiffSFTMod) { maxDiffSFTMod = tmpDiffSFTMod; } if (tmpDiffSFTMod > overallMaxDiffSFTMod) { overallMaxDiffSFTMod = tmpDiffSFTMod; } if (sftMod > maxMod) { maxMod = sftMod; jMaxMod = j; } if (fastSFTMod > fastMaxMod) { fastMaxMod = fastSFTMod; jFastMaxMod = j; } } } if (fabs(maxMod) > smallMod) { diffAtMaxPower = fabs((maxMod - fastMaxMod)/maxMod); } #ifdef PRINT_MAXSFTPOWER fprintf(stdout,"maxSFTMod, testNumber %i, SFT %i, bin %i = %g \n",testNumber,i,jMaxMod,maxMod); fprintf(stdout,"maxFastSFTMod, testNumber %i, SFT %i, bin %i = %g \n",testNumber,i,jFastMaxMod,fastMaxMod); fflush(stdout); #endif #ifdef PRINT_MAXDIFFSFTPOWER fprintf(stdout,"maxDiffSFTMod, testNumber %i, SFT %i, bin %i = %g \n",testNumber,i,jMaxDiff,maxDiffSFTMod); fflush(stdout); #endif #ifdef PRINT_DIFFATMAXPOWER fprintf(stdout,"diffAtMaxPower, testNumber %i, SFT %i, bins %i and %i = %g \n",testNumber,i,jMaxMod,jFastMaxMod,diffAtMaxPower); fflush(stdout); #endif #ifdef PRINT_ERRORATMAXPOWER if (diffAtMaxPower > epsDiffMod) { fprintf(stdout,"diffAtMaxPower, testNumber %i, SFT %i, bins %i and %i = %g exceeded epsDiffMod = %g \n",testNumber,i,jMaxMod,jFastMaxMod,diffAtMaxPower,epsDiffMod); fflush(stdout); /* break; */ /* only report 1 error per test */ } if (jMaxMod != jFastMaxMod) { fprintf(stdout,"MaxPower occurred in different bins: testNumber %i, SFT %i, bins %i and %i\n",testNumber,i,jMaxMod,jFastMaxMod); fflush(stdout); /* break; */ /* only report 1 error per test */ } #endif if (jMaxMod != jFastMaxMod) { binErrorCount++; /* 10/12/04 gam; count up bin errors; if too ABORT at bottom of code */ } if ( diffAtMaxPower > epsDiffMod ) { ABORT( status, GENERATEPULSARSIGNALTESTC_EMOD, GENERATEPULSARSIGNALTESTC_MSGEMOD); } /* 10/12/04 gam; turn on test above and add test below */ if ( fabs(((REAL8)(jMaxMod - jFastMaxMod))) > 1.1 ) { ABORT( status, GENERATEPULSARSIGNALTESTC_EBIN, GENERATEPULSARSIGNALTESTC_MSGEBIN); } } /* END for(i = 0; i < numSFTs; i++) */ #ifdef PRINT_OVERALLMAXDIFFSFTPOWER fprintf(stdout,"overallMaxDiffSFTMod, testNumber = %i, SFT %i, bin %i = %g \n",testNumber,iOverallMaxDiffSFTMod,jOverallMaxDiff,overallMaxDiffSFTMod); fflush(stdout); #endif /* 09/07/05 gam; Initialize fastOutputSFTs since only 2*Dterms bins are changed by LALFastGeneratePulsarSFTs */ for (i = 0; i < numSFTs; i++) { for(j=0;j<nBinsSFT;j++) { fastOutputSFTs->data[i].data->data[j] = 0.0; } } XLALDestroySFTVector( outputSFTs); LALFree(signalvec->data->data); LALFree(signalvec->data); LALFree(signalvec); } /* END for(iFreq=0;iFreq<nBinsSGNL;iFreq++) */ /* ***************************************************/ /* */ /* END SECTION: LOOP OVER FREQUENCIES */ /* */ /* ***************************************************/ } /* END for(jDeriv=0;jDeriv<numFreqDerivIncludingNoSpinDown;jDeriv++) */ /* ********************************************************/ /* */ /* END SECTION: LOOP OVER SPINDOWN */ /* */ /* ********************************************************/ } /* END for(iSky=0;iSky<numSkyPosTotal;iSky++) */ /* ********************************************************/ /* */ /* END SECTION: LOOP OVER SKY POSITIONS */ /* */ /* ********************************************************/ /* 10/12/04 gam; check if too many bin errors */ epsBinErrorRate = 0.20; /* 10/12/04 gam; maximum allowed bin errors */ if ( (((REAL4)binErrorCount)/((REAL4)testNumber)) > epsBinErrorRate ) { ABORT( status, GENERATEPULSARSIGNALTESTC_EBINS, GENERATEPULSARSIGNALTESTC_MSGEBINS); } #ifdef INCLUDE_RANDVAL_MISMATCH LALDestroyRandomParams(status->statusPtr, &randPar); CHECKSTATUSPTR (status); #endif /* fprintf(stdout,"Total number of tests completed = %i. \n", testNumber); fflush(stdout); */ LALFree(pSFTParams); if (numSpinDown > 0) { LALDDestroyVector(status->statusPtr, &(pPulsarSignalParams->pulsar.spindown)); CHECKSTATUSPTR (status); } LALFree(pPulsarSignalParams); /* deallocate memory for structs needed by LALComputeSkyAndZeroPsiAMResponse and LALFastGeneratePulsarSFTs */ XLALDestroySFTVector( fastOutputSFTs); LALFree(pSkyConstAndZeroPsiAMResponse->fCrossZeroPsi); LALFree(pSkyConstAndZeroPsiAMResponse->fPlusZeroPsi); LALFree(pSkyConstAndZeroPsiAMResponse->skyConst); LALFree(pSkyConstAndZeroPsiAMResponse); /* 02/02/05 gam; if NOT pSFTandSignalParams->resTrig > 0 should not create trigArg etc... */ if (pSFTandSignalParams->resTrig > 0) { LALFree(pSFTandSignalParams->trigArg); LALFree(pSFTandSignalParams->sinVal); LALFree(pSFTandSignalParams->cosVal); } LALFree(pSFTandSignalParams); /* deallocate skyPosData */ for(i=0;i<numSkyPosTotal;i++) { LALFree(skyPosData[i]); } LALFree(skyPosData); if (numSpinDown > 0) { /* deallocate freqDerivData */ for(i=0;i<numFreqDerivTotal;i++) { LALFree(freqDerivData[i]); } LALFree(freqDerivData); } LALFree(timeStamps->data); LALFree(timeStamps); XLALDestroyEphemerisData(edat); CHECKSTATUSPTR (status); DETATCHSTATUSPTR (status); }
/*--------------main function---------------*/ int main(int argc, char **argv){ const CHAR *fn = __func__; InputParams XLAL_INIT_DECL(inputs); REAL8 srate = 16384.0; /*sample rate defaulted to 16384 */ /* read in command line input args */ ReadInput( &inputs, argc, argv ); LALStatus XLAL_INIT_DECL(status); EphemerisData *edat; if ( (edat = InitEphemeris ( inputs.ephemType, inputs.ephemDir)) == NULL ){ XLALPrintError ( "%s: Failed to init ephemeris data\n", fn ); XLAL_ERROR ( XLAL_EFUNC ); } /*init detector info */ LALDetector *site; if ( ( site = XLALGetSiteInfo ( inputs.det )) == NULL ){ XLALPrintError("%s: Failed to get site-info for detector '%s'\n", fn, inputs.det ); XLAL_ERROR ( XLAL_EFUNC ); } if( inputs.geocentre ){ /* set site to the geocentre */ site->location[0] = 0.0; site->location[1] = 0.0; site->location[2] = 0.0; } struct dirent **pulsars; INT4 n=scandir(inputs.pulsarDir, &pulsars, 0, alphasort); if ( n < 0){ XLALPrintError("scandir failed\n"); XLAL_ERROR(XLAL_EIO); } UINT4 numpulsars = (UINT4)n; UINT4 h=0; CHAR parname[256]; PulsarParameters *pulparams[numpulsars]; for(h=2; h<numpulsars; h++){ if(strstr(pulsars[h]->d_name,".par") == NULL){ free(pulsars[h]); continue; } else{ sprintf(parname,"%s/%s", inputs.pulsarDir, pulsars[h]->d_name); fprintf(stderr, "%s\n", parname); FILE *inject; if (( inject = fopen ( parname, "r" )) == NULL ){ fprintf(stderr,"Error opening file: %s\n", parname); XLAL_ERROR ( XLAL_EIO ); } pulparams[h] = XLALReadTEMPOParFile( parname ); fclose( inject ); } } LIGOTimeGPS epoch; UINT4 ndata; epoch.gpsSeconds = inputs.epoch; epoch.gpsNanoSeconds = 0; ndata = inputs.frDur; REAL8TimeSeries *series=NULL; CHAR out_file[256]; sprintf(out_file, "%s-%s-%d-%d.gwf", inputs.det, inputs.outStr, epoch.gpsSeconds, ndata ); LALFrameH *outFrame = NULL; if ((outFrame = XLALFrameNew( &epoch, (REAL8)ndata, inputs.channel, 1, 0, 0 )) == NULL) { LogPrintf(LOG_CRITICAL, "%s : XLALFrameNew() filed with error = %d.\n", fn, xlalErrno); XLAL_ERROR( XLAL_EFAILED); } if ((series = XLALCreateREAL8TimeSeries( inputs.channel, &epoch, 0., 1./srate,&lalSecondUnit, (int)(ndata*srate) )) == NULL) { XLAL_ERROR( XLAL_EFUNC ); } UINT4 counter=0; for (counter = 0; counter < series->data->length; counter++) series->data->data[counter] = 0; /*** Read Pulsar Data ***/ for (h=0; h < numpulsars; h++){ if(strstr(pulsars[h]->d_name,".par")==NULL){ free(pulsars[h]); continue; } else{ PulsarSignalParams XLAL_INIT_DECL(params); /* set signal generation barycenter delay look-up table step size */ params.dtDelayBy2 = 10.; /* generate table every 10 seconds */ if (( params.pulsar.spindown = XLALCreateREAL8Vector(1)) == NULL ){ XLALPrintError("Out of memory"); XLAL_ERROR ( XLAL_EFUNC ); } INT4 dtpos = 0; if ( PulsarCheckParam(pulparams[h], "POSEPOCH") ) dtpos = epoch.gpsSeconds - (INT4)PulsarGetREAL8Param(pulparams[h], "POSEPOCH"); else dtpos = epoch.gpsSeconds - (INT4)PulsarGetREAL8Param(pulparams[h], "PEPOCH"); REAL8 ra = 0., dec = 0.; if ( PulsarCheckParam( pulparams[h], "RAJ" ) ) { ra = PulsarGetREAL8Param( pulparams[h], "RAJ" ); } else if ( PulsarCheckParam( pulparams[h], "RA" ) ){ ra = PulsarGetREAL8Param( pulparams[h], "RA" ); } else{ XLALPrintError("No right ascension found"); XLAL_ERROR ( XLAL_EFUNC ); } if ( PulsarCheckParam( pulparams[h], "DECJ" ) ) { dec = PulsarGetREAL8Param( pulparams[h], "DECJ" ); } else if ( PulsarCheckParam( pulparams[h], "DEC" ) ){ dec = PulsarGetREAL8Param( pulparams[h], "DEC" ); } else{ XLALPrintError("No declination found"); XLAL_ERROR ( XLAL_EFUNC ); } params.pulsar.position.latitude = dec + (REAL8)dtpos * PulsarGetREAL8ParamOrZero(pulparams[h], "PMDEC"); params.pulsar.position.longitude = ra + (REAL8)dtpos * PulsarGetREAL8ParamOrZero(pulparams[h], "PMRA") / cos(params.pulsar.position.latitude); params.pulsar.position.system = COORDINATESYSTEM_EQUATORIAL; REAL8Vector *fs = PulsarGetREAL8VectorParam(pulparams[h], "F"); if ( fs->length == 0 ){ XLALPrintError("No frequencies found"); XLAL_ERROR ( XLAL_EFUNC ); } params.pulsar.f0 = 2.*fs->data[0]; if ( fs->length > 1 ){ params.pulsar.spindown->data[0] = 2.*fs->data[1]; } if (( XLALGPSSetREAL8(&(params.pulsar.refTime), PulsarGetREAL8Param(pulparams[h], "PEPOCH")) ) == NULL ) XLAL_ERROR ( XLAL_EFUNC ); params.pulsar.psi = PulsarGetREAL8ParamOrZero(pulparams[h], "PSI"); params.pulsar.phi0 = PulsarGetREAL8ParamOrZero(pulparams[h], "PHI0"); REAL8 cosiota = PulsarGetREAL8ParamOrZero(pulparams[h], "COSIOTA"); REAL8 h0 = PulsarGetREAL8ParamOrZero(pulparams[h], "H0"); params.pulsar.aPlus = 0.5 * h0 * (1. + cosiota * cosiota ); params.pulsar.aCross = h0 * cosiota; /*Add binary later if needed!*/ params.site = site; params.ephemerides = edat; params.startTimeGPS = epoch; params.duration = ndata; params.samplingRate = srate; params.fHeterodyne = 0.; REAL4TimeSeries *TSeries = NULL; LALGeneratePulsarSignal( &status, &TSeries, ¶ms ); if (status.statusCode){ fprintf(stderr, "LAL Routine failed!\n"); XLAL_ERROR (XLAL_EFAILED); } UINT4 i; for (i=0; i < TSeries->data->length; i++) series->data->data[i] += TSeries->data->data[i]; XLALDestroyREAL4TimeSeries(TSeries); XLALDestroyREAL8Vector(params.pulsar.spindown); } } if (XLALFrameAddREAL8TimeSeriesProcData(outFrame,series)){ LogPrintf(LOG_CRITICAL, "%s : XLALFrameAddREAL8TimeSeries() failed with error = %d.\n",fn,xlalErrno); XLAL_ERROR(XLAL_EFAILED); } CHAR OUTFILE[256]; sprintf(OUTFILE, "%s/%s", inputs.outDir, out_file); if ( XLALFrameWrite(outFrame, OUTFILE)){ LogPrintf(LOG_CRITICAL, "%s : XLALFrameWrite() failed with error = %d.\n", fn, xlalErrno); XLAL_ERROR( XLAL_EFAILED ); } XLALFrameFree(outFrame); XLALDestroyREAL8TimeSeries( series ); return 0; }
int main( int argc, char *argv[]){ static LALStatus status; static LALDetector detector; static LIGOTimeGPSVector timeV; static REAL8Cart3CoorVector velV; static REAL8Vector timeDiffV; static REAL8Vector foft; static PulsarSignalParams params; static SFTParams sftParams; static UCHARPeakGram pg1; static COMPLEX8SFTData1 sft1; static REAL8PeriodoPSD periPSD; REAL4TimeSeries *signalTseries = NULL; SFTVector *inputSFTs = NULL; SFTVector *outputSFTs = NULL; /* data about injected signal */ static PulsarData pulsarInject; /* the template */ static HoughTemplate pulsarTemplate, pulsarTemplate1; /*FILE *fpOUT = NULL; output file pointer */ FILE *fpLog = NULL; /* log file pointer */ CHAR *logstr=NULL; /* log string containing user input variables */ CHAR *fnamelog=NULL; /* name of log file */ INT4 nfSizeCylinder; EphemerisData *edat = NULL; INT4 mObsCoh, numberCount; REAL8 sftBand; REAL8 timeBase, deltaF, normalizeThr, threshold; UINT4 sftlength; INT4 sftFminBin; REAL8 fHeterodyne; REAL8 tSamplingRate; /* grid spacings */ REAL8 deltaTheta; INT4 mmP, mmT; /* for loop over mismatched templates */ /* user input variables */ INT4 uvar_ifo, uvar_blocksRngMed; REAL8 uvar_peakThreshold; REAL8 uvar_alpha, uvar_delta, uvar_h0, uvar_f0; REAL8 uvar_psi, uvar_phi0, uvar_fdot, uvar_cosiota; CHAR *uvar_earthEphemeris=NULL; CHAR *uvar_sunEphemeris=NULL; CHAR *uvar_sftDir=NULL; CHAR *uvar_fnameout=NULL; /* set up the default parameters */ nfSizeCylinder = NFSIZE; /* set other user input variables */ uvar_peakThreshold = THRESHOLD; uvar_ifo = IFO; uvar_blocksRngMed = BLOCKSRNGMED; /* set default pulsar parameters */ uvar_h0 = H0; uvar_alpha = ALPHA; uvar_delta = DELTA; uvar_f0 = F0; uvar_fdot = FDOT; uvar_psi = PSI; uvar_cosiota = COSIOTA; uvar_phi0 = PHI0; /* now set the default filenames */ uvar_earthEphemeris = (CHAR *)LALMalloc(512*sizeof(CHAR)); strcpy(uvar_earthEphemeris,EARTHEPHEMERIS); uvar_sunEphemeris = (CHAR *)LALMalloc(512*sizeof(CHAR)); strcpy(uvar_sunEphemeris,SUNEPHEMERIS); uvar_sftDir = (CHAR *)LALMalloc(512*sizeof(CHAR)); strcpy(uvar_sftDir,SFTDIRECTORY); uvar_fnameout = (CHAR *)LALMalloc(512*sizeof(CHAR)); strcpy(uvar_fnameout, FILEOUT); /* register user input variables */ XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_ifo, "ifo", INT4, 'i', OPTIONAL, "Detector GEO(1) LLO(2) LHO(3)" ) == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_blocksRngMed, "blocksRngMed", INT4, 'w', OPTIONAL, "RngMed block size") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_peakThreshold, "peakThreshold", REAL8, 't', OPTIONAL, "Peak selection threshold") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_earthEphemeris, "earthEphemeris", STRING, 'E', OPTIONAL, "Earth Ephemeris file") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_sunEphemeris, "sunEphemeris", STRING, 'S', OPTIONAL, "Sun Ephemeris file") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_sftDir, "sftDir", STRING, 'D', OPTIONAL, "SFT Directory") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_fnameout, "fnameout", STRING, 'o', OPTIONAL, "Output file prefix") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_alpha, "alpha", REAL8, 'r', OPTIONAL, "Right ascension") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_delta, "delta", REAL8, 'l', OPTIONAL, "Declination") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_h0, "h0", REAL8, 'm', OPTIONAL, "h0 to inject") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_f0, "f0", REAL8, 'f', OPTIONAL, "Start search frequency") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_psi, "psi", REAL8, 'p', OPTIONAL, "Polarization angle") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_phi0, "phi0", REAL8, 'P', OPTIONAL, "Initial phase") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_cosiota, "cosiota", REAL8, 'c', OPTIONAL, "Cosine of iota") == XLAL_SUCCESS, XLAL_EFUNC); XLAL_CHECK_MAIN( XLALRegisterNamedUvar( &uvar_fdot, "fdot", REAL8, 'd', OPTIONAL, "Spindown parameter") == XLAL_SUCCESS, XLAL_EFUNC); /* read all command line variables */ BOOLEAN should_exit = 0; XLAL_CHECK_MAIN( XLALUserVarReadAllInput(&should_exit, argc, argv, lalAppsVCSInfoList) == XLAL_SUCCESS, XLAL_EFUNC); if (should_exit) exit(1); /* write the log file */ fnamelog = (CHAR *)LALMalloc( 512*sizeof(CHAR)); strcpy(fnamelog, uvar_fnameout); strcat(fnamelog, "_log"); /* open the log file for writing */ if ((fpLog = fopen(fnamelog, "w")) == NULL) { fprintf(stderr, "Unable to open file %s for writing\n", fnamelog); LALFree(fnamelog); exit(1); } /* get the log string */ XLAL_CHECK_MAIN( ( logstr = XLALUserVarGetLog(UVAR_LOGFMT_CFGFILE) ) != NULL, XLAL_EFUNC); fprintf( fpLog, "## Log file for HoughMismatch\n\n"); fprintf( fpLog, "# User Input:\n"); fprintf( fpLog, "#-------------------------------------------\n"); fprintf( fpLog, "%s", logstr); LALFree(logstr); /* append an ident-string defining the exact CVS-version of the code used */ { CHAR command[1024] = ""; fprintf (fpLog, "\n\n# CVS-versions of executable:\n"); fprintf (fpLog, "# -----------------------------------------\n"); fclose (fpLog); sprintf (command, "ident %s | sort -u >> %s", argv[0], fnamelog); /* we don't check this. If it fails, we assume that */ /* one of the system-commands was not available, and */ /* therefore the CVS-versions will not be logged */ if ( system(command) ) fprintf (stderr, "\nsystem('%s') returned non-zero status!\n\n", command ); LALFree(fnamelog); } /* set peak selection threshold */ SUB( LALRngMedBias( &status, &normalizeThr, uvar_blocksRngMed ), &status ); threshold = uvar_peakThreshold/normalizeThr; /* set detector */ if (uvar_ifo ==1) detector=lalCachedDetectors[LALDetectorIndexGEO600DIFF]; if (uvar_ifo ==2) detector=lalCachedDetectors[LALDetectorIndexLLODIFF]; if (uvar_ifo ==3) detector=lalCachedDetectors[LALDetectorIndexLHODIFF]; /* copy user input values */ pulsarInject.f0 = uvar_f0; pulsarInject.latitude = uvar_delta; pulsarInject.longitude = uvar_alpha; pulsarInject.aPlus = 0.5 * uvar_h0 * ( 1.0 + uvar_cosiota * uvar_cosiota ); pulsarInject.aCross = uvar_h0 * uvar_cosiota; pulsarInject.psi = uvar_psi; pulsarInject.phi0 = uvar_phi0; pulsarInject.spindown.length = 1; pulsarInject.spindown.data = NULL; pulsarInject.spindown.data = (REAL8 *)LALMalloc(sizeof(REAL8)); pulsarInject.spindown.data[0] = uvar_fdot; /* copy these values also to the pulsar template */ /* template is complately matched at this point */ pulsarTemplate.f0 = uvar_f0; pulsarTemplate.latitude = uvar_delta; pulsarTemplate.longitude = uvar_alpha; pulsarTemplate.spindown.length = 1; pulsarTemplate.spindown.data = NULL; pulsarTemplate.spindown.data = (REAL8 *)LALMalloc(sizeof(REAL8)); pulsarTemplate.spindown.data[0] = uvar_fdot; /* allocate memory for mismatched spindown template */ pulsarTemplate1.spindown.length = 1; pulsarTemplate1.spindown.data = NULL; pulsarTemplate1.spindown.data = (REAL8 *)LALMalloc(sizeof(REAL8)); /* read sfts */ { CHAR *tempDir; tempDir = (CHAR *)LALMalloc(512*sizeof(CHAR)); strcpy(tempDir, uvar_sftDir); strcat(tempDir, "/*SFT*.*"); sftBand = 0.5; SUB( LALReadSFTfiles ( &status, &inputSFTs, uvar_f0 - sftBand, uvar_f0 + sftBand, nfSizeCylinder + uvar_blocksRngMed , tempDir), &status); LALFree(tempDir); } /* get sft parameters */ mObsCoh = inputSFTs->length; sftlength = inputSFTs->data->data->length; deltaF = inputSFTs->data->deltaF; timeBase = 1.0/deltaF; sftFminBin = floor( timeBase * inputSFTs->data->f0 + 0.5); fHeterodyne = sftFminBin*deltaF; tSamplingRate = 2.0*deltaF*(sftlength -1.); /* create timestamp vector */ timeV.length = mObsCoh; timeV.data = NULL; timeV.data = (LIGOTimeGPS *)LALMalloc(mObsCoh*sizeof(LIGOTimeGPS)); /* read timestamps */ { INT4 i; SFTtype *sft= NULL; sft = inputSFTs->data; for (i=0; i < mObsCoh; i++){ timeV.data[i].gpsSeconds = sft->epoch.gpsSeconds; timeV.data[i].gpsNanoSeconds = sft->epoch.gpsNanoSeconds; ++sft; } } /* compute the time difference relative to startTime for all SFT */ timeDiffV.length = mObsCoh; timeDiffV.data = NULL; timeDiffV.data = (REAL8 *)LALMalloc(mObsCoh*sizeof(REAL8)); { REAL8 t0, ts, tn, midTimeBase; INT4 j; midTimeBase=0.5*timeBase; ts = timeV.data[0].gpsSeconds; tn = timeV.data[0].gpsNanoSeconds * 1.00E-9; t0=ts+tn; timeDiffV.data[0] = midTimeBase; for (j=1; j< mObsCoh; ++j){ ts = timeV.data[j].gpsSeconds; tn = timeV.data[j].gpsNanoSeconds * 1.00E-9; timeDiffV.data[j] = ts+tn -t0+midTimeBase; } } /* compute detector velocity for those time stamps */ velV.length = mObsCoh; velV.data = NULL; velV.data = (REAL8Cart3Coor *)LALMalloc(mObsCoh*sizeof(REAL8Cart3Coor)); { VelocityPar velPar; REAL8 vel[3]; UINT4 j; velPar.detector = detector; velPar.tBase = timeBase; velPar.vTol = 0.0; /* irrelevant */ velPar.edat = NULL; /* read in ephemeris data */ XLAL_CHECK_MAIN( ( edat = XLALInitBarycenter( uvar_earthEphemeris, uvar_sunEphemeris ) ) != NULL, XLAL_EFUNC); velPar.edat = edat; /* calculate detector velocity */ for(j=0; j< velV.length; ++j){ velPar.startTime.gpsSeconds = timeV.data[j].gpsSeconds; velPar.startTime.gpsNanoSeconds = timeV.data[j].gpsNanoSeconds; SUB( LALAvgDetectorVel ( &status, vel, &velPar), &status ); velV.data[j].x= vel[0]; velV.data[j].y= vel[1]; velV.data[j].z= vel[2]; } } /* set grid spacings */ { deltaTheta = 1.0 / ( VTOT * uvar_f0 * timeBase ); /* currently unused: REAL8 deltaFdot = deltaF / timeBase; */ } /* allocate memory for f(t) pattern */ foft.length = mObsCoh; foft.data = NULL; foft.data = (REAL8 *)LALMalloc(mObsCoh*sizeof(REAL8)); /* allocate memory for Hough peripsd structure */ periPSD.periodogram.length = sftlength; periPSD.periodogram.data = NULL; periPSD.periodogram.data = (REAL8 *)LALMalloc(sftlength* sizeof(REAL8)); periPSD.psd.length = sftlength; periPSD.psd.data = NULL; periPSD.psd.data = (REAL8 *)LALMalloc(sftlength* sizeof(REAL8)); /* allocate memory for peakgram */ pg1.length = sftlength; pg1.data = NULL; pg1.data = (UCHAR *)LALMalloc(sftlength* sizeof(UCHAR)); /* generate signal and add to input sfts */ /* parameters for output sfts */ sftParams.Tsft = timeBase; sftParams.timestamps = &(timeV); sftParams.noiseSFTs = inputSFTs; /* signal generation parameters */ params.orbit.asini = 0 /* isolated pulsar */; /* params.transferFunction = NULL; */ params.site = &(detector); params.ephemerides = edat; params.startTimeGPS.gpsSeconds = timeV.data[0].gpsSeconds; /* start time of output time series */ params.startTimeGPS.gpsNanoSeconds = timeV.data[0].gpsNanoSeconds; /* start time of output time series */ params.duration = timeDiffV.data[mObsCoh-1] + 0.5 * timeBase; /* length of time series in seconds */ params.samplingRate = tSamplingRate; params.fHeterodyne = fHeterodyne; /* reference time for frequency and spindown is first timestamp */ params.pulsar.refTime.gpsSeconds = timeV.data[0].gpsSeconds; params.pulsar.refTime.gpsNanoSeconds = timeV.data[0].gpsNanoSeconds; params.pulsar.position.longitude = pulsarInject.longitude; params.pulsar.position.latitude = pulsarInject.latitude ; params.pulsar.position.system = COORDINATESYSTEM_EQUATORIAL; params.pulsar.psi = pulsarInject.psi; params.pulsar.aPlus = pulsarInject.aPlus; params.pulsar.aCross = pulsarInject.aCross; params.pulsar.phi0 = pulsarInject.phi0; params.pulsar.f0 = pulsarInject.f0; params.pulsar.spindown = &pulsarInject.spindown ; SUB( LALGeneratePulsarSignal(&status, &signalTseries, ¶ms ), &status); SUB( LALSignalToSFTs(&status, &outputSFTs, signalTseries, &sftParams), &status); /* fill in elements of sft structure sft1 used in peak selection */ sft1.length = sftlength; sft1.fminBinIndex = sftFminBin; sft1.timeBase = timeBase; /* loop over mismatched templates */ for (mmT = -2; mmT <= 2; mmT++) { for (mmP = -2; mmP <= 2; mmP++) { INT4 mmFactor; /* displace the template */ mmFactor = 1.0; pulsarTemplate1.f0 = pulsarTemplate.f0 /*+ mmFactor * mm * deltaF*/; pulsarTemplate1.latitude = pulsarTemplate.latitude + mmFactor * mmT * deltaTheta; pulsarTemplate1.longitude = pulsarTemplate.longitude + mmFactor * mmP * deltaTheta; pulsarTemplate1.spindown.data[0] = pulsarTemplate.spindown.data[0] /*+ mmFactor * mm * deltaFdot*/; numberCount = 0; /* now calculate the number count for the template */ INT4 j; for (j=0; j < mObsCoh; j++) { INT4 ind; sft1.epoch.gpsSeconds = timeV.data[j].gpsSeconds; sft1.epoch.gpsNanoSeconds = timeV.data[j].gpsNanoSeconds; sft1.data = outputSFTs->data[j].data->data; SUB( COMPLEX8SFT2Periodogram1(&status, &periPSD.periodogram, &sft1), &status ); SUB( LALPeriodo2PSDrng( &status, &periPSD.psd, &periPSD.periodogram, &uvar_blocksRngMed), &status ); SUB( LALSelectPeakColorNoise(&status,&pg1,&threshold,&periPSD), &status); SUB( ComputeFoft(&status, &foft, &pulsarTemplate1, &timeDiffV, &velV, timeBase), &status); ind = floor( foft.data[j]*timeBase - sftFminBin + 0.5); numberCount += pg1.data[ind]; } /* print the number count */ fprintf(stdout, "%d %d %d\n", mmT, mmP, numberCount); } } /* free structures created by signal generation routines */ LALFree(signalTseries->data->data); LALFree(signalTseries->data); LALFree(signalTseries); signalTseries =NULL; XLALDestroySFTVector( outputSFTs); /* destroy input sfts */ XLALDestroySFTVector( inputSFTs); /* free other structures */ LALFree(foft.data); LALFree(pulsarInject.spindown.data); LALFree(pulsarTemplate.spindown.data); LALFree(pulsarTemplate1.spindown.data); LALFree(timeV.data); LALFree(timeDiffV.data); LALFree(velV.data); XLALDestroyEphemerisData(edat); LALFree(periPSD.periodogram.data); LALFree(periPSD.psd.data); LALFree(pg1.data); XLALDestroyUserVars(); LALCheckMemoryLeaks(); INFO( DRIVEHOUGHCOLOR_MSGENORM ); return DRIVEHOUGHCOLOR_ENORM; }